My master thesis

Scalable solving of boundary element methods utilizing the Green cross approximation method and GPUs

This project is maintained by Bennet Carstensen

Documentation

Algorithms and functions to perform the Green cross approximation method. More...

Data Structures

struct  _greencross
 Main container object for performing the Green cross approximation method. More...
 

Typedefs

typedef struct _greencross greencross
 greencross is just an abbreviation for the struct _greencross. More...
 
typedef greencrosspgreencross
 Pointer to a greencross object. More...
 
typedef const greencrosspcgreencross
 Pointer to a constant greencross object. More...
 

Functions

INLINE_PREFIX uint get_max_lvl_clusterbasis (pclusterbasis cb, const uint level)
 Retrieve the maximum level of a cluster basis. More...
 
INLINE_PREFIX void set_on_lvl_info_clusterbasis (pclusterbasis cb, const uint level, uint *on_level)
 Retrieve the number of cluster basis for each level from a given root cluster basis. More...
 
INLINE_PREFIX void set_per_level_clusterbasis (pclusterbasis cb, const uint level, uint *lvl_idx, pclusterbasis **cb_per_level)
 Retrieve the cluster bases of a given root cluster basis divided into the level they are on. More...
 
INLINE_PREFIX pclustergeometry build_clustergeometry_greencross (const void *data, const uint dim, uint **idx)
 Creats a clustergeometry object from a triangular mesh. More...
 
HEADER_PREFIX void init_greencross (pgreencross gc, uint dim)
 Initilize components related to the dimension of the problem and zero initializes all other components of the given greencross object gc. More...
 
HEADER_PREFIX void uninit_greencross (pgreencross gc)
 Uninitilize a greencross object. More...
 
HEADER_PREFIX pgreencross new_laplace2d_greencross (pcurve2d c2d, uint res, void *eta, uint q, uint m)
 Prepares a greencross object for approximating the integral equation \( \int\limits_{\Omega} -\frac{1}{2 \pi} \log{\left \Vert x - y \right \Vert_2} u(y) dy = f(x) \), where \(\Omega \subset \mathbb{R}^2 \). More...
 
HEADER_PREFIX void del_greencross (pgreencross gc)
 Free memory and set pointer to NULL of the corresponding _greencross greencross object gc. More...
 
HEADER_PREFIX void nearfield_greencross (pcgreencross gc, const uint rsize, const uint *ridx, const uint csize, const uint *cidx, pamatrix G)
 Constructs the matrix resulting from the Galerkin discretization of a variational formulation described in gc. More...
 
HEADER_PREFIX void fill_green_left_greencross (pcgreencross gc, pccluster t, pamatrix A)
 
HEADER_PREFIX void fill_green_right_greencross (pcgreencross gc, pccluster t, pccluster s, pamatrix B)
 
HEADER_PREFIX void green_cross_approximation (pcgreencross gc, ph2matrix H2)
 Constructs a Green cross approximation of a problem given by a greencross, represented as an \(\mathcal{H}^{2}\)-Matrix [1]. More...
 

Variables

const uint greencross_min_dim = 2
 
const uint greencross_max_dim = 3
 

Detailed Description

Algorithms and functions to perform the Green cross approximation method.

Typedef Documentation

◆ greencross

typedef struct _greencross greencross

greencross is just an abbreviation for the struct _greencross.

◆ pcgreencross

typedef const greencross* pcgreencross

Pointer to a constant greencross object.

◆ pgreencross

Pointer to a greencross object.

Function Documentation

◆ build_clustergeometry_greencross()

pclustergeometry build_clustergeometry_greencross ( const void *  data,
const uint  dim,
uint **  idx 
)

Creats a clustergeometry object from a triangular mesh.

The geometry is taken from data interpreted as different objects depending on the dimension of the problem dim.

Attention
For dim = 2, data is expected to be an object of type curve2d.
Parameters
[in]dataPointer to an object containing the informations to construct the cluster geometry.
[in]dimDimension of the cluster geometry.
[out]idxIndex set which will be allocated by this method to index each triangle in data. The entries of idx will be 0, 1, ..., N-1, where N is equal to the number of triangles.
Returns
A valid cluster geometry object which can be used to construct a cluster tree along with the index set idx.

◆ del_greencross()

HEADER_PREFIX void del_greencross ( pgreencross  gc)

Free memory and set pointer to NULL of the corresponding _greencross greencross object gc.

Parameters
gcThe greencross object which needs to be deleted.

◆ fill_green_left_greencross()

HEADER_PREFIX void fill_green_left_greencross ( pcgreencross  gc,
pccluster  t,
pamatrix  A 
)

Constructs the first factor resulting of using Green's representation formula as a degenerated approximation of the Galerkin discretization [1].

Constructs the matrix \(A_{t} = (A_{t+} \ A_{t-}) \in \mathbb{R^{\hat{t} \times K}}\), where \(\hat{t}\) and K are given by t and gc respectively, and

\[ \begin{align*} a_{t+, i \nu} &:= \sqrt{w_{\nu}} \int\limits_{\Omega} \Phi_{i}(x) g(x, z_{\nu}) dx, \\ a_{t-, i \nu} &:= \delta_{t} \sqrt{w_{\nu}} \int\limits_{\Omega} \varphi_{i}(x) \frac{\partial g}{\partial n_{\iota}}(x, z_{\nu}) dx, \end{align*} \]

where \(w_{\nu}, \delta{t}, \Omega, \varphi_i, n_{\iota}, z_{\nu}\) and the quadrature rules can be obtained or calculated from gc and t.

Parameters
[in]gcThe greencross contains the problem description and several data to approximate the integral.
[in]tCluster containing the indices of the base functions needed to discretize the integral and mapping informations for the green quadrature points.
[out]AMatrix object which will be filled with \(A_{t}\).

◆ fill_green_right_greencross()

HEADER_PREFIX void fill_green_right_greencross ( pcgreencross  gc,
pccluster  t,
pccluster  s,
pamatrix  B 
)

Constructs the second factor resulting of using Green's representation formula as a degenerated approximation of the Galerkin discretization [1].

Constructs the matrix \(B_{ts} = (B_{ts+} \ B_{ts-}) \in \mathbb{R^{\hat{s} \times K}}\), where \(\hat{s}\) and K are given by s and gc respectively, and

\[ \begin{align*} b_{t+, j \nu} &:= \sqrt{w_{\nu}} \int\limits_{\Omega} \varphi_{j}(y) \frac{\partial g}{\partial n_{\iota}}(z_{\nu}, y) dy, \\ b_{t-, j \nu} &:= -\frac{\sqrt{w_{\nu}}}{\delta_{t}}\int\limits_{\Omega} \Phi_{j}(y) g(z_{\nu}, y) dy, \\ \end{align*} \]

where \(w_{\nu}, \delta{t}, \Omega, \varphi_i, n_{\iota}, z_{\nu}\) and the quadrature rules can obtained or calculated from gc, t and s.

Parameters
[in]gcThe greencross object contains the problem description and several data to approximate the integral.
[in]tCluster containing mapping informations for the green quadrature points.
[in]sCluster containing the indices of the base functions needed to discretize the integral.
[out]BMatrix object which will be filled with \(A_{t}\).

◆ get_max_lvl_clusterbasis()

uint get_max_lvl_clusterbasis ( pclusterbasis  cb,
const uint  level 
)

Retrieve the maximum level of a cluster basis.

Attention
level should be 0 since it's only used insid the function to keep track of the level of the current viewed cluster basis.
Parameters
cbThe cluster basis of which we search the maximum level.
levelThe level of the current viewed cluster basis.

◆ green_cross_approximation()

HEADER_PREFIX void green_cross_approximation ( pcgreencross  gc,
ph2matrix  H2 
)

Constructs a Green cross approximation of a problem given by a greencross, represented as an \(\mathcal{H}^{2}\)-Matrix [1].

Parameters
[in]gcThe greencross object contains the problem description and several data to perform the approximation.
[out]H2\(\mathcal{H}^{2}\)-Matrix object which will be filled according to the Green cross approximation method.

◆ init_greencross()

HEADER_PREFIX void init_greencross ( pgreencross  gc,
uint  dim 
)

Initilize components related to the dimension of the problem and zero initializes all other components of the given greencross object gc.

Note
Should be called at the beginning of the initialization of each greencross object.
Should always be matched with a call of uninit_greencross.
Parameters
[in,out]gcThe greencross object which needs to be initialized.
[in]dimThe dimension of the problem.

◆ nearfield_greencross()

HEADER_PREFIX void nearfield_greencross ( pcgreencross  gc,
const uint  rsize,
const uint *  ridx,
const uint  csize,
const uint *  cidx,
pamatrix  G 
)

Constructs the matrix resulting from the Galerkin discretization of a variational formulation described in gc.

The arrays ridx and cidx serve as the indices for the set of base functions described through the geometry of the problem greencross::geom. Further choices like quadrature for solving the integral are taken from gc.

Parameters
[in]gcObject containing the problem describtion like base functions and Quadratur choices.
[in]rsizeNumber of base functions used from greencross::geom for the row of the Galerkin discretization matrix.
[in]ridxIndices of the base functions used from greencross::geom for the row of the Galerkin discretization matrix.
[in]csizeNumber of base functions used from greencross::geom for the column of the Galerkin discretization matrix.
[in]cidxIndices of the base functions used from greencross::geom for the column of the Galerkin discretization matrix.
[out]GMatrix which will contain the Galerkin discretization.

◆ new_laplace2d_greencross()

HEADER_PREFIX pgreencross new_laplace2d_greencross ( pcurve2d  c2d,
uint  res,
void *  eta,
uint  q,
uint  m 
)

Prepares a greencross object for approximating the integral equation \( \int\limits_{\Omega} -\frac{1}{2 \pi} \log{\left \Vert x - y \right \Vert_2} u(y) dy = f(x) \), where \(\Omega \subset \mathbb{R}^2 \).

Given a subspace \(\Omega \subset \mathbb{R}^2 \) by c2d, an object to approximate the integral equation \( \int\limits_{\Omega} -\frac{1}{2 \pi} \log{\left \Vert x - y \right \Vert_2} u(y) dy = f(x) \) is given with an approximation order m, a quadrature order q, necessary data for the admissibility condition in the hierarchical clustering eta and a maximal size of leaf clusters res.

Parameters
[in]c2dGeometry describing the subspace \(\Omega \subset \mathbb{R}^2 \).
[in]resMaximal size of leaf clusters.
[in]etaNecessary data for the admissibility condition in the hierarchical clustering.
[in]qQuadratur order.
[in]mApproximation order.
Returns
A greencross object ready for approximating the integral equation \( \int\limits_{\Omega} -\frac{1}{2 \pi} \log{\left \Vert x - y \right \Vert_2} u(y) dy = f(x) \).

◆ set_on_lvl_info_clusterbasis()

void set_on_lvl_info_clusterbasis ( pclusterbasis  cb,
const uint  level,
uint *  on_level 
)

Retrieve the number of cluster basis for each level from a given root cluster basis.

After execution, each entry of on_level with index i will hold the number of cluster bases whose level are i.

Attention
level should be 0 since it's only used insid the function to keep track of the level of the current viewed cluster basis.
on_level should be an zero initialized array of length get_max_lvl_clusterbasis(cb, level) + 1.
Parameters
cbThe root cluster basis.
levelThe level of the current viewed cluster basis.
on_levelAn array whose entries will hold the number of cluster bases for each level.

◆ set_per_level_clusterbasis()

void set_per_level_clusterbasis ( pclusterbasis  cb,
const uint  level,
uint *  lvl_idx,
pclusterbasis **  cb_per_level 
)

Retrieve the cluster bases of a given root cluster basis divided into the level they are on.

After execution, each entry of cb_per_level will hold the pointers to all cluster bases whose level is i, where i is the index of the entry inside cb_per_level.

Attention
level should be 0 since it's only used insid the function to keep track of the level of the current viewed cluster basis.
lvl_idx should be an zero initialized array of length get_max_lvl_clusterbasis(cb, level) + 1.
cb_per_level should be an initialized array of length get_max_lvl_clusterbasis(cb, level) + 1, where each entry with index should be an initialized array of length on_level[i], where on_level is retrieved by set_on_lvl_info_clusterbasis(cb, 0, on_level).
Parameters
cbThe root cluster basis.
levelThe level of the current viewed cluster basis.
lvl_idxThe number of cluster bases already found on each level respectively.
cb_per_levelAn arrays whose entries will hold pointers to all cluster bases of the corresponding level.

◆ uninit_greencross()

HEADER_PREFIX void uninit_greencross ( pgreencross  gc)

Uninitilize a greencross object.

Invalidates member pointers, freeing corresponding storage if appropriate, and prepares the object for deletion.

Note
Should always be matched with a call of init_greencross.
Parameters
[out]gcThe greencross object which needs to be initialized.

Variable Documentation

◆ greencross_max_dim

const uint greencross_max_dim = 3

◆ greencross_min_dim

const uint greencross_min_dim = 2