H2Lib  3.0
Data Structures | Typedefs | Enumerations | Functions
bem2d.h File Reference
#include <stdio.h>
#include "settings.h"
#include "amatrix.h"
#include "gaussquad.h"
#include "factorizations.h"
#include "cluster.h"
#include "clustergeometry.h"
#include "hmatrix.h"
#include "clusterbasis.h"
#include "h2matrix.h"
#include "hcoarsen.h"
#include "h2update.h"
#include "aca.h"
#include "curve2d.h"
#include "singquad1d.h"

Go to the source code of this file.

Data Structures

struct  _bem2d
 Main container object for computation of BEM related matrices and vectors. More...
 
struct  _kernelbem2d
 Substructure containing callback functions to different types of kernel evaluations. More...
 

Typedefs

typedef struct _bem2d bem2d
 
typedef bem2dpbem2d
 
typedef const bem2dpcbem2d
 
typedef struct _aprxbem2d aprxbem2d
 
typedef aprxbem2dpaprxbem2d
 
typedef const aprxbem2dpcaprxbem2d
 
typedef struct _parbem2d parbem2d
 
typedef parbem2dpparbem2d
 
typedef const parbem2dpcparbem2d
 
typedef struct _kernelbem2d kernelbem2d
 
typedef kernelbem2dpkernelbem2d
 
typedef const kernelbem2dpckernelbem2d
 
typedef void(* quadpoints2d) (pcbem2d bem, const real bmin[2], const real bmax[2], const real delta, real(**Z)[2], real(**N)[2])
 Callback function type for parameterizations. More...
 
typedef enum _basisfunctionbem2d basisfunctionbem2d
 
typedef field(* boundary_func2d) (const real *x, const real *n)
 

Enumerations

enum  _basisfunctionbem2d { BASIS_NONE_BEM2D = 0, BASIS_CONSTANT_BEM2D = 'c', BASIS_LINEAR_BEM2D = 'l' }
 Possible types of basis functions. More...
 

Functions

pbem2d new_bem2d (pccurve2d gr)
 Main constructor for bem2d objects. More...
 
void del_bem2d (pbem2d bem)
 Destructor for bem2d objects. More...
 
pclustergeometry build_bem2d_const_clustergeometry (pcbem2d bem, uint **idx)
 Creates a clustergeometry object for a BEM-Problem using piecewise constant basis functions. More...
 
pclustergeometry build_bem2d_linear_clustergeometry (pcbem2d bem, uint **idx)
 Creates a clustergeometry object for a BEM-Problem using linear basis functions. More...
 
pclustergeometry build_bem2d_clustergeometry (pcbem2d bem, uint **idx, basisfunctionbem2d basis)
 Creates a clustergeometry object for a BEM-Problem using the type of basis functions specified by basis. More...
 
pcluster build_bem2d_cluster (pcbem2d bem, uint clf, basisfunctionbem2d basis)
 Creates a clustertree for specified basis functions. More...
 
void setup_hmatrix_recomp_bem2d (pbem2d bem, bool recomp, real accur_recomp, bool coarsen, real accur_coarsen)
 Initialize the bem object for on the fly recompression techniques. More...
 
void setup_hmatrix_aprx_inter_row_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
 This function initializes the bem2d object for approximating a matrix with tensorinterpolation within the row cluster. More...
 
void setup_hmatrix_aprx_inter_col_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
 This function initializes the bem2d object for approximating a matrix with tensorinterpolation within the column cluster. More...
 
void setup_hmatrix_aprx_inter_mixed_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m)
 This function initializes the bem2d object for approximating a matrix with tensorinterpolation within one of row or column cluster. More...
 
void setup_hmatrix_aprx_green_row_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints2d quadpoints)
 creating hmatrix approximation using green's method with row cluster . More...
 
void setup_hmatrix_aprx_green_col_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints2d quadpoints)
 creating hmatrix approximation using green's method with column cluster . More...
 
void setup_hmatrix_aprx_green_mixed_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, quadpoints2d quadpoints)
 creating hmatrix approximation using green's method with one of row or column cluster . More...
 
void setup_hmatrix_aprx_greenhybrid_row_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
 creating hmatrix approximation using green's method with row cluster connected with recompression technique using ACA. More...
 
void setup_hmatrix_aprx_greenhybrid_col_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
 creating hmatrix approximation using green's method with column cluster connected with recompression technique using ACA. More...
 
void setup_hmatrix_aprx_greenhybrid_mixed_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
 creating hmatrix approximation using green's method with row or column cluster connected with recompression technique using ACA. More...
 
void setup_hmatrix_aprx_aca_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, real accur)
 Approximate matrix block with ACA using full pivoting. More...
 
void setup_hmatrix_aprx_paca_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, real accur)
 Approximate matrix block with ACA using partial pivoting. More...
 
void setup_hmatrix_aprx_hca_bem2d (pbem2d bem, pccluster rc, pccluster cc, pcblock tree, uint m, real accur)
 Approximate matrix block with hybrid cross approximation using Interpolation and ACA with full pivoting. More...
 
void setup_h2matrix_recomp_bem2d (pbem2d bem, bool hiercomp, real accur_hiercomp)
 Enables hierarchical recompression for hmatrices. More...
 
void setup_h2matrix_aprx_inter_bem2d (pbem2d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m)
 Initialize the bem2d object for approximating a h2matrix with tensorinterpolation. More...
 
void setup_h2matrix_aprx_greenhybrid_bem2d (pbem2d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
 Initialize the bem2d object for approximating a h2matrix with green's method and ACA based recompression. More...
 
void setup_h2matrix_aprx_greenhybrid_ortho_bem2d (pbem2d bem, pcclusterbasis rb, pcclusterbasis cb, pcblock tree, uint m, uint l, real delta, real accur, quadpoints2d quadpoints)
 Initialize the bem2d object for approximating a h2matrix with green's method and ACA based recompression. The resulting clusterbasis will be orthogonal. More...
 
void assemble_bem2d_hmatrix (pbem2d bem, pblock b, phmatrix G)
 Fills an hmatrix with a predefined approximation technique. More...
 
void assemblecoarsen_bem2d_hmatrix (pbem2d bem, pblock b, phmatrix G)
 Fills an hmatrix with a predefined approximation technique using coarsening strategy. More...
 
void assemble_bem2d_h2matrix_row_clusterbasis (pcbem2d bem, pclusterbasis rb)
 This function computes the matrix entries for the nested clusterbasis $ V_t \, , t \in \mathcal T_{\mathcal I} $. More...
 
void assemble_bem2d_h2matrix_col_clusterbasis (pcbem2d bem, pclusterbasis cb)
 This function computes the matrix entries for the nested clusterbasis $ W_s \, , s \in \mathcal T_{\mathcal J} $. More...
 
void assemble_bem2d_h2matrix (pbem2d bem, pblock b, ph2matrix G)
 Fills an h2matrix with a predefined approximation technique. More...
 
void assemblehiercomp_bem2d_h2matrix (pbem2d bem, pblock b, ph2matrix G)
 Fills an h2matrix with a predefined approximation technique using hierarchical recompression. More...
 
void assemble_bem2d_lagrange_const_amatrix (const uint *idx, pcavector px, pcavector py, pcbem2d bem, pamatrix V)
 This function will integrate Lagrange polynomials on the boundary domain using piecewise constant basis function and store the results in a matrix V. More...
 
void assemble_bem2d_dn_lagrange_const_amatrix (const uint *idx, pcavector px, pcavector py, pcbem2d bem, pamatrix V)
 This function will integrate the normal derivatives of the Lagrange polynomials on the boundary domain using piecewise constant basis function and store the results in a matrix V. More...
 
void assemble_bem2d_lagrange_amatrix (const real(*X)[2], pcavector px, pcavector py, pamatrix V)
 This function will evaluate the Lagrange polynomials in a given point set X and store the results in a matrix V. More...
 
void projectl2_bem2d_const_avector (pbem2d bem, field(*rhs)(const real *x, const real *n), pavector f)
 Computes the $ L_2 $-projection of a given function using piecewise constant basis functions. More...
 
void build_bem2d_rect_quadpoints (pcbem2d bem, const real a[2], const real b[2], const real delta, real(**Z)[2], real(**N)[2])
 Generating quadrature points, weights and normal vectors on a rectangular parameterization. More...
 
void build_bem2d_circle_quadpoints (pcbem2d bem, const real a[2], const real b[2], const real delta, real(**Z)[2], real(**N)[2])
 Generating quadrature points, weights and normal vectors on a circular parameterization. More...
 
prkmatrix build_bem2d_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 rkmatrix. More...
 
pamatrix build_bem2d_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 amatrix. More...
 

Detailed Description

Author
Sven Christophersen
Date
2011