H2Lib  3.0
bem3d.h
Go to the documentation of this file.
1 /* ------------------------------------------------------------
2  * This is the file "bem3d.h" of the H2Lib package.
3  * All rights reserved, Sven Christophersen 2011
4  * ------------------------------------------------------------ */
5 
12 #ifndef BEM3D_H_
13 #define BEM3D_H_
14 
15 /* C STD LIBRARY */
16 #include <stdio.h>
17 #ifdef USE_OPENMP
18 #include <omp.h>
19 #endif
20 /* CORE 0 */
21 #ifdef USE_SIMD
22 #include "simd.h"
23 #endif
24 #include "basic.h"
25 /* CORE 1 */
26 #include "amatrix.h"
27 #include "gaussquad.h"
28 #include "factorizations.h"
29 #include "krylov.h"
30 #include "sparsematrix.h"
31 /* CORE 2 */
32 #include "cluster.h"
33 #include "clustergeometry.h"
34 #include "block.h"
35 #include "hmatrix.h"
36 #include "clusterbasis.h"
37 #include "h2matrix.h"
38 #include "dcluster.h"
39 #include "dblock.h"
40 #include "dclusterbasis.h"
41 #include "dh2matrix.h"
42 #include "krylovsolvers.h"
43 /* CORE 3 */
44 #include "hcoarsen.h"
45 #include "h2update.h"
46 #include "aca.h"
47 /* SIMPLE */
48 /* BEM */
49 #include "macrosurface3d.h"
50 #include "surface3d.h"
51 #include "singquad2d.h"
52 
67 /* ------------------------------------------------------------
68  * Type definitions
69  * ------------------------------------------------------------ */
70 
75 typedef struct _bem3d bem3d;
79 typedef bem3d *pbem3d;
83 typedef const bem3d* pcbem3d;
84 
92 typedef struct _aprxbem3d aprxbem3d;
100 typedef const aprxbem3d *pcaprxbem3d;
101 
107 typedef struct _parbem3d parbem3d;
115 typedef const parbem3d *pcparbem3d;
116 
122 typedef struct _kernelbem3d kernelbem3d;
123 
128 
132 typedef const kernelbem3d *pckernelbem3d;
133 
153 typedef void (*quadpoints3d)(pcbem3d bem, const real bmin[3],
154  const real bmax[3], const real delta, real (**Z)[3], real (**N)[3]);
155 
176 };
177 
182 
197 typedef field (*boundary_func3d)(const real *x, const real *n, void *data);
198 
210 typedef field (*kernel_func3d)(const real *x, const real *y, const real *nx,
211  const real *ny, void *data);
212 
230 typedef field (*kernel_wave_func3d)(const real *x, const real *y,
231  const real *nx, const real *ny, pcreal dir, void *data);
232 
233 #ifdef USE_SIMD
234 typedef void (*kernel_simd_func3d)(const vreal *x, const vreal *y,
235  const vreal *nx, const vreal *ny, void *data, vreal *res_re, vreal *res_im);
236 
237 typedef void (*kernel_simd_wave_func3d)(const vreal *x, const vreal *y,
238  const vreal *nx, const vreal *ny, pcreal dir, void *data, vreal *res_re,
239  vreal *res_im);
240 #endif
241 
245 typedef struct _listnode listnode;
246 
251 
255 typedef struct _tri_list tri_list;
256 
261 
265 typedef struct _vert_list vert_list;
266 
271 
290 struct _bem3d {
298 
308 
313 
318 
328 
333 
334 
339 
345 
357 
380  void (*nearfield)(const uint *ridx, const uint *cidx, pcbem3d bem,
381  bool ntrans, pamatrix N);
382 
405  void (*nearfield_far)(const uint *ridx, const uint *cidx, pcbem3d bem,
406  bool ntrans, pamatrix N);
407 
428  void (*farfield_rk)(pccluster rc, uint rname, pccluster cc, uint cname,
429  pcbem3d bem, prkmatrix R);
430 
453  void (*farfield_u)(uint rname, uint cname, uint bname, pcbem3d bem);
454 
477  void (*farfield_wave_u)(uint rname, uint cname, uint bname, pcbem3d bem);
478 
496  void (*leaf_row)(uint rname, pcbem3d bem);
497 
515  void (*leaf_wave_row)(uint rname, pcbem3d bem);
516 
534  void (*leaf_col)(uint cname, pcbem3d bem);
535 
554  void (*leaf_wave_col)(uint cname, pcbem3d bem);
555 
573  void (*transfer_row)(uint rname, pcbem3d bem);
574 
599  void (*transfer_wave_row)(uint rname, pcbem3d bem);
600 
633  void (*transfer_wave_wave_row)(uint rname, pcbem3d bem);
634 
652  void (*transfer_col)(uint cname, pcbem3d bem);
653 
678  void (*transfer_wave_col)(uint cname, pcbem3d bem);
679 
712  void (*transfer_wave_wave_col)(uint cname, pcbem3d bem);
713 
722 
731 
749 };
750 
768 struct _kernelbem3d {
791  void (*fundamental)(pcbem3d bem, const real (*X)[3], const real (*Y)[3],
792  pamatrix V);
793 
821  void (*fundamental_wave)(pcbem3d bem, const real (*X)[3], const real (*Y)[3],
822  pcreal dir, pamatrix V);
823 
850  void (*dny_fundamental)(pcbem3d bem, const real (*X)[3], const real (*Y)[3],
851  const real (*NY)[3], pamatrix V);
852 
881  void (*dnx_dny_fundamental)(pcbem3d bem, const real (*X)[3],
882  const real (*NX)[3], const real (*Y)[3], const real (*NY)[3], pamatrix V);
883 
911  void (*kernel_row)(const uint *idx, const real (*Z)[3], pcbem3d bem,
912  pamatrix A);
913 
942  void (*kernel_col)(const uint *idx, const real (*Z)[3], pcbem3d bem,
943  pamatrix B);
944 
976  void (*dnz_kernel_row)(const uint *idx, const real (*Z)[3],
977  const real (*NZ)[3], pcbem3d bem, pamatrix A);
978 
1010  void (*dnz_kernel_col)(const uint *idx, const real (*Z)[3],
1011  const real (*NZ)[3], pcbem3d bem, pamatrix B);
1012 
1036  void (*fundamental_row)(const uint *idx, const real (*Z)[3], pcbem3d bem,
1037  pamatrix A);
1038 
1063  void (*fundamental_col)(const uint *idx, const real (*Z)[3], pcbem3d bem,
1064  pamatrix B);
1065 
1093  void (*dnz_fundamental_row)(const uint *idx, const real (*Z)[3],
1094  const real (*NZ)[3], pcbem3d bem, pamatrix A);
1095 
1123  void (*dnz_fundamental_col)(const uint *idx, const real (*Z)[3],
1124  const real (*NZ)[3], pcbem3d bem, pamatrix B);
1125 
1160  void (*lagrange_row)(const uint *idx, pcrealavector px, pcrealavector py,
1161  pcrealavector pz, pcbem3d bem, pamatrix V);
1162 
1201  void (*lagrange_wave_row)(const uint *idx, pcrealavector px, pcrealavector py,
1202  pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V);
1203 
1238  void (*lagrange_col)(const uint *idx, pcrealavector px, pcrealavector py,
1239  pcrealavector pz, pcbem3d bem, pamatrix W);
1240 
1279  void (*lagrange_wave_col)(const uint *idx, pcrealavector px, pcrealavector py,
1280  pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix W);
1281 };
1282 
1290 struct _listnode {
1295 };
1296 
1304 struct _tri_list {
1314 };
1315 
1323 struct _vert_list {
1328 };
1329 
1330 /* ------------------------------------------------------------
1331  * Constructors and destructors
1332  * ------------------------------------------------------------ */
1333 
1345 
1351 HEADER_PREFIX void
1353 
1364 new_tri_list(ptri_list next);
1365 
1371 HEADER_PREFIX void
1373 
1392 
1399 HEADER_PREFIX void
1400 del_bem3d(pbem3d bem);
1401 
1402 /* ------------------------------------------------------------
1403  * Methods to build clustertrees
1404  * ------------------------------------------------------------ */
1405 
1425 
1446 
1472 
1495 
1496 /* ------------------------------------------------------------
1497  * Nearfield quadrature routines
1498  * ------------------------------------------------------------ */
1499 
1523 HEADER_PREFIX void
1524 assemble_cc_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1525  bool ntrans, pamatrix N, kernel_func3d kernel);
1526 #ifdef USE_SIMD
1527 HEADER_PREFIX void
1528 assemble_cc_simd_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1529  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1530 #endif
1531 
1558 HEADER_PREFIX void
1559 assemble_cc_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1560  bool ntrans, pamatrix N, kernel_func3d kernel);
1561 #ifdef USE_SIMD
1562 HEADER_PREFIX void
1563 assemble_cc_simd_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1564  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1565 #endif
1566 
1591 HEADER_PREFIX void
1592 assemble_cl_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1593  bool ntrans, pamatrix N, kernel_func3d kernel);
1594 #ifdef USE_SIMD
1595 HEADER_PREFIX void
1596 assemble_cl_simd_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1597  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1598 #endif
1599 
1627 HEADER_PREFIX void
1628 assemble_cl_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1629  bool ntrans, pamatrix N, kernel_func3d kernel);
1630 #ifdef USE_SIMD
1631 HEADER_PREFIX void
1632 assemble_cl_simd_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1633  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1634 #endif
1635 
1660 HEADER_PREFIX void
1661 assemble_lc_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1662  bool ntrans, pamatrix N, kernel_func3d kernel);
1663 #ifdef USE_SIMD
1664 HEADER_PREFIX void
1665 assemble_lc_simd_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1666  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1667 #endif
1668 
1696 HEADER_PREFIX void
1697 assemble_lc_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1698  bool ntrans, pamatrix N, kernel_func3d kernel);
1699 #ifdef USE_SIMD
1700 HEADER_PREFIX void
1701 assemble_lc_simd_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1702  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1703 #endif
1704 
1728 HEADER_PREFIX void
1729 assemble_ll_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1730  bool ntrans, pamatrix N, kernel_func3d kernel);
1731 #ifdef USE_SIMD
1732 HEADER_PREFIX void
1733 assemble_ll_simd_near_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1734  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1735 #endif
1736 
1763 HEADER_PREFIX void
1764 assemble_ll_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1765  bool ntrans, pamatrix N, kernel_func3d kernel);
1766 #ifdef USE_SIMD
1767 HEADER_PREFIX void
1768 assemble_ll_simd_far_bem3d(const uint * ridx, const uint * cidx, pcbem3d bem,
1769  bool ntrans, pamatrix N, kernel_simd_func3d kernel);
1770 #endif
1771 
1772 /* ------------------------------------------------------------
1773  * Evaluate kernel function at given points
1774  * ------------------------------------------------------------ */
1775 
1798 HEADER_PREFIX void
1799 fill_bem3d(pcbem3d bem, const real (*X)[3], const real (*Y)[3],
1800  const real (*NX)[3], const real (*NY)[3], pamatrix V, kernel_func3d kernel);
1801 
1829 HEADER_PREFIX void
1830 fill_wave_bem3d(pcbem3d bem, const real (*X)[3], const real (*Y)[3],
1831  const real (*NX)[3], const real (*NY)[3], pamatrix V, pcreal dir,
1832  kernel_wave_func3d kernel);
1833 
1834 /* ------------------------------------------------------------
1835  * Compute single integrals with kernel functions
1836  * ------------------------------------------------------------ */
1837 
1864 HEADER_PREFIX void
1865 fill_row_c_bem3d(const uint * idx, const real (*Z)[3], pcbem3d bem, pamatrix V,
1866  kernel_func3d kernel);
1867 #ifdef USE_SIMD
1868 HEADER_PREFIX void
1869 fill_row_simd_c_bem3d(const uint * idx, const real (*Z)[3], pcbem3d bem,
1870  pamatrix V, kernel_simd_func3d kernel);
1871 #endif
1872 
1899 HEADER_PREFIX void
1900 fill_col_c_bem3d(const uint * idx, const real (*Z)[3], pcbem3d bem, pamatrix V,
1901  kernel_func3d kernel);
1902 #ifdef USE_SIMD
1903 HEADER_PREFIX void
1904 fill_col_simd_c_bem3d(const uint * idx, const real (*Z)[3], pcbem3d bem,
1905  pamatrix V, kernel_simd_func3d kernel);
1906 #endif
1907 
1934 HEADER_PREFIX void
1935 fill_row_l_bem3d(const uint * idx, const real (*Z)[3], pcbem3d bem, pamatrix V,
1936  kernel_func3d kernel);
1937 
1964 HEADER_PREFIX void
1965 fill_col_l_bem3d(const uint * idx, const real (*Z)[3], pcbem3d bem, pamatrix V,
1966  kernel_func3d kernel);
1967 
1996 HEADER_PREFIX void
1997 fill_dnz_row_c_bem3d(const uint * idx, const real (*Z)[3], const real (*N)[3],
1998  pcbem3d bem, pamatrix V, kernel_func3d kernel);
1999 #ifdef USE_SIMD
2000 HEADER_PREFIX void
2001 fill_dnz_row_simd_c_bem3d(const uint * idx, const real (*Z)[3],
2002  const real (*N)[3], pcbem3d bem, pamatrix V, kernel_simd_func3d kernel);
2003 #endif
2004 
2033 HEADER_PREFIX void
2034 fill_dnz_col_c_bem3d(const uint * idx, const real (*Z)[3], const real (*N)[3],
2035  pcbem3d bem, pamatrix V, kernel_func3d kernel);
2036 #ifdef USE_SIMD
2037 HEADER_PREFIX void
2038 fill_dnz_col_simd_c_bem3d(const uint * idx, const real (*Z)[3],
2039  const real (*N)[3], pcbem3d bem, pamatrix V, kernel_simd_func3d kernel);
2040 #endif
2041 
2070 HEADER_PREFIX void
2071 fill_dnz_row_l_bem3d(const uint * idx, const real (*Z)[3], const real (*N)[3],
2072  pcbem3d bem, pamatrix V, kernel_func3d kernel);
2073 
2102 HEADER_PREFIX void
2103 fill_dnz_col_l_bem3d(const uint * idx, const real (*Z)[3], const real (*N)[3],
2104  pcbem3d bem, pamatrix V, kernel_func3d kernel);
2105 
2106 /* ------------------------------------------------------------
2107  * Initializer functions for h-matrix approximations
2108  * ------------------------------------------------------------ */
2109 
2140 HEADER_PREFIX void
2141 setup_hmatrix_recomp_bem3d(pbem3d bem, bool recomp, real accur_recomp,
2142  bool coarsen, real accur_coarsen);
2143 
2144 /* ------------------------------------------------------------
2145  * Interpolation
2146  * ------------------------------------------------------------ */
2147 
2177 HEADER_PREFIX void
2179  pcblock tree, uint m);
2180 
2210 HEADER_PREFIX void
2212  pcblock tree, uint m);
2213 
2237 HEADER_PREFIX void
2239  pcblock tree, uint m);
2240 
2241 /* ------------------------------------------------------------
2242  * Green
2243  * ------------------------------------------------------------ */
2244 
2305 HEADER_PREFIX void
2307  pcblock tree, uint m, uint l, real delta, quadpoints3d quadpoints);
2308 
2371 HEADER_PREFIX void
2373  pcblock tree, uint m, uint l, real delta, quadpoints3d quadpoints);
2374 
2411 HEADER_PREFIX void
2413  pcblock tree, uint m, uint l, real delta, quadpoints3d quadpoints);
2414 
2415 /* ------------------------------------------------------------
2416  * Greenhybrid
2417  * ------------------------------------------------------------ */
2418 
2476 HEADER_PREFIX void
2478  pcblock tree, uint m, uint l, real delta, real accur,
2479  quadpoints3d quadpoints);
2480 
2538 HEADER_PREFIX void
2540  pcblock tree, uint m, uint l, real delta, real accur,
2541  quadpoints3d quadpoints);
2542 
2582 HEADER_PREFIX void
2584  pccluster cc, pcblock tree, uint m, uint l, real delta, real accur,
2585  quadpoints3d quadpoints);
2586 
2587 /* ------------------------------------------------------------
2588  * ACA
2589  * ------------------------------------------------------------ */
2590 
2607 HEADER_PREFIX void
2609  pcblock tree, real accur);
2610 
2627 HEADER_PREFIX void
2629  pcblock tree, real accur);
2630 
2631 /* ------------------------------------------------------------
2632  * HCA
2633  * ------------------------------------------------------------ */
2634 
2693 HEADER_PREFIX void
2695  pcblock tree, uint m, real accur);
2696 
2697 /* ------------------------------------------------------------
2698  * Initializer functions for H2-matrix approximations
2699  * ------------------------------------------------------------ */
2700 
2717 HEADER_PREFIX void
2718 setup_h2matrix_recomp_bem3d(pbem3d bem, bool hiercomp, real accur_hiercomp);
2719 
2772 HEADER_PREFIX void
2774  pcclusterbasis cb, pcblock tree, uint m);
2775 
2851 HEADER_PREFIX void
2853  pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur,
2854  quadpoints3d quadpoints);
2855 
2931 HEADER_PREFIX void
2933  pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur,
2934  quadpoints3d quadpoints);
2935 
2936 /* ------------------------------------------------------------
2937  * Initializer functions for DH2-matrix approximations
2938  * ------------------------------------------------------------ */
2939 
2959 HEADER_PREFIX void
2961  pcdclusterbasis cb, pcdblock tree, uint m);
2962 
2984 HEADER_PREFIX void
2986  pcdclusterbasis cb, pcdblock tree, uint m);
2987 
3012 HEADER_PREFIX void
3014  pcdclusterbasis cb, pcdblock tree, uint m, ptruncmode tm, real eps);
3015 
3016 /* ------------------------------------------------------------
3017  * Fill amatrix
3018  * ------------------------------------------------------------ */
3019 
3028 HEADER_PREFIX void
3030 
3031 /* ------------------------------------------------------------
3032  * Fill hmatrix
3033  * ------------------------------------------------------------ */
3034 
3056 HEADER_PREFIX void
3058 
3087 HEADER_PREFIX void
3089 
3104 HEADER_PREFIX void
3106 
3127 HEADER_PREFIX void
3129 
3130 /* ------------------------------------------------------------
3131  * Fill H2-matrix
3132  * ------------------------------------------------------------ */
3133 
3159 HEADER_PREFIX void
3161 
3187 HEADER_PREFIX void
3189 
3215 HEADER_PREFIX void
3217 
3230 HEADER_PREFIX void
3232 
3257 HEADER_PREFIX void
3259 
3291 HEADER_PREFIX void
3293 
3322 HEADER_PREFIX void
3324 
3354 HEADER_PREFIX void
3356 
3376 HEADER_PREFIX void
3378  pdclusteroperator ro);
3379 
3399 HEADER_PREFIX void
3401  pdclusteroperator co);
3402 
3429 HEADER_PREFIX void
3432  pdclusteroperator bco, pcdblock broot);
3433 
3459 HEADER_PREFIX void
3461 
3474 HEADER_PREFIX void
3476 
3501 HEADER_PREFIX void
3503 
3504 /* ------------------------------------------------------------
3505  * Lagrange polynomials
3506  * ------------------------------------------------------------ */
3507 
3542 HEADER_PREFIX void
3544  pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V);
3545 
3583 HEADER_PREFIX void
3585  pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V);
3586 
3621 HEADER_PREFIX void
3623  pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V);
3624 
3662 HEADER_PREFIX void
3664  pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V);
3665 
3700 HEADER_PREFIX void
3702  pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V);
3703 
3742 HEADER_PREFIX void
3744  pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V);
3745 
3780 HEADER_PREFIX void
3782  pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V);
3783 
3822 HEADER_PREFIX void
3824  pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V);
3825 
3858 HEADER_PREFIX void
3860  pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V);
3861 
3896 HEADER_PREFIX void
3898  pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V);
3899 
3900 /* ------------------------------------------------------------
3901  * Some useful functions
3902  * ------------------------------------------------------------ */
3903 
3913 normL2_bem3d(pbem3d bem, boundary_func3d f, void *data);
3914 
3931 
3950 normL2diff_c_bem3d(pbem3d bem, pavector x, boundary_func3d f, void *data);
3951 
3968 
3987 normL2diff_l_bem3d(pbem3d bem, pavector x, boundary_func3d f, void *data);
3988 
4009 HEADER_PREFIX void
4011 
4032 HEADER_PREFIX void
4034 
4057 HEADER_PREFIX void
4059 
4091 HEADER_PREFIX void
4093 
4101 HEADER_PREFIX void
4103 
4138 void build_bem3d_cube_quadpoints(pcbem3d bem, const real a[3], const real b[3],
4139  const real delta, real (**Z)[3], real (**N)[3]);
4140 
4164 build_bem3d_rkmatrix(pccluster row, pccluster col, void* data);
4165 
4183 build_bem3d_amatrix(pccluster row, pccluster col, void* data);
4184 
4197 HEADER_PREFIX void
4199  psparsematrix *C2);
4200 
4203 #endif /* BEM3D_H_ */
4204 
field(* boundary_func3d)(const real *x, const real *n, void *data)
Definition: bem3d.h:197
pamatrix build_bem3d_amatrix(pccluster row, pccluster col, void *data)
For a block defined by row and col this function will compute the full submatrix and return an amatri...
real normL2_l_bem3d(pbem3d bem, pavector x)
Compute the -norm of some given function in terms of a piecewise linear basis .
Substructure containing callback functions to different types of kernel evaluations.
Definition: bem3d.h:768
void assemble_bem3d_dh2matrix_ortho_row_dclusterbasis(pcbem3d bem, pdclusterbasis rb, pdclusteroperator ro)
Function for filling a directional row cluster basis with orthogonal matrices.
void assemblehiercomp_bem3d_h2matrix(pbem3d bem, pblock b, ph2matrix G)
Fills an h2matrix with a predefined approximation technique using hierarchical recompression.
void setup_hmatrix_aprx_hca_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, real accur)
Approximate matrix block with hybrid cross approximation using Interpolation and ACA with full pivoti...
const kernelbem3d * pckernelbem3d
Definition: bem3d.h:132
pbem3d new_bem3d(pcsurface3d gr, basisfunctionbem3d row_basis, basisfunctionbem3d col_basis)
Main constructor for bem3d objects.
void setup_hmatrix_aprx_greenhybrid_row_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints3d quadpoints)
creating hmatrix approximation using green&#39;s method with row cluster connected with recompression tec...
void assemble_lc_far_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise linear basis functions for Tes...
void projectL2_bem3d_c_avector(pbem3d bem, boundary_func3d f, pavector w, void *data)
Computes the -projection of a given function using piecewise constant basis functions over the bounda...
void setup_hmatrix_aprx_green_col_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints3d quadpoints)
creating hmatrix approximation using green&#39;s method with column cluster .
void assemble_bem3d_dh2matrix_row_dclusterbasis(pcbem3d bem, pdclusterbasis rb)
This function computes the matrix entries for the nested dclusterbasis .
Enum value to represent piecewise linear basis functions.
Definition: bem3d.h:175
Definition: realavector.h:48
void integrate_bem3d_c_avector(pbem3d bem, boundary_func3d f, pavector w, void *data)
Computes the integral of a given function using piecewise constant basis functions over the boundary ...
void integrate_bem3d_l_avector(pbem3d bem, boundary_func3d f, pavector w, void *data)
Computes the integral of a given function using piecewise linear basis functions over the boundary of...
tri_list * ptri_list
Definition: bem3d.h:260
pvert_list new_vert_list(pvert_list next)
Create a new list to store a number of vertex indices.
basisfunctionbem3d col_basis
Type of basis function for dirichlet-data.
Definition: bem3d.h:317
void(* nearfield_far)(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N)
Computes nearfield entries of Galkerin matrices for disjoint clusters.
Definition: bem3d.h:405
Definition: avector.h:39
void setup_hmatrix_aprx_inter_col_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
This function initializes the bem3d object for approximating a matrix with tensorinterpolation within...
void(* farfield_u)(uint rname, uint cname, uint bname, pcbem3d bem)
Computes coupling matrix for uniform matrices.
Definition: bem3d.h:453
pkernelbem3d kernels
A collection of callback functions for different types of kernel evaluations.
Definition: bem3d.h:748
void assemble_bem3d_lagrange_l_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V)
This function will integrate Lagrange polynomials on the boundary domain using linear basis function ...
void setup_dh2matrix_aprx_inter_bem3d(pbem3d bem, pcdclusterbasis rb, pcdclusterbasis cb, pcdblock tree, uint m)
Initialize the bem3d object for an interpolation based -matrix approximation.
void fill_wave_bem3d(pcbem3d bem, const real(*X)[3], const real(*Y)[3], const real(*NX)[3], const real(*NY)[3], pamatrix V, pcreal dir, kernel_wave_func3d kernel)
Evaluate a modified kernel function at some points and for some direction .
void(* transfer_col)(uint cname, pcbem3d bem)
Computes the transfer matrices for a cluster .
Definition: bem3d.h:652
pclustergeometry build_bem3d_clustergeometry(pcbem3d bem, uint **idx, basisfunctionbem3d basis)
Creates a clustergeometry object for a BEM-Problem using the type of basis functions specified by bas...
void fill_dnz_col_c_bem3d(const uint *idx, const real(*Z)[3], const real(*N)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a normal derivative of a kernel function with respect to on the bounda...
void assemble_bem3d_lagrange_wave_c_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V)
This function will integrate modified Lagrange polynomials on the boundary domain using piecewise con...
void build_bem3d_curl_sparsematrix(pcbem3d bem, psparsematrix *C0, psparsematrix *C1, psparsematrix *C2)
Set up the surface curl operator.
struct _parbem3d parbem3d
Definition: bem3d.h:107
void assemble_ll_near_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise linear basis functions for bot...
void assemble_bem3d_lagrange_c_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V)
This function will integrate Lagrange polynomials on the boundary domain using piecewise constant bas...
Simple singly connected list to efficiently store a list of vertices.
Definition: bem3d.h:1323
void assemble_bem3d_amatrix(pbem3d bem, pamatrix G)
Assemble a dense matrix for some boundary integral operator given by bem.
Representation of -matrices.
Definition: h2matrix.h:48
field(* kernel_wave_func3d)(const real *x, const real *y, const real *nx, const real *ny, pcreal dir, void *data)
Evaluate a modified fundamental solution or its normal derivatives at points x and y...
Definition: bem3d.h:230
real * mass
Mass-matrix of reference basis functions.
Definition: bem3d.h:327
void projectL2_bem3d_l_avector(pbem3d bem, boundary_func3d f, pavector w, void *data)
Computes the -projection of a given function using linear basis functions over the boundary of a doma...
basisfunctionbem3d row_basis
Type of basis function for neumann-data.
Definition: bem3d.h:312
void setup_hmatrix_aprx_greenhybrid_mixed_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints3d quadpoints)
creating hmatrix approximation using green&#39;s method with row or column cluster connected with recompr...
void assemble_cc_far_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise constant basis functions for b...
Representation of -matrices.
Definition: hmatrix.h:49
void assemble_bem3d_dh2matrix_recomp_both_dclusterbasis(pcbem3d bem, pdclusterbasis rb, pdclusteroperator bro, pdclusterbasis cb, pdclusteroperator bco, pcdblock broot)
Function for filling both directional cluster bases with orthogonal matrices.
Representation of a cluster basis.
Definition: clusterbasis.h:40
void assemble_ll_far_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise linear basis functions for bot...
void assemble_bem3d_lagrange_amatrix(const real(*X)[3], pcrealavector px, pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V)
This function will evaluate the Lagrange polynomials in a given point set X and store the results in ...
void fill_bem3d(pcbem3d bem, const real(*X)[3], const real(*Y)[3], const real(*NX)[3], const real(*NY)[3], pamatrix V, kernel_func3d kernel)
Evaluate a kernel function at some points and .
unsigned uint
Unsigned integer type.
Definition: settings.h:70
void fill_col_l_bem3d(const uint *idx, const real(*Z)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a kernel function on the boundary domain in the second argument using p...
void assemble_bem3d_dn_lagrange_l_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V)
This function will integrate the normal derivatives of the Lagrange polynomials on the boundary domai...
pcsurface3d gr
Polyedric 3 dimensional surface mesh.
Definition: bem3d.h:297
void(* quadpoints3d)(pcbem3d bem, const real bmin[3], const real bmax[3], const real delta, real(**Z)[3], real(**N)[3])
Callback function type for parameterizations.
Definition: bem3d.h:153
void assemble_bem3d_h2matrix(pbem3d bem, ph2matrix G)
Fills a h2matrix with a predefined approximation technique.
ptri_list new_tri_list(ptri_list next)
Create a new list to store a number of triangles indices.
void setup_h2matrix_recomp_bem3d(pbem3d bem, bool hiercomp, real accur_hiercomp)
Enables hierarchical recompression for hmatrices.
void assemble_lc_near_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise linear basis functions for Tes...
bem3d * pbem3d
Definition: bem3d.h:79
double _Complex field
Field type.
Definition: settings.h:171
ptri_list next
Pointer to the next _tri_list element.
Definition: bem3d.h:1313
void(* transfer_row)(uint rname, pcbem3d bem)
Computes the transfer matrices for a cluster .
Definition: bem3d.h:573
pvert_list next
Pointer to the next _vert_list element.
Definition: bem3d.h:1327
pparbem3d par
Some helpers to make parallel computations possible.
Definition: bem3d.h:730
Representation of a clustergeometry object.
Definition: clustergeometry.h:45
void assemble_bem3d_dn_lagrange_wave_c_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V)
This function will integrate the normal derivatives of the modified Lagrange polynomials on the bound...
void assemble_bem3d_dn_lagrange_c_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcbem3d bem, pamatrix V)
This function will integrate the normal derivatives of the Lagrange polynomials on the boundary domai...
void assemble_bem3d_dh2matrix_ortho_col_dclusterbasis(pcbem3d bem, pdclusterbasis cb, pdclusteroperator co)
Function for filling a directional column cluster basis with orthogonal matrices. ...
void setup_h2matrix_aprx_inter_bem3d(pbem3d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m)
Initialize the bem3d object for approximating a h2matrix with tensorinterpolation.
vert_list * pvert_list
Definition: bem3d.h:270
_basisfunctionbem3d
Possible types of basis functions.
Definition: bem3d.h:163
Tree structure representing a -matrix.
Definition: dh2matrix.h:50
Dummy value, should never be used.
Definition: bem3d.h:167
parbem3d * pparbem3d
Definition: bem3d.h:111
void fill_dnz_row_c_bem3d(const uint *idx, const real(*Z)[3], const real(*N)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a normal derivative of a kernel function with respect to on the bounda...
psingquad2d sq
Quadrature rules used within BEM computation.
Definition: bem3d.h:307
void(* transfer_wave_row)(uint rname, pcbem3d bem)
Computes the transfer matrices for a cluster and directions . The father cluster is using the wave...
Definition: bem3d.h:599
void setup_vertex_to_triangle_map_bem3d(pbem3d bem)
initializes the field bem->v2t when using linear basis functions
real normL2_c_bem3d(pbem3d bem, pavector x)
Compute the -norm of some given function in terms of a piecewise constant basis .
Representation of a directional cluster operator.
Definition: dclusteroperator.h:33
void setup_hmatrix_aprx_green_mixed_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints3d quadpoints)
creating hmatrix approximation using green&#39;s method with one of row or column cluster ...
simple singly linked list to store unsigned inter values.
Definition: bem3d.h:1290
real normL2diff_c_bem3d(pbem3d bem, pavector x, boundary_func3d f, void *data)
Compute the difference norm of some given function in terms of a piecewise constant basis and some...
This struct collects all type of quadrature formulas needed by the computation of matrix entries with...
Definition: singquad2d.h:38
void setup_hmatrix_aprx_inter_row_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
This function initializes the bem3d object for approximating a matrix with tensorinterpolation within...
field(* kernel_func3d)(const real *x, const real *y, const real *nx, const real *ny, void *data)
Evaluate a fundamental solution or its normal derivatives at points x and y.
Definition: bem3d.h:210
void fill_row_l_bem3d(const uint *idx, const real(*Z)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a kernel function on the boundary domain in the first argument using pi...
Directional block tree.
Definition: dblock.h:45
void(* transfer_wave_wave_col)(uint cname, pcbem3d bem)
Computes the transfer matrices for a cluster and directions . Both, the father cluster and sons a...
Definition: bem3d.h:712
uint t
Index of a the triangle to be stored.
Definition: bem3d.h:1306
void fill_col_c_bem3d(const uint *idx, const real(*Z)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a kernel function on the boundary domain in the second argument using p...
Enum value to represent piecewise constant basis functions.
Definition: bem3d.h:171
void assemble_bem3d_lagrange_wave_l_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V)
This function will integrate modified Lagrange polynomials on the boundary domain using linear basis ...
void assemble_bem3d_farfield_dh2matrix(pbem3d bem, pdh2matrix G)
Fills a dh2matrix with a predefined approximation technique.
Define different strategies used by various truncation and compression algorithms for hmatrices and h...
Definition: truncation.h:43
void assemble_cl_near_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise constant basis functions for T...
void assemble_bem3d_farfield_hmatrix(pbem3d bem, pblock b, phmatrix G)
Fills the farfield blocks of a hmatrix with a predefined approximation technique. ...
Simple singly connected list to efficiently store a list of triangles.
Definition: bem3d.h:1304
field alpha
Boundary integral operator + mass matrix.
Definition: bem3d.h:332
void(* transfer_wave_col)(uint cname, pcbem3d bem)
Computes the transfer matrices for a cluster and directions . The father cluster is using the wave...
Definition: bem3d.h:678
void assemblecoarsen_bem3d_hmatrix(pbem3d bem, pblock b, phmatrix G)
Fills a hmatrix with a predefined approximation technique using coarsening strategy.
void fill_dnz_col_l_bem3d(const uint *idx, const real(*Z)[3], const real(*N)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a normal derivative of a kernel function with respect to on the bounda...
void fill_row_c_bem3d(const uint *idx, const real(*Z)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a kernel function on the boundary domain in the first argument using pi...
void(* farfield_wave_u)(uint rname, uint cname, uint bname, pcbem3d bem)
Computes coupling matrix for duniform matrices for some direction .
Definition: bem3d.h:477
void(* transfer_wave_wave_row)(uint rname, pcbem3d bem)
Computes the transfer matrices for a cluster and directions . Both, the father cluster and sons a...
Definition: bem3d.h:633
void setup_dh2matrix_aprx_inter_recomp_bem3d(pbem3d bem, pcdclusterbasis rb, pcdclusterbasis cb, pcdblock tree, uint m, ptruncmode tm, real eps)
Initialize the bem3d object for an interpolation based -matrix approximation with orthogonal directio...
Representation of a triangle surface mesh.
Definition: surface3d.h:45
pclustergeometry build_bem3d_const_clustergeometry(pcbem3d bem, uint **idx)
Creates a clustergeometry object for a BEM-Problem using piecewise constant basis functions...
void assemble_bem3d_farfield_h2matrix(pbem3d bem, ph2matrix G)
Fills a h2matrix with a predefined approximation technique.
void assemble_bem3d_lagrange_wave_amatrix(const real(*X)[3], pcrealavector px, pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V)
This function will evaluate the modified Lagrange polynomials in a given point set X and store the re...
uint data
Data element to be stored. here the number of a triangle.
Definition: bem3d.h:1292
void assemble_bem3d_nearfield_h2matrix(pbem3d bem, ph2matrix G)
Fills the nearfield part of a h2matrix.
pcluster build_bem3d_cluster(pcbem3d bem, uint clf, basisfunctionbem3d basis)
Creates a clustertree for specified basis functions.
#define HEADER_PREFIX
Prefix for function declarations.
Definition: settings.h:43
double real
real floating point type.
Definition: settings.h:97
void assemble_bem3d_h2matrix_row_clusterbasis(pcbem3d bem, pclusterbasis rb)
This function computes the matrix entries for the nested clusterbasis .
void del_vert_list(pvert_list vl)
Recursively deletes a vert_list.
void setup_hmatrix_aprx_greenhybrid_col_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints3d quadpoints)
creating hmatrix approximation using green&#39;s method with column cluster connected with recompression ...
void(* leaf_col)(uint cname, pcbem3d bem)
Computes the the matrix for a leaf cluster .
Definition: bem3d.h:534
void(* leaf_wave_row)(uint rname, pcbem3d bem)
Computes the the matrix for a leaf cluster and a direction .
Definition: bem3d.h:515
kernelbem3d * pkernelbem3d
Definition: bem3d.h:127
field kernel_const
A constant factor extracted from the kernel function to speed up quadrature.
Definition: bem3d.h:344
pclustergeometry build_bem3d_linear_clustergeometry(pcbem3d bem, uint **idx)
Creates a clustergeometry object for a BEM-Problem using linear basis functions.
pvert_list vl
a List of vertices whose basis functions are connected with this triangle and are needed within the c...
Definition: bem3d.h:1311
void(* farfield_rk)(pccluster rc, uint rname, pccluster cc, uint cname, pcbem3d bem, prkmatrix R)
Computes rank-k-approximations of a given block.
Definition: bem3d.h:428
enum _basisfunctionbem3d basisfunctionbem3d
Definition: bem3d.h:181
void assemble_bem3d_dn_lagrange_wave_l_amatrix(const uint *idx, pcrealavector px, pcrealavector py, pcrealavector pz, pcreal dir, pcbem3d bem, pamatrix V)
This function will integrate the normal derivatives of the modified Lagrange polynomials on the bound...
void assemble_bem3d_h2matrix_col_clusterbasis(pcbem3d bem, pclusterbasis cb)
This function computes the matrix entries for the nested clusterbasis .
void assemble_bem3d_dh2matrix(pbem3d bem, pdh2matrix G)
Fills a dh2matrix with a predefined approximation technique.
void assemble_cl_far_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise constant basis functions for T...
void assemble_bem3d_hmatrix(pbem3d bem, pblock b, phmatrix G)
Fills a hmatrix with a predefined approximation technique.
void setup_h2matrix_aprx_greenhybrid_bem3d(pbem3d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur, quadpoints3d quadpoints)
Initialize the bem3d object for approximating a h2matrix with green&#39;s method and ACA based recompress...
uint v
Index of the vertex to be stored.
Definition: bem3d.h:1325
Representation of block trees.
Definition: block.h:48
plistnode * v2t
This field describes the mapping from the vertices of a geometry to the triangles they belong to...
Definition: bem3d.h:356
void(* leaf_row)(uint rname, pcbem3d bem)
Computes the the matrix for a leaf cluster .
Definition: bem3d.h:496
void setup_hmatrix_aprx_green_row_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints3d quadpoints)
creating hmatrix approximation using green&#39;s method with row cluster .
void del_tri_list(ptri_list tl)
Recursively deletes a tri_list.
Representation of cluster trees.
Definition: cluster.h:40
prkmatrix build_bem3d_rkmatrix(pccluster row, pccluster col, void *data)
For a block defined by row and col this function will build a rank-k-approximation and return the rkm...
void(* nearfield)(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N)
Computes nearfield entries of Galerkin matrices.
Definition: bem3d.h:380
Representation of a low-rank matrix in factorized form .
Definition: rkmatrix.h:36
void del_bem3d(pbem3d bem)
Destructor for bem3d objects.
void setup_hmatrix_recomp_bem3d(pbem3d bem, bool recomp, real accur_recomp, bool coarsen, real accur_coarsen)
Initialize the bem object for on the fly recompression techniques.
real normL2_bem3d(pbem3d bem, boundary_func3d f, void *data)
Compute the -norm of some given function f.
const bem3d * pcbem3d
Definition: bem3d.h:83
Main container object for computation of BEM related matrices and vectors.
Definition: bem3d.h:290
void setup_h2matrix_aprx_greenhybrid_ortho_bem3d(pbem3d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur, quadpoints3d quadpoints)
Initialize the bem3d object for approximating a h2matrix with green&#39;s method and ACA based recompress...
Representation of a sparse matrix in compressed row format.
Definition: sparsematrix.h:42
void assemble_cc_near_bem3d(const uint *ridx, const uint *cidx, pcbem3d bem, bool ntrans, pamatrix N, kernel_func3d kernel)
Compute general entries of a boundary integral operator with piecewise constant basis functions for b...
Representation of a matrix as an array in column-major order.
Definition: amatrix.h:43
void setup_dh2matrix_aprx_inter_ortho_bem3d(pbem3d bem, pcdclusterbasis rb, pcdclusterbasis cb, pcdblock tree, uint m)
Initialize the bem3d object for an interpolation based -matrix approximation with orthogonal directio...
void assemble_bem3d_dh2matrix_col_dclusterbasis(pcbem3d bem, pdclusterbasis cb)
This function computes the matrix entries for the nested dclusterbasis .
void fill_dnz_row_l_bem3d(const uint *idx, const real(*Z)[3], const real(*N)[3], pcbem3d bem, pamatrix V, kernel_func3d kernel)
This function will integrate a normal derivative of a kernel function with respect to on the bounda...
const parbem3d * pcparbem3d
Definition: bem3d.h:115
void(* leaf_wave_col)(uint cname, pcbem3d bem)
Computes the the matrix for a leaf cluster and a direction .
Definition: bem3d.h:554
void setup_hmatrix_aprx_inter_mixed_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
This function initializes the bem3d object for approximating a matrix with tensorinterpolation within...
void build_bem3d_cube_quadpoints(pcbem3d bem, const real a[3], const real b[3], const real delta, real(**Z)[3], real(**N)[3])
Generating quadrature points, weights and normal vectors on a cube parameterization.
aprxbem3d * paprxbem3d
Definition: bem3d.h:96
void assemble_bem3d_nearfield_hmatrix(pbem3d bem, pblock b, phmatrix G)
Fills the nearfield blocks of a hmatrix.
paprxbem3d aprx
A collection of necessary data structures for approximating matrices with various techniques...
Definition: bem3d.h:721
struct _aprxbem3d aprxbem3d
Definition: bem3d.h:92
void assemble_bem3d_nearfield_dh2matrix(pbem3d bem, pdh2matrix G)
Fills the nearfield part of a dh2matrix.
plistnode next
Pointer to the next _listnode element.
Definition: bem3d.h:1294
void setup_hmatrix_aprx_aca_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, real accur)
Approximate matrix block with ACA using full pivoting.
Representation of a directional cluster basis.
Definition: dclusterbasis.h:45
const aprxbem3d * pcaprxbem3d
Definition: bem3d.h:100
field k
Wavenumber for Helmholtz type problems, possibly complex valued.
Definition: bem3d.h:338
listnode * plistnode
Definition: bem3d.h:250
real normL2diff_l_bem3d(pbem3d bem, pavector x, boundary_func3d f, void *data)
Compute the difference norm of some given function in terms of a piecewise linear basis and some o...
void setup_hmatrix_aprx_paca_bem3d(pbem3d bem, pccluster rc, pccluster cc, pcblock tree, real accur)
Approximate matrix block with ACA using partial pivoting.
const real * pcreal
Pointer to constant real array.
Definition: settings.h:148