Scalable solving of boundary element methods utilizing the Green cross approximation method and GPUs
This project is maintained by Bennet Carstensen
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 greencross * | pgreencross |
Pointer to a greencross object. More... | |
typedef const greencross * | pcgreencross |
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 |
Algorithms and functions to perform the Green cross approximation method.
typedef struct _greencross greencross |
greencross is just an abbreviation for the struct _greencross.
typedef const greencross* pcgreencross |
Pointer to a constant greencross object.
typedef greencross* pgreencross |
Pointer to a greencross object.
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
.
dim
= 2, data
is expected to be an object of type curve2d.[in] | data | Pointer to an object containing the informations to construct the cluster geometry. |
[in] | dim | Dimension of the cluster geometry. |
[out] | idx | Index 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. |
idx
. HEADER_PREFIX void del_greencross | ( | pgreencross | gc | ) |
Free memory and set pointer to NULL of the corresponding _greencross greencross object gc
.
gc | The greencross object which needs to be deleted. |
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
.
[in] | gc | The greencross contains the problem description and several data to approximate the integral. |
[in] | t | Cluster containing the indices of the base functions needed to discretize the integral and mapping informations for the green quadrature points. |
[out] | A | Matrix object which will be filled with \(A_{t}\). |
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
.
[in] | gc | The greencross object contains the problem description and several data to approximate the integral. |
[in] | t | Cluster containing mapping informations for the green quadrature points. |
[in] | s | Cluster containing the indices of the base functions needed to discretize the integral. |
[out] | B | Matrix object which will be filled with \(A_{t}\). |
uint get_max_lvl_clusterbasis | ( | pclusterbasis | cb, |
const uint | level | ||
) |
Retrieve the maximum level of a cluster basis.
level
should be 0 since it's only used insid the function to keep track of the level of the current viewed cluster basis.cb | The cluster basis of which we search the maximum level. |
level | The level of the current viewed cluster basis. |
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].
[in] | gc | The 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. |
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
.
[in,out] | gc | The greencross object which needs to be initialized. |
[in] | dim | The dimension of the problem. |
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
.
[in] | gc | Object containing the problem describtion like base functions and Quadratur choices. |
[in] | rsize | Number of base functions used from greencross::geom for the row of the Galerkin discretization matrix. |
[in] | ridx | Indices of the base functions used from greencross::geom for the row of the Galerkin discretization matrix. |
[in] | csize | Number of base functions used from greencross::geom for the column of the Galerkin discretization matrix. |
[in] | cidx | Indices of the base functions used from greencross::geom for the column of the Galerkin discretization matrix. |
[out] | G | Matrix which will contain the Galerkin discretization. |
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
.
[in] | c2d | Geometry describing the subspace \(\Omega \subset \mathbb{R}^2 \). |
[in] | res | Maximal size of leaf clusters. |
[in] | eta | Necessary data for the admissibility condition in the hierarchical clustering. |
[in] | q | Quadratur order. |
[in] | m | Approximation order. |
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.
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.cb | The root cluster basis. |
level | The level of the current viewed cluster basis. |
on_level | An array whose entries will hold the number of cluster bases for each level. |
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.
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).cb | The root cluster basis. |
level | The level of the current viewed cluster basis. |
lvl_idx | The number of cluster bases already found on each level respectively. |
cb_per_level | An arrays whose entries will hold pointers to all cluster bases of the corresponding level. |
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.
[out] | gc | The greencross object which needs to be initialized. |
const uint greencross_max_dim = 3 |
const uint greencross_min_dim = 2 |