nlfem
Functions
integration.cpp File Reference

Go to the source code of this file.

Functions

int integrate (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_weakSingular (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_fullyContained (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_retriangulate (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_retriangulate_unysmm (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_exact (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_baryCenter (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_subSuperSetBalls (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int integrate_baryCenterRT (const ElementType &aT, const ElementType &bT, const QuadratureType &quadRule, const MeshType &mesh, const ConfigurationType &conf, bool is_firstbfslayer, double *termLocal, double *termNonloc, double *termLocalPrime, double *termNonlocPrime)
 
int method_subSuperSetBalls (const double *x, const ElementType &T, const MeshType &mesh)
 This truncation method returns the number of vertices of triangle T which interact with x_center, i.e. have a l2 distance smaller than delta. More...
 
int method_baryCenter (const double *x_center, const ElementType &T, const MeshType &mesh, double *reTriangle_list, int is_placePointOnCap)
 
int method_retriangulate1D (const double *xCenter, const ElementType &T, const MeshType &mesh, double *reTriangleList, int isPlacePointOnCap)
 
int method_retriangulate (const double *xCenter, const ElementType &T, const MeshType &mesh, double *reTriangleList, int isPlacePointOnCap)
 
int method_retriangulateInfty (const double *xCenter, const ElementType &T, const MeshType &mesh, double *reTriangleList, int isPlacePointOnCap)
 
int method_exact (const double *xCenter, const ElementType &T, const MeshType &mesh, double *reTriangleList, double *capsList, double *capsWeights, int *prtnCaps)
 

Detailed Description

Author
Manuel Klar

Definition in file integration.cpp.

Function Documentation

◆ integrate()

int integrate ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

All integration routines share a common signature. The central difference between different routines is the way they handle the truncation of the domain of integration. In almost all cases, the truncation of the inner triangle bT is performed for each quadrature point in the outer triangle aT. In addition for singular kernels special care is needed. Due to the necessity of transformations of the domain of integration those integration routines are tied to a kernel.

The integration routine changes the data in *termLocal and *termNonloc to. The arrays have to be zero-initialized.

termLocal = int_aT phiA(x) phiB(x) int_bT ker(x,y) dy dx,\n
termNonloc = int_aT phiA(x) int_bT phiB(y) ker'(y,x) dy dx.
termLocalPrime = int_aT  int_bT phiA(y) phiB(y) ker'(x,y) dy dx,\n
termNonlocPrime = int_aT phiA(x) int_bT phiB(y) ker(y,x) dy dx,

where ker'(x,y) = ker(y,x). Please note that the nonlocal term has to be subtracted, while the local term has to be added to the stiffness matrix.

Parameters
aTTriangle of the outer integral.
bTTriangle of the inner integral.
quadRuleQuadrature rule.
meshMesh.
confConfiguration.
is_firstbfslayerSwitch to tell whether the integration is happening in the first layer of the breadth first search. This variable is true only if the kernel is singular. In that case the integrals between aT and its immediate neighbours have to be handled with special care.
termLocalThis term contains the local part of the integral
termNonlocThis term contains the nonlocal part of the integral

Definition at line 33 of file integration.cpp.

◆ integrate_baryCenter()

int integrate_baryCenter ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routines uses method_baryCenter() to truncate the inner domain bT. See integrate() for general information about the integration routines.

Definition at line 877 of file integration.cpp.

◆ integrate_baryCenterRT()

int integrate_baryCenterRT ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routines uses method_retriangulate() to truncate the outer domain bT. See integrate() for general information about the integration routines.

Definition at line 1121 of file integration.cpp.

◆ integrate_exact()

int integrate_exact ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routines uses method_exact() to truncate the inner domain bT. See integrate() for general information about the integration routines.

Definition at line 672 of file integration.cpp.

◆ integrate_fullyContained()

int integrate_fullyContained ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routine performs no truncation of any domain. It can be applied to integrate in cases where no truncation is needed. See integrate() for general information about the integration routines.

Definition at line 256 of file integration.cpp.

◆ integrate_retriangulate()

int integrate_retriangulate ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routines uses method_retriangulate() to truncate the inner domain bT. See integrate() for general information about the integration routines.

Definition at line 311 of file integration.cpp.

◆ integrate_retriangulate_unysmm()

int integrate_retriangulate_unysmm ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routines uses method_retriangulate() to truncate the inner domain bT. See integrate() for general information about the integration routines. The approx balls are not symmetrified.

Definition at line 508 of file integration.cpp.

◆ integrate_subSuperSetBalls()

int integrate_subSuperSetBalls ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routines uses method_subSuperSetBalls() to truncate the inner domain bT. See integrate() for general information about the integration routines.

Definition at line 988 of file integration.cpp.

◆ integrate_weakSingular()

int integrate_weakSingular ( const ElementType aT,
const ElementType bT,
const QuadratureType &  quadRule,
const MeshType mesh,
const ConfigurationType conf,
bool  is_firstbfslayer,
double *  termLocal,
double *  termNonloc,
double *  termLocalPrime,
double *  termNonlocPrime 
)

This integration routine handles the weak singularity close to the origin. Due to the transformation of the set aT x bT truncations are not possible. However, the truncation is necessary only close to the boundary of the interaction set. There, this function is used in combination with integration routines which truncate the domain.

Definition at line 48 of file integration.cpp.

◆ method_baryCenter()

int method_baryCenter ( const double *  x_center,
const ElementType T,
const MeshType mesh,
double *  reTriangle_list,
int  is_placePointOnCap 
)

This truncation method checks whether the distance of the point point x_center (of aT) to the barycenter of bT is smaller than delta. The function writes the list of triangles which are obtained from the truncation into reTriangle_list. In this case the list is either untouched or contains bT itself (there is no retriangulation).

Parameters
x_centerPhysical quadrature point. This point is obtained by mapping a quadrature point of the reference element to the triangle aT of the outer domain.
TTriangle of the inner integral.
meshMesh
reTriangle_listList of triangles which are obtained from the truncation.
is_placePointOnCap(Unused)
Returns
0 if there is no interaction, -1 otherwise.

Definition at line 1228 of file integration.cpp.

◆ method_exact()

int method_exact ( const double *  xCenter,
const ElementType T,
const MeshType mesh,
double *  reTriangleList,
double *  capsList,
double *  capsWeights,
int *  prtnCaps 
)

This truncation method retriangulates the triangle bT depending on the distance of xCenter to bT w.r.t. an exact quadrature rule for caps. The function writes the list of triangles which are obtained from the trunaction into reTriangle_list, and additionally provides one point in the cap center with the corresponding weight.

Parameters
x_centerPhysical quadrature point. This point is obtained by mapping a quadrature point of the reference element to the triangle aT of the outer domain.
TTriangle of the inner integral.
meshMesh
reTriangle_listList of triangles which are obtained from the retriangulation.
capsListList of cap center points
capsWeightsCorresponding weights of the quadrature
prtnCapsNumber of caps.
Returns

Definition at line 1754 of file integration.cpp.

◆ method_retriangulate()

int method_retriangulate ( const double *  xCenter,
const ElementType T,
const MeshType mesh,
double *  reTriangleList,
int  isPlacePointOnCap 
)

This truncation method retriangulates the triangle bT depending on the distance of xCenter to bT w.r.t. the L2-norm ball. The function writes the list of triangles which are obtained from the truncation into reTriangle_list.

Parameters
x_centerPhysical quadrature point. This point is obtained by mapping a quadrature point of the reference element to the triangle aT of the outer domain.
TTriangle of the inner integral.
meshMesh
reTriangle_listList of triangles which are obtained from the retriangulation.
is_placePointOnCapSwitch for with caps integration.
Returns
0 if there is no interaction, -1 otherwise.

Definition at line 1379 of file integration.cpp.

◆ method_retriangulate1D()

int method_retriangulate1D ( const double *  xCenter,
const ElementType T,
const MeshType mesh,
double *  reTriangleList,
int  isPlacePointOnCap 
)

This method allows a truncation on 1D domains. In 1D the truncated interaction set is always an interval.

Parameters
x_centerPhysical quadrature point.
TA physical element of the inner integral.
meshMesh
reTriangle_listSubinterval which is obtained from the truncation.
(unused)
Returns
0 if there is no interaction, -1 otherwise.

Definition at line 1359 of file integration.cpp.

◆ method_retriangulateInfty()

int method_retriangulateInfty ( const double *  xCenter,
const ElementType T,
const MeshType mesh,
double *  reTriangleList,
int  isPlacePointOnCap 
)

This truncation method retriangulates the triangle bT depending on the distance of xCenter to bT w.r.t. the L-infinity-norm ball. The function writes the list of triangles which are obtained from the trunaction into reTriangle_list.

Parameters
x_centerPhysical quadrature point. This point is obtained by mapping a quadrature point of the reference element to the triangle aT of the outer domain.
TTriangle of the inner integral.
meshMesh
reTriangle_listList of triangles which are obtained from the retriangulation.
is_placePointOnCapSwitch for with caps integration.
Returns
0 if there is no interaction, -1 otherwise.

Definition at line 1611 of file integration.cpp.

◆ method_subSuperSetBalls()

int method_subSuperSetBalls ( const double *  x_center,
const ElementType T,
const MeshType mesh 
)

This truncation method returns the number of vertices of triangle T which interact with x_center, i.e. have a l2 distance smaller than delta.

Parameters
x_centerPhysical quadrature point. This point is obtained by mapping a quadrature point of the reference element to the triangle aT of the outer domain.
TTriangle of the inner integral.
meshMesh
Returns
0 if there is no interaction. 1,2 or 3 otherwise.

Definition at line 1219 of file integration.cpp.