H2Lib  3.0
singquad2d.h
Go to the documentation of this file.
1 /* ------------------------------------------------------------
2  This is the file "singquad2d.h" of the H2Lib package.
3  All rights reserved, Sven Christophersen 2011-2013
4  ------------------------------------------------------------ */
5 
10 #ifndef SINGQUAD2D_H_
11 #define SINGQUAD2D_H_
12 
13 /* C STD LIBRARY */
14 #include <assert.h>
15 /* CORE 0 */
16 #ifdef USE_SIMD
17 #include "simd.h"
18 #endif
19 #include "settings.h"
20 /* CORE 1 */
21 #include "gaussquad.h"
22 /* CORE 2 */
23 /* CORE 3 */
24 /* SIMPLE */
25 /* PARTICLES */
26 /* BEM */
27 #include "surface3d.h"
28 
38 struct _singquad2d {
119 
120 #ifdef USE_TRIQUADPOINTS
121 
124  real *tri_x;
125 
129  real *tri_y;
130 
134  real *tri_z;
135 #endif
136 
144 };
145 
151 typedef struct _singquad2d singquad2d;
152 
157 
161 typedef const singquad2d *pcsingquad2d;
162 
163 /* ------------------------------------------------------------
164  Constructors and destructors
165  ------------------------------------------------------------ */
166 
181 
188 HEADER_PREFIX void
190 
191 /* ------------------------------------------------------------
192  Weighting quadrature rules
193  ------------------------------------------------------------ */
194 
210 HEADER_PREFIX void
212 
228 HEADER_PREFIX void
230 
246 HEADER_PREFIX void
248 
264 HEADER_PREFIX void
266 
267 /* ------------------------------------------------------------
268  Select quadrature rule
269  ------------------------------------------------------------ */
270 
280 fast_select_quadrature(uint (*geo_t)[3], uint t, uint s);
281 
282 
302 select_quadrature_singquad2d(pcsingquad2d sq, const uint *tv, const uint *sv,
303  uint *tp, uint *sp, real **x, real **y, real **w, uint *n, real *base);
304 
307 #endif /* SINGQUAD2D_H_ */
real * y_edge
Y-component of quadrature points for singular integrals on domains sharing a common edge...
Definition: singquad2d.h:58
real * x_id
X-component of quadrature points for singular integrals on same domain.
Definition: singquad2d.h:40
real base_vert
Constant offset for singular integrals on domains sharing a common vertex.
Definition: singquad2d.h:93
uint nmax
maximal number of quadrature points.
Definition: singquad2d.h:143
uint select_quadrature_singquad2d(pcsingquad2d sq, const uint *tv, const uint *sv, uint *tp, uint *sp, real **x, real **y, real **w, uint *n, real *base)
This function is designed to select the correct quadrature rule for a current pair of triangles...
uint fast_select_quadrature(uint(*geo_t)[3], uint t, uint s)
Determine the number of common vertices of a pair of triangles.
real * y_vert
Y-component of quadrature points for singular integrals on domains sharing a common vertex...
Definition: singquad2d.h:83
void weight_basisfunc_cl_singquad2d(real *x, real *y, real *w, uint nq)
Weighting a quadrature rule for a double integral with a combination of piecewise constant and linear...
real base_single
Constant offset for a single triangle.
Definition: singquad2d.h:116
real * w_dist
Quadrature weights for singular integrals on distant domains.
Definition: singquad2d.h:104
real * y_id
Y-component of quadrature points for singular integrals on same domain.
Definition: singquad2d.h:42
unsigned uint
Unsigned integer type.
Definition: settings.h:70
uint q
Order of basic quadrature rule for single and regular double integrals.
Definition: singquad2d.h:139
singquad2d * psingquad2d
Definition: singquad2d.h:156
const singquad2d * pcsingquad2d
Definition: singquad2d.h:161
real base_dist
Constant offset for singular integrals on distant domains.
Definition: singquad2d.h:106
uint n_edge
Number of quadrature points for singular integrals on domains sharing a common edge.
Definition: singquad2d.h:73
void weight_basisfunc_l_singquad2d(real *x, real *y, real *w, uint nq)
Weighting a quadrature rule for a single integral with linear basis functions.
real base_edge
Constant offset for singular integrals on domains sharing a common edge.
Definition: singquad2d.h:68
real * w_id
Quadrature weights for singular integrals on same domain.
Definition: singquad2d.h:44
real * x_single
X-component of quadrature points for a single triangle.
Definition: singquad2d.h:110
This struct collects all type of quadrature formulas needed by the computation of matrix entries with...
Definition: singquad2d.h:38
uint n_dist
Number of quadrature points for singular integrals on distant domains.
Definition: singquad2d.h:108
void del_singquad2d(psingquad2d sq)
Destructor for singquad2d objects.
Representation of a triangle surface mesh.
Definition: surface3d.h:45
real * y_dist
Y-component of quadrature points for singular integrals on distant domains.
Definition: singquad2d.h:102
#define HEADER_PREFIX
Prefix for function declarations.
Definition: settings.h:43
real * w_single
Quadrature weights for a single triangle.
Definition: singquad2d.h:114
double real
real floating point type.
Definition: settings.h:97
uint n_vert
Number of quadrature points for singular integrals on domains sharing a common vertex.
Definition: singquad2d.h:98
real * x_edge
X-component of quadrature points for singular integrals on domains sharing a common edge...
Definition: singquad2d.h:53
uint q2
Order of basic quadrature rule for singular double integrals.
Definition: singquad2d.h:141
uint n_single
Number of quadrature points for a single triangle.
Definition: singquad2d.h:118
uint n_id
Number of quadrature points for singular integrals on same domain.
Definition: singquad2d.h:48
real * w_vert
Quadrature weights for singular integrals on domains sharing a common vertex.
Definition: singquad2d.h:88
void weight_basisfunc_lc_singquad2d(real *x, real *y, real *w, uint nq)
Weighting a quadrature rule for a double integral with a combination of piecewise linear and constant...
real * w_edge
Quadrature weights for singular integrals on domains sharing a common edge.
Definition: singquad2d.h:63
real * y_single
Y-component of quadrature points for a single triangle.
Definition: singquad2d.h:112
real * x_vert
X-component of quadrature points for singular integrals on domains sharing a common vertex...
Definition: singquad2d.h:78
real base_id
Constant offset for singular integrals on same domain.
Definition: singquad2d.h:46
void weight_basisfunc_ll_singquad2d(real *x, real *y, real *w, uint nq)
Weighting a quadrature rule for a double integral with both linear basis functions.
psingquad2d build_singquad2d(pcsurface3d gr, uint q, uint q2)
Creates a new singquad2d object containing all necessary quadrature rules for singular integrals aris...
real * x_dist
X-component of quadrature points for singular integrals on distant domains.
Definition: singquad2d.h:100