H2Lib  3.0
Data Fields
_bem2d Struct Reference

Main container object for computation of BEM related matrices and vectors. More...

#include <bem2d.h>

Data Fields

pccurve2d gr
 Polygonal two dimensional mesh. More...
 
psingquad1d sq
 Quadrature rules used within BEM computation. More...
 
uint N_neumann
 Number of degrees of freedom for neumann data.
 
basisfunctionbem2d basis_neumann
 Type of basis function for neumann-data.
 
uint N_dirichlet
 Number of degrees of freedom for dirichlet data.
 
basisfunctionbem2d basis_dirichlet
 Type of basis function for dirichlet-data.
 
realmass
 Mass-matrix of reference basis functions. More...
 
field alpha
 Boundary integral operator + $\alpha$ mass matrix.
 
uintv2t
 vertex to triangle map. More...
 
uintv2ts
 Vertex to triangle map stating indices. More...
 
void(* nearfield )(const uint *ridx, const uint *cidx, pcbem2d bem, bool ntrans, pamatrix N)
 Computes nearfield entries of Galerkin matrices. More...
 
void(* farfield_rk )(pccluster rc, uint rname, pccluster cc, uint cname, pcbem2d bem, prkmatrix R)
 Computes rank-k-approximations of a given block. More...
 
void(* farfield_u )(uint rname, uint cname, pcbem2d bem, puniform U)
 Computes coupling matrix $ S_b $ for uniform matrices. More...
 
void(* leaf_row )(pcbem2d bem, pclusterbasis rb, uint rname)
 Computes the the matrix $ V_t $ for a leaf cluster $ t \in \mathcal L_{\mathcal I} $ . More...
 
void(* leaf_col )(pcbem2d bem, pclusterbasis cb, uint cname)
 Computes the the matrix $ W_s $ for a leaf cluster $ s \in \mathcal L_{\mathcal J} $ . More...
 
void(* transfer_row )(pcbem2d bem, pclusterbasis rb, uint rname)
 Computes the transfermatrices $ E_{t'} $ for a cluster $ t \in \mathcal T_{\mathcal I} \, , t' \in \operatorname{sons}(t) $ . More...
 
void(* transfer_col )(pcbem2d bem, pclusterbasis cb, uint cname)
 Computes the transfermatrices $ F_{s'} $ for a cluster $ s \in \mathcal T_{\mathcal J} \, , s' \in \operatorname{sons}(s) $ . More...
 
paprxbem2d aprx
 A collection of necessary data structures for approximating matrices with various techniques. More...
 
pparbem2d par
 Some helpers to make parallel computations possible. More...
 
pkernelbem2d kernels
 A collection of callback functions for different types of kernel evaluations. More...
 

Detailed Description

Main container object for computation of BEM related matrices and vectors.

This struct defines the essential part of all calculations concerning BEM. Just creating a simple bem2d object is not intended, though possible. Therefore call a problem specific constructor such as new_slp_laplace_bem2d or new_dlp_laplace_bem2d . This will initialize all needed callback functions to build up at least dense matrices. To be able to create h- or h2matrix approximations one needs to call one of the set_hmatrix_aprx* or set_h2matrix_aprx* depending on the type of your matrix, to fully initialize the bem2d object. Afterwards you can build up your desired matrix approximation with a few simple functions calls, delivered by the bem2d object, using your favored approximation technique. Calling another of these " set " functions will reinitialize the bem2d object and you can build another approximation with a different technique.

Field Documentation

paprxbem2d aprx

A collection of necessary data structures for approximating matrices with various techniques.

Necessary Objects for approximation matrices with different approximation schemes are stored within aprx.

void(* farfield_rk) (pccluster rc, uint rname, pccluster cc, uint cname, pcbem2d bem, prkmatrix R)

Computes rank-k-approximations of a given block.

This callback function computes a rank-k-approximation of a matrix block given by row cluster rc and column cluster cc. The Matrices A and B are stored in a rank-k-matrix R.

Parameters
rcCurrent row cluster .
rnameEnumerated number of current row cluster rc.
ccCurrent column cluster
cnameEnumerated number of current column cluster cc.
bemBEM-object containing additional information for computation of the matrix entries.
RRank-k-matrix, which takes up the resulting matrices A and B . The rank k is defined by the underlying approximation technique.
void(* farfield_u) (uint rname, uint cname, pcbem2d bem, puniform U)

Computes coupling matrix $ S_b $ for uniform matrices.

While using h2-matrices the coupling matrices $ S_b $ are computed by this callback function. Coupling matrices are stored in the struct uniform defining a certain block of the whole h2-matrix .

Parameters
rnameEnumerated number of current row cluster clusterbasis U->rb.
cnameEnumerated number of current column clusterbasis U->cb.
bemBEM-object containing additional information for computation of the matrix entries.
UAn Object of type uniform containing the coupling matrix $ S_b $ , the row clusterbasis U->rb and the column clusterbasis U->cb .

Polygonal two dimensional mesh.

Polygonal boundary mesh used for the BEM-problem.

See also
curve2d.
pkernelbem2d kernels

A collection of callback functions for different types of kernel evaluations.

This struct defines a set of callback functions that are based on the fundamental solution of the underlying problem. They are used to build up approximations that will be used within farfield_rk , farfield_u , leaf_row , leaf_col , transfer_row and transfer_col. In Order to receive correct result, the callback functions defined in kernels have to be set correctly to the given problem by a constructor. These callback functions also depend on the utilized basis functions for either neumann- and dirichlet-data. As an example the constructors new_slp_laplace_bem2d and new_dlp_laplace_bem2d should serve.

See also
kernelbem2d
void(* leaf_col) (pcbem2d bem, pclusterbasis cb, uint cname)

Computes the the matrix $ W_s $ for a leaf cluster $ s \in \mathcal L_{\mathcal J} $ .

leaf_col is used for building up a clusterbasis for the set $ \mathcal L_{\mathcal J} $ for h2matrices. For each $ s \in \mathcal L_{\mathcal J} $ leaf_col will construct the matrix $ W_s $ . In case of nested clusterbasis transfer_col has also to be set to a valid function constructing suitable transfer matrices $ F_s $ .

Parameters
bemBEM-object containing additional information for computation of the matrix entries.
cbcurrent column clusterbasis. The Matrix $ W_s $ will be saved into cb->V .
cnameEnumerated number of current column clusterbasis cb.
void(* leaf_row) (pcbem2d bem, pclusterbasis rb, uint rname)

Computes the the matrix $ V_t $ for a leaf cluster $ t \in \mathcal L_{\mathcal I} $ .

leaf_row is used for building up a clusterbasis for the set $ \mathcal L_{\mathcal I} $ for h2matrices. For each $ t \in \mathcal L_{\mathcal I}$ leaf_row will construct the matrix $ V_t $ . In case of nested clusterbasis transfer_row has also to be set to a valid function constructing suitable transfer matrices $ E_t $.

Parameters
bemBEM-object containing additional information for computation of the matrix entries.
rbcurrent row clusterbasis. The Matrix $ V_t $ will be saved into rb->V .
rnameEnumerated number of current row clusterbasis rb.
real* mass

Mass-matrix of reference basis functions.

Mass-matrix for all combinations of reference basis functions $ \widehat\varphi_i $ and $ \widehat\psi_j $ . Entries are computed as:

\[ M_{ij} = \int_\Gamma \widehat\varphi_i(x) \cdot \widehat\psi_j(x) \, \mathrm d x \]

void(* nearfield) (const uint *ridx, const uint *cidx, pcbem2d bem, bool ntrans, pamatrix N)

Computes nearfield entries of Galerkin matrices.

This callback function computes 'nearfield' entries of the underlying problem and integral-operator.

Parameters
ridxDefines the indices of row boundary elements used. If ridx equals NULL, then Elements 0, ..., N->rows-1 are used.
cidxDefines the indices of column boundary elements used. If cidx equals NULL, then Elements 0, ..., N->cols-1 are used.
bemContains additional Information for the computation of the matrix entries.
ntransIs a boolean flag to indicates the way of storing the entries in matrix N . If ntrans == true the matrix entries will be stored in a transposed way.
NFor ntrans == false the matrix entries are computed as:

\[ N_{ij} = \int_\Gamma \int_\Gamma \varphi_i(x) \, g(x,y) \, \psi_j(y) \, \mathrm d y \, \mathrm d x \, ,\]

otherwise they are stored as

\[ N_{ji} = \int_\Gamma \int_\Gamma \varphi_i(x) \, g(x,y) \, \psi_j(y) \, \mathrm d y \, \mathrm d x \, .\]

pparbem2d par

Some helpers to make parallel computations possible.

For the sake of parallelism we need some auxiliary data structure, which are all collected within the struct par. Most commonly used is this struct for enumerated non linear data structures such als cluster- or blocktrees.

Quadrature rules used within BEM computation.

Singular quadrature formulas used within boundary-integral-Operator computation.

See also
_singquad1d.
void(* transfer_col) (pcbem2d bem, pclusterbasis cb, uint cname)

Computes the transfermatrices $ F_{s'} $ for a cluster $ s \in \mathcal T_{\mathcal J} \, , s' \in \operatorname{sons}(s) $ .

transfer_col will build up the transfer matrices $ F_{s'} $ , for all $ s \in \mathcal T_{\mathcal J}$ and $ s' \in \operatorname{sons}(s) $ . For a cluster $ s \in \mathcal T_{\mathcal J} $ transfer_col will compute the matrices $ F_{s_1}, \ldots, F_{t_\sigma} $ , $ \# \operatorname{sons}(s) = \sigma $ . The matrix $ F_{s_1} $ is stored in cb->son[0]->E, $ F_{s_2} $ is stored in cb->son[1]->E, $ \ldots F_{t_\sigma} $ is stored in cb->son[ $\sigma-1$]->E .

Parameters
bemBEM-object containing additional information for computation of the matrix entries.
cbCurrent column clusterbasis. The Matrices $ F_{s'}, s' \in \operatorname{sons}(s) $ will be saved into cb->son[s'-1]->E .
cnameEnumerated number of current column clusterbasis cb.
void(* transfer_row) (pcbem2d bem, pclusterbasis rb, uint rname)

Computes the transfermatrices $ E_{t'} $ for a cluster $ t \in \mathcal T_{\mathcal I} \, , t' \in \operatorname{sons}(t) $ .

transfer_row will build up the transfer matrices $ E_{t'} $ , for all $ t \in \mathcal T_{\mathcal I}$ and $ t' \in \operatorname{sons}(t) $ . For a cluster $ t \in \mathcal T_{\mathcal I} $ transfer_row will compute the matrices $ E_{t_1}, \ldots, E_{t_\tau} $ , $ \# \operatorname{sons}(t) = \tau $ . The matrix $ E_{t_1} $ is stored in rb->son[0]->E, $ E_{t_2} $ is stored in rb->son[1]->E, $ \ldots E_{t_\tau} $ is stored in rb->son[ $\tau-1$]->E .

Parameters
bemBEM-object containing additional information for computation of the matrix entries.
rbCurrent row clusterbasis. The Matrices $ E_{t'}, t' \in \operatorname{sons}(t) $ will be saved into rb->son[t'-1]->E .
rnameEnumerated number of current row clusterbasis rb.
uint* v2t

vertex to triangle map.

A list of triangle indices corresponding to each vertex index. A vertex $ v $ is part of all triangles ranging from $ \texttt{v2t[v2ts[v]]} \ldots \texttt{v2t[v2ts[v+1]-1]}$

uint* v2ts

Vertex to triangle map stating indices.

Array of starting indices for each vertex needed by v2t. v2ts has length vertices + 1.


The documentation for this struct was generated from the following file: