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

greencross.h
Go to the documentation of this file.
1 /* ------------------------------------------------------------
2  * This is the file "greencross.h" of this master thesis.
3  * All rights reserved, Bennet Carstensen 2017
4  * ------------------------------------------------------------ */
5 
13 #ifndef GREENCROSS_H
14 #define GREENCROSS_H
15 
16 #include "curve2d.h"
17 #include "h2matrix.h"
18 #include "surface3d.h"
19 
26 const uint greencross_min_dim = 2;
27 const uint greencross_max_dim = 3;
28 
31 typedef struct _greencross greencross;
32 
35 
37 typedef const greencross *pcgreencross;
38 
42 {
44  void *geom;
45 
47  uint dim;
48 
50  uint n;
51 
53  uint *idx;
54 
56  uint nq;
57 
59  preal zq;
60 
62  preal wq;
63 
65  void *sq;
66 
68  uint ng;
69 
71  preal zg;
72 
74  preal wg;
75 
77  uint K;
78 
79  ph2matrix h2;
80 
82  real (*kernel_function) (const field* x, const field* y, const uint dim);
83 
86  real (*pdx_kernel_function) (const field* x,
87  const field* y,
88  const uint dim,
89  const uint i);
90 
93  real (*pdy_kernel_function) (const field* x,
94  const field* y,
95  const uint dim,
96  const uint i);
97 };
98 
109 HEADER_PREFIX void
111 
120 HEADER_PREFIX void
122 
146 HEADER_PREFIX pgreencross
147 new_laplace2d_greencross(pcurve2d c2d, uint res, void *eta, uint q, uint m);
148 
153 HEADER_PREFIX void
155 
179 HEADER_PREFIX void
181  const uint rsize,
182  const uint *ridx,
183  const uint csize,
184  const uint *cidx,
185  pamatrix G);
186 
212 HEADER_PREFIX void
213 fill_green_left_greencross(pcgreencross gc, pccluster t, pamatrix A);
214 
241 HEADER_PREFIX void
243  pccluster t,
244  pccluster s,
245  pamatrix B);
246 
260 HEADER_PREFIX void
261 green_cross_approximation(pcgreencross gc, ph2matrix H2);
262 
263 
266 #endif // GREENCROSS_H
greencross * pgreencross
Pointer to a greencross object.
Definition: greencross.h:34
uint ng
Approximation order / quadrature order for Green's formula.
Definition: greencross.h:68
real(* pdy_kernel_function)(const field *x, const field *y, const uint dim, const uint i)
Partial derivative of the kernel function describing the problem in the i-th component of y...
Definition: greencross.h:93
Main container object for performing the Green cross approximation method.
Definition: greencross.h:41
const uint greencross_min_dim
Definition: greencross.h:26
preal wg
Quadrature weights for Green's formula.
Definition: greencross.h:74
ph2matrix h2
Definition: greencross.h:79
uint dim
Dimension of the problem.
Definition: greencross.h:47
preal zg
Quadrature points in for Green's formula.
Definition: greencross.h:71
uint nq
Quadrature order for regular integrals.
Definition: greencross.h:56
HEADER_PREFIX void init_greencross(pgreencross gc, uint dim)
Initilize components related to the dimension of the problem and zero initializes all other component...
Definition: greencross.c:243
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 describ...
Definition: greencross.c:415
real(* pdx_kernel_function)(const field *x, const field *y, const uint dim, const uint i)
Partial derivative of the kernel function describing the problem in the i-th component of x...
Definition: greencross.h:86
const greencross * pcgreencross
Pointer to a constant greencross object.
Definition: greencross.h:37
HEADER_PREFIX void del_greencross(pgreencross gc)
Free memory and set pointer to NULL of the corresponding _greencross greencross object gc...
Definition: greencross.c:403
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 -Matrix ...
Definition: greencross.c:781
uint * idx
Index set.
Definition: greencross.h:53
preal zq
Quadrature points in for regular integrals.
Definition: greencross.h:59
uint K
Local rank of approximated submatrices.
Definition: greencross.h:77
void * geom
Geometry of the problem.
Definition: greencross.h:44
real(* kernel_function)(const field *x, const field *y, const uint dim)
Kernel function describing the problem.
Definition: greencross.h:82
uint n
Number of basis functions.
Definition: greencross.h:50
HEADER_PREFIX void fill_green_right_greencross(pcgreencross gc, pccluster t, pccluster s, pamatrix B)
Definition: greencross.c:650
void * sq
Quadrature for singular integrals based on q, zq and @ wq.
Definition: greencross.h:65
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 , where .
Definition: greencross.c:348
HEADER_PREFIX void fill_green_left_greencross(pcgreencross gc, pccluster t, pamatrix A)
Definition: greencross.c:522
preal wq
Quadrature weights for regular integrals.
Definition: greencross.h:62
const uint greencross_max_dim
Definition: greencross.h:27
HEADER_PREFIX void uninit_greencross(pgreencross gc)
Uninitilize a greencross object.
Definition: greencross.c:266