H2Lib  3.0
bem2d.h
Go to the documentation of this file.
1 
2 /* ------------------------------------------------------------
3  * This is the file "bem2d.h" of the H2Lib package.
4  * All rights reserved, Sven Christophersen 2011
5  * ------------------------------------------------------------ */
6 
13 #ifndef BEM2D_H_
14 #define BEM2D_H_
15 
16 /* C STD LIBRARY */
17 #include <stdio.h>
18 /* CORE 0 */
19 #include "settings.h"
20 /* CORE 1 */
21 #include "amatrix.h"
22 #include "gaussquad.h"
23 #include "factorizations.h"
24 /* CORE 2 */
25 #include "cluster.h"
26 #include "clustergeometry.h"
27 #include "hmatrix.h"
28 #include "clusterbasis.h"
29 #include "h2matrix.h"
30 /* CORE 3 */
31 #include "hcoarsen.h"
32 #include "h2update.h"
33 #include "aca.h"
34 /* SIMPLE */
35 /* PARTICLES */
36 /* BEM */
37 #include "curve2d.h"
38 #include "singquad1d.h"
39 
54 /* ------------------------------------------------------------ *
55  * Type definitions *
56  * ------------------------------------------------------------ */
57 
62 typedef struct _bem2d bem2d;
66 typedef bem2d *pbem2d;
70 typedef const bem2d* pcbem2d;
71 
79 typedef struct _aprxbem2d aprxbem2d;
87 typedef const aprxbem2d *pcaprxbem2d;
88 
94 typedef struct _parbem2d parbem2d;
102 typedef const parbem2d *pcparbem2d;
103 
109 typedef struct _kernelbem2d kernelbem2d;
110 
115 
119 typedef const kernelbem2d *pckernelbem2d;
120 
140 typedef void (*quadpoints2d)(pcbem2d bem, const real bmin[2],
141  const real bmax[2], const real delta, real (**Z)[2], real (**N)[2]);
142 
165 };
166 
171 
185 typedef field (*boundary_func2d)(const real *x, const real *n);
186 
205 struct _bem2d {
213 
223 
228 
233 
238 
243 
253 
258 
267 
275 
298 void (*nearfield)(const uint *ridx, const uint *cidx, pcbem2d bem, bool ntrans,
299  pamatrix N);
300 
321 void (*farfield_rk)(pccluster rc, uint rname, pccluster cc, uint cname,
322  pcbem2d bem, prkmatrix R);
323 
344 void (*farfield_u)(uint rname, uint cname, pcbem2d bem, puniform U);
345 
365 void (*leaf_row)(pcbem2d bem, pclusterbasis rb, uint rname);
366 
386 void (*leaf_col)(pcbem2d bem, pclusterbasis cb, uint cname);
387 
408 void (*transfer_row)(pcbem2d bem, pclusterbasis rb, uint rname);
409 
430 void (*transfer_col)(pcbem2d bem, pclusterbasis cb, uint cname);
431 
440 
449 
467 };
468 
486 struct _kernelbem2d {
507 void (*fundamental)(const real (*X)[2], const real (*Y)[2], pamatrix V);
508 
533 void (*dny_fundamental)(const real (*X)[2], const real (*Y)[2],
534  const real (*NY)[2], pamatrix V);
535 
562 void (*dnx_dny_fundamental)(const real (*X)[2], const real (*NX)[2],
563  const real (*Y)[2], const real (*NY)[2], pamatrix V);
564 
592 void (*kernel_row)(const uint *idx, const real (*Z)[2], pcbem2d bem, pamatrix A);
593 
622 void (*kernel_col)(const uint *idx, const real (*Z)[2], pcbem2d bem, pamatrix B);
623 
655 void (*dnz_kernel_row)(const uint *idx, const real (*Z)[2], const real (*NZ)[2],
656  pcbem2d bem, pamatrix A);
657 
689 void (*dnz_kernel_col)(const uint *idx, const real (*Z)[2], const real (*NZ)[2],
690  pcbem2d bem, pamatrix B);
691 
715 void (*fundamental_row)(const uint *idx, const real (*Z)[2], pcbem2d bem,
716  pamatrix A);
717 
742 void (*fundamental_col)(const uint *idx, const real (*Z)[2], pcbem2d bem,
743  pamatrix B);
744 
772 void (*dnz_fundamental_row)(const uint *idx, const real (*Z)[2],
773  const real (*NZ)[2], pcbem2d bem, pamatrix A);
774 
802 void (*dnz_fundamental_col)(const uint *idx, const real (*Z)[2],
803  const real (*NZ)[2], pcbem2d bem, pamatrix B);
804 
835 void (*lagrange_row)(const uint *idx, pcavector px, pcavector py, pcbem2d bem,
836  pamatrix V);
837 
868 void (*lagrange_col)(const uint *idx, pcavector px, pcavector py, pcbem2d bem,
869  pamatrix W);
870 };
871 
872 /* ------------------------------------------------------------
873  * Constructors and destructors
874  * ------------------------------------------------------------ */
875 
888 
895 HEADER_PREFIX void del_bem2d(pbem2d bem);
896 
897 /* ------------------------------------------------------------
898  * Methods to build cluster trees
899  * ------------------------------------------------------------ */
900 
919  uint **idx);
920 
940  uint **idx);
941 
966  uint **idx, basisfunctionbem2d basis);
967 
989  basisfunctionbem2d basis);
990 
991 /* ------------------------------------------------------------
992  * Initializer functions for h-matrix approximations
993  * ------------------------------------------------------------ */
994 
1025 HEADER_PREFIX void setup_hmatrix_recomp_bem2d(pbem2d bem, bool recomp,
1026  real accur_recomp, bool coarsen, real accur_coarsen);
1027 
1028 /* ------------------------------------------------------------
1029  * Interpolation
1030  * ------------------------------------------------------------ */
1031 
1062  pccluster cc, pcblock tree, uint m);
1063 
1094  pccluster cc, pcblock tree, uint m);
1095 
1120  pccluster rc, pccluster cc, pcblock tree, uint m);
1121 
1122 /* ------------------------------------------------------------
1123  * Green
1124  * ------------------------------------------------------------ */
1125 
1187  pccluster cc, pcblock tree, uint m, uint l, real delta,
1188  quadpoints2d quadpoints);
1189 
1253  pccluster cc, pcblock tree, uint m, uint l, real delta,
1254  quadpoints2d quadpoints);
1255 
1293  pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta,
1294  quadpoints2d quadpoints);
1295 
1296 /* ------------------------------------------------------------
1297  * Greenhybrid
1298  * ------------------------------------------------------------ */
1299 
1358  pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta,
1359  real accur, quadpoints2d quadpoints);
1360 
1419  pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta,
1420  real accur, quadpoints2d quadpoints);
1421 
1462  pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta,
1463  real accur, quadpoints2d quadpoints);
1464 
1465 /* ------------------------------------------------------------
1466  * ACA
1467  * ------------------------------------------------------------ */
1468 
1486  pccluster cc, pcblock tree, real accur);
1487 
1505  pccluster cc, pcblock tree, real accur);
1506 
1507 /* ------------------------------------------------------------
1508  * HCA
1509  * ------------------------------------------------------------ */
1510 
1570  pccluster cc, pcblock tree, uint m, real accur);
1571 
1572 /* ------------------------------------------------------------
1573  * Initializer functions for h2-matrix approximations
1574  * ------------------------------------------------------------ */
1575 
1592 HEADER_PREFIX void setup_h2matrix_recomp_bem2d(pbem2d bem, bool hiercomp,
1593  real accur_hiercomp);
1594 
1648  pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m);
1649 
1726  pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l,
1727  real delta, real accur, quadpoints2d quadpoints);
1728 
1805  pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l,
1806  real delta, real accur, quadpoints2d quadpoints);
1807 
1808 /* ------------------------------------------------------------
1809  * Fill hmatrix
1810  * ------------------------------------------------------------ */
1811 
1834 
1864  phmatrix G);
1865 
1866 /* ------------------------------------------------------------
1867  * Fill h2-matrix
1868  * ------------------------------------------------------------ */
1869 
1896  pclusterbasis rb);
1897 
1924  pclusterbasis cb);
1925 
1954 
1987  ph2matrix G);
1988 
1989 /* ------------------------------------------------------------
1990  * Lagrange polynomials
1991  * ------------------------------------------------------------ */
1992 
2024  pcavector px, pcavector py, pcbem2d bem, pamatrix V);
2025 
2057  pcavector px, pcavector py, pcbem2d bem, pamatrix V);
2058 
2086  pcavector px, pcavector py, pamatrix V);
2087 
2088 /* ------------------------------------------------------------
2089  * some useful functions
2090  * ------------------------------------------------------------ */
2091 
2114  field (*rhs)(const real *x, const real *n), pavector f);
2115 
2150  const real b[2], const real delta, real (**Z)[2], real (**N)[2]);
2151 
2186  const real b[2], const real delta, real (**Z)[2], real (**N)[2]);
2187 
2211  void* data);
2212 
2230  void* data);
2231 
2234 #endif /* BEM2D_H_ */
void setup_hmatrix_aprx_inter_row_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
This function initializes the bem2d object for approximating a matrix with tensorinterpolation within...
void assemblecoarsen_bem2d_hmatrix(pbem2d bem, pblock b, phmatrix G)
Fills an hmatrix with a predefined approximation technique using coarsening strategy.
void setup_h2matrix_aprx_greenhybrid_bem2d(pbem2d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
Initialize the bem2d object for approximating a h2matrix with green&#39;s method and ACA based recompress...
real * mass
Mass-matrix of reference basis functions.
Definition: bem2d.h:252
field(* boundary_func2d)(const real *x, const real *n)
Definition: bem2d.h:185
void assemble_bem2d_h2matrix_col_clusterbasis(pcbem2d bem, pclusterbasis cb)
This function computes the matrix entries for the nested clusterbasis .
void assemble_bem2d_dn_lagrange_const_amatrix(const uint *idx, pcavector px, pcavector py, pcbem2d bem, pamatrix V)
This function will integrate the normal derivatives of the Lagrange polynomials on the boundary domai...
void(* leaf_row)(pcbem2d bem, pclusterbasis rb, uint rname)
Computes the the matrix for a leaf cluster .
Definition: bem2d.h:365
pclustergeometry build_bem2d_const_clustergeometry(pcbem2d bem, uint **idx)
Creates a clustergeometry object for a BEM-Problem using piecewise constant basis functions...
pbem2d new_bem2d(pccurve2d gr)
Main constructor for bem2d objects.
Definition: avector.h:39
const aprxbem2d * pcaprxbem2d
Definition: bem2d.h:87
This struct collects all type of quadrature formulas needed by the computation of matrix entries with...
Definition: singquad1d.h:39
void build_bem2d_rect_quadpoints(pcbem2d bem, const real a[2], const real b[2], const real delta, real(**Z)[2], real(**N)[2])
Generating quadrature points, weights and normal vectors on a rectangular parameterization.
void setup_hmatrix_aprx_aca_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, real accur)
Approximate matrix block with ACA using full pivoting.
void(* quadpoints2d)(pcbem2d bem, const real bmin[2], const real bmax[2], const real delta, real(**Z)[2], real(**N)[2])
Callback function type for parameterizations.
Definition: bem2d.h:140
uint N_neumann
Number of degrees of freedom for neumann data.
Definition: bem2d.h:227
void(* farfield_rk)(pccluster rc, uint rname, pccluster cc, uint cname, pcbem2d bem, prkmatrix R)
Computes rank-k-approximations of a given block.
Definition: bem2d.h:321
Representation of a polygon in 2D.
Definition: curve2d.h:46
Representation of -matrices.
Definition: h2matrix.h:48
void assemble_bem2d_hmatrix(pbem2d bem, pblock b, phmatrix G)
Fills an hmatrix with a predefined approximation technique.
pclustergeometry build_bem2d_clustergeometry(pcbem2d bem, uint **idx, basisfunctionbem2d basis)
Creates a clustergeometry object for a BEM-Problem using the type of basis functions specified by bas...
void assemblehiercomp_bem2d_h2matrix(pbem2d bem, pblock b, ph2matrix G)
Fills an h2matrix with a predefined approximation technique using hierarchical recompression.
Representation of -matrices.
Definition: hmatrix.h:49
Representation of a cluster basis.
Definition: clusterbasis.h:40
enum _basisfunctionbem2d basisfunctionbem2d
Definition: bem2d.h:170
Enum value to represent piecewise constant basis functions.
Definition: bem2d.h:160
const parbem2d * pcparbem2d
Definition: bem2d.h:102
uint N_dirichlet
Number of degrees of freedom for dirichlet data.
Definition: bem2d.h:237
void(* transfer_row)(pcbem2d bem, pclusterbasis rb, uint rname)
Computes the transfermatrices for a cluster .
Definition: bem2d.h:408
const kernelbem2d * pckernelbem2d
Definition: bem2d.h:119
struct _parbem2d parbem2d
Definition: bem2d.h:94
void setup_hmatrix_aprx_greenhybrid_col_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
creating hmatrix approximation using green&#39;s method with column cluster connected with recompression ...
unsigned uint
Unsigned integer type.
Definition: settings.h:70
void(* transfer_col)(pcbem2d bem, pclusterbasis cb, uint cname)
Computes the transfermatrices for a cluster .
Definition: bem2d.h:430
pamatrix build_bem2d_amatrix(pccluster row, pccluster col, void *data)
For a block defined by row and col this function will compute the full submatrix and return an amatri...
void projectl2_bem2d_const_avector(pbem2d bem, field(*rhs)(const real *x, const real *n), pavector f)
Computes the -projection of a given function using piecewise constant basis functions.
double _Complex field
Field type.
Definition: settings.h:171
void build_bem2d_circle_quadpoints(pcbem2d bem, const real a[2], const real b[2], const real delta, real(**Z)[2], real(**N)[2])
Generating quadrature points, weights and normal vectors on a circular parameterization.
basisfunctionbem2d basis_dirichlet
Type of basis function for dirichlet-data.
Definition: bem2d.h:242
void setup_hmatrix_aprx_greenhybrid_mixed_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
creating hmatrix approximation using green&#39;s method with row or column cluster connected with recompr...
Representation of a clustergeometry object.
Definition: clustergeometry.h:45
void setup_hmatrix_recomp_bem2d(pbem2d bem, bool recomp, real accur_recomp, bool coarsen, real accur_coarsen)
Initialize the bem object for on the fly recompression techniques.
uint * v2ts
Vertex to triangle map stating indices.
Definition: bem2d.h:274
Main container object for computation of BEM related matrices and vectors.
Definition: bem2d.h:205
aprxbem2d * paprxbem2d
Definition: bem2d.h:83
pparbem2d par
Some helpers to make parallel computations possible.
Definition: bem2d.h:448
void setup_h2matrix_recomp_bem2d(pbem2d bem, bool hiercomp, real accur_hiercomp)
Enables hierarchical recompression for hmatrices.
basisfunctionbem2d basis_neumann
Type of basis function for neumann-data.
Definition: bem2d.h:232
struct _aprxbem2d aprxbem2d
Definition: bem2d.h:79
void setup_h2matrix_aprx_inter_bem2d(pbem2d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m)
Initialize the bem2d object for approximating a h2matrix with tensorinterpolation.
pcluster build_bem2d_cluster(pcbem2d bem, uint clf, basisfunctionbem2d basis)
Creates a clustertree for specified basis functions.
void(* leaf_col)(pcbem2d bem, pclusterbasis cb, uint cname)
Computes the the matrix for a leaf cluster .
Definition: bem2d.h:386
psingquad1d sq
Quadrature rules used within BEM computation.
Definition: bem2d.h:222
void setup_h2matrix_aprx_greenhybrid_ortho_bem2d(pbem2d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
Initialize the bem2d object for approximating a h2matrix with green&#39;s method and ACA based recompress...
pccurve2d gr
Polygonal two dimensional mesh.
Definition: bem2d.h:212
Substructure containing callback functions to different types of kernel evaluations.
Definition: bem2d.h:486
_basisfunctionbem2d
Possible types of basis functions.
Definition: bem2d.h:152
Enum value to represent piecewise linear basis functions.
Definition: bem2d.h:164
prkmatrix build_bem2d_rkmatrix(pccluster row, pccluster col, void *data)
For a block defined by row and col this function will build a rank-k-approximation and return the rkm...
field alpha
Boundary integral operator + mass matrix.
Definition: bem2d.h:257
paprxbem2d aprx
A collection of necessary data structures for approximating matrices with various techniques...
Definition: bem2d.h:439
void setup_hmatrix_aprx_paca_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, real accur)
Approximate matrix block with ACA using partial pivoting.
#define HEADER_PREFIX
Prefix for function declarations.
Definition: settings.h:43
void assemble_bem2d_lagrange_const_amatrix(const uint *idx, pcavector px, pcavector py, pcbem2d bem, pamatrix V)
This function will integrate Lagrange polynomials on the boundary domain using piecewise constant bas...
double real
real floating point type.
Definition: settings.h:97
bem2d * pbem2d
Definition: bem2d.h:66
kernelbem2d * pkernelbem2d
Definition: bem2d.h:114
void(* nearfield)(const uint *ridx, const uint *cidx, pcbem2d bem, bool ntrans, pamatrix N)
Computes nearfield entries of Galerkin matrices.
Definition: bem2d.h:298
void setup_hmatrix_aprx_hca_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, real accur)
Approximate matrix block with hybrid cross approximation using Interpolation and ACA with full pivoti...
Representation of block trees.
Definition: block.h:48
void assemble_bem2d_h2matrix(pbem2d bem, pblock b, ph2matrix G)
Fills an h2matrix with a predefined approximation technique.
Dummy value, should never be used.
Definition: bem2d.h:156
Representation of cluster trees.
Definition: cluster.h:40
void assemble_bem2d_h2matrix_row_clusterbasis(pcbem2d bem, pclusterbasis rb)
This function computes the matrix entries for the nested clusterbasis .
Representation of an admissible block for -matrices.
Definition: uniform.h:69
Representation of a low-rank matrix in factorized form .
Definition: rkmatrix.h:36
void setup_hmatrix_aprx_green_col_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints2d quadpoints)
creating hmatrix approximation using green&#39;s method with column cluster .
void setup_hmatrix_aprx_greenhybrid_row_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
creating hmatrix approximation using green&#39;s method with row cluster connected with recompression tec...
Representation of a matrix as an array in column-major order.
Definition: amatrix.h:43
void(* farfield_u)(uint rname, uint cname, pcbem2d bem, puniform U)
Computes coupling matrix for uniform matrices.
Definition: bem2d.h:344
void setup_hmatrix_aprx_inter_mixed_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
This function initializes the bem2d object for approximating a matrix with tensorinterpolation within...
uint * v2t
vertex to triangle map.
Definition: bem2d.h:266
void assemble_bem2d_lagrange_amatrix(const real(*X)[2], pcavector px, pcavector py, pamatrix V)
This function will evaluate the Lagrange polynomials in a given point set X and store the results in ...
void setup_hmatrix_aprx_green_mixed_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints2d quadpoints)
creating hmatrix approximation using green&#39;s method with one of row or column cluster ...
void del_bem2d(pbem2d bem)
Destructor for bem2d objects.
const bem2d * pcbem2d
Definition: bem2d.h:70
void setup_hmatrix_aprx_green_row_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints2d quadpoints)
creating hmatrix approximation using green&#39;s method with row cluster .
pclustergeometry build_bem2d_linear_clustergeometry(pcbem2d bem, uint **idx)
Creates a clustergeometry object for a BEM-Problem using linear basis functions.
void setup_hmatrix_aprx_inter_col_bem2d(pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
This function initializes the bem2d object for approximating a matrix with tensorinterpolation within...
parbem2d * pparbem2d
Definition: bem2d.h:98
pkernelbem2d kernels
A collection of callback functions for different types of kernel evaluations.
Definition: bem2d.h:466