nlfem
|
Go to the source code of this file.
Functions | |
void | lookup_configuration (ConfigurationType &conf, int verbose=0) |
void | par_evaluateMass (double *vd, const double *ud, long *Elements, const long *ElementLabels, const double *Verts, const long *VertexLabels, int K_Omega, int J, int nP, double *P, const double *dx, const int dim, const int outdim, const int is_DiscontinuousGalerkin) |
Evaluate the mass matrix v = Mu. More... | |
void | par_assemble (const string compute, const string path_spAd, const string path_fd, const int K_Omega, const int K, const long *ptrTriangles, const long *ptrLabelTriangles, const double *ptrVerts, const long *ptrLabelVerts, const int nE, const int nE_Omega, const int nV, const int nV_Omega, const double *Px, const int nPx, const double *dx, const double *Py, const int nPy, const double *dy, const double sqdelta, const int is_DiscontinuousGalerkin, const int is_NeumannBoundary, const string str_model_kernel, const string str_model_f, const string str_integration_method_remote, const string str_integration_method_close, const int is_PlacePointOnCap, const int dim, const int outdim, const long *ptrTheta_indices, const long *ptrTheta_indptr, const double *ptrTheta_data, const long nTheta, const double *Pg, const int degree, const double *dg, double maxDiameter, double fractional_s, int is_fullConnectedComponentSearch, int verbose) |
Parallel assembly of nonlocal operator using a finite element approach. This function is a wrapper for the functions par_system() and par_forcing(). It allows to call the C++ functions from Cython. More... | |
void | par_system (MeshType &mesh, QuadratureType &quadRule, ConfigurationType &conf) |
Multi threaded assembly of a nonlocal operator using a finite element approach. Kernel functions can be defined in model see model_kernel() for more information. More... | |
void | par_forcing (MeshType &mesh, QuadratureType &quadRule, ConfigurationType &conf) |
Multi threaded assembly of the forcing term. Forcing functions can be defined in model see model_f() for more information. The resulting vector is written to the disk. More... | |
Variables | |
map< string, void(*)(const double *x, const MeshType &mesh, double *forcing_out)> | lookup_f |
map< string, void(*)(const double *x, const ElementType &aT, const double *y, const ElementType &bT, const MeshType &mesh, double *kernel_val)> | lookup_kernel |
map< string, int(*)(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)> | lookup_integrate |
Contains the main assembly algorithm for nonlocal stiffnes matrix and forcing function.
Definition in file Cassemble.cpp.
void lookup_configuration | ( | ConfigurationType & | conf, |
int | verbose = 0 |
||
) |
This function looks up the configuration.
conf | ConfigurationType |
verbose |
Definition at line 101 of file Cassemble.cpp.
void par_assemble | ( | string | compute, |
string | path_spAd, | ||
string | path_fd, | ||
int | K_Omega, | ||
int | K, | ||
const long * | ptrTriangles, | ||
const long * | ptrLabelTriangles, | ||
const double * | ptrVerts, | ||
const long * | ptrLabelVerts, | ||
int | nE, | ||
int | nE_Omega, | ||
int | nV, | ||
int | nV_Omega, | ||
const double * | Px, | ||
int | nPx, | ||
const double * | dx, | ||
const double * | Py, | ||
int | nPy, | ||
const double * | dy, | ||
double | sqdelta, | ||
int | is_DiscontinuousGalerkin, | ||
int | is_NeumannBoundary, | ||
string | str_model_kernel, | ||
string | str_model_f, | ||
string | str_integration_method_remote, | ||
string | str_integration_method_close, | ||
int | is_PlacePointOnCap, | ||
int | dim, | ||
int | outdim, | ||
const long * | ptrTheta_indices, | ||
const long * | ptrTheta_indptr, | ||
const double * | ptrTheta_data, | ||
long | nTheta = 0 , |
||
const double * | Pg = nullptr , |
||
int | degree = 0 , |
||
const double * | dg = nullptr , |
||
double | maxDiameter = 0.0 , |
||
double | fractional_s = -1.0 , |
||
int | is_fullConnectedComponentSearch = 0 , |
||
int | verbose = 0 |
||
) |
Parallel assembly of nonlocal operator using a finite element approach. This function is a wrapper for the functions par_system() and par_forcing(). It allows to call the C++ functions from Cython.
Any 2-dimensional array is handed over by a pointer only. The expected shape of the input arrays is given in the parameter list where the underlying data is expected to be in C-contiguous (i.e. row-major) order. We denote the dimension of the domain by d. See par_system() and par_forcing() for mor details on the assembly.
compute | string "system", "forcing", "systemforcing". This string determines whether the stiffnes matrix, load or both should be computed. |
path_spAd | Save path for sparse matrix (armadillo binary) |
path_fd | Save path for forcing (armadillo binary) |
K_Omega | Number of rows of the stiffness matrix A. Example: If you want to solve a scalar equation using continuous Galerkin basis functions then K_Omega is equal to the number of basisfunctions (nV_Omega) which are not part of the Dirichlet boundary. If your problem is not scalar as for example in peridynamics then K_Omega = nV_Omega*outdim. |
K | Number of Columns of the stiffness matrix A |
ptrTriangles | (nE, d+1) Pointer to elements. Label 1 for elements in Omega, Label 2 for elements in OmegaI. |
ptrLabelTriangles | (nE,) Pointer to element labels |
ptrVerts | (L, d) Pointer to vertices |
ptrLabelVerts | Pointer to vertex labels |
nE | Number of elements |
nE_Omega | Number of elements in Omega |
nV | Number of vertices |
nV_Omega | Number of vertices in Omega |
Px | (nPx, d) Pointer to quadrature points for the outer integral |
nPx | Number of quadrature points in the outer integral |
dx | (nPx,) Pointer to quadrature weights of the outer integral |
Py | Number of quadrature points in the inner integral |
nPy | (nPy, d) Pointer to quadrature points for the inner integral |
dy | (nPx,) Pointer to quadrature weights of the inner integral |
sqdelta | Squared delta |
is_DiscontinuousGalerkin | Switch for discontinuous Galerkin |
is_NeumannBoundary | Switch of Neumann Boundary Conditions |
str_model_kernel | Name of kernel |
str_model_f | Name of right hand side |
str_integration_method_remote | Name of integration method |
is_PlacePointOnCap | Switch for withcaps parameter in retriangulation |
dim | Dimension of the domain |
outdim | Dimension of the solution space (e.g. 1 for scalar problems, dim for linear elasticity) |
ptrTheta_indices | Varying coefficients Theta, indices of CSR-format (optional) |
ptrTheta_indptr | Varying coefficients Theta, index pointer of CSR-format (optional) |
ptrTheta_data | Varying coefficients Theta, data of CSR-format (optional) |
nTheta | Number of rows, that is of subdomains |
Pg | (nPg, dim^2) Quadrature points for tensor Gauss quadrature (optional, needed for singular kernels). |
nPg | (tensorGaussDegree^dim,) Number of quadrature points. |
dg | (nPg,) Weights for tensor Gauss quadrature. |
maxDiameter | Maximal diameter of finite elements (optional). Might increase speed of retriangulation if provided. |
fractional_s | Degree of fractional kernel (d=2 only). Required for the correct choice of the integration routine. Default is s=-1.0, which corresponds to a kernel with no singularity. |
is_fullConnectedComponentSearch | If equal to 1 the breadth first traversal is not truncated. Default is 0. |
verbose | switch for verbose mode. |
Definition at line 388 of file Cassemble.cpp.
void par_evaluateMass | ( | double * | vd, |
const double * | ud, | ||
long * | Elements, | ||
const long * | ElementLabels, | ||
const double * | Verts, | ||
const long * | VertexLabels, | ||
int | K_Omega, | ||
int | J, | ||
int | nP, | ||
double * | P, | ||
const double * | dx, | ||
int | dim, | ||
int | outdim, | ||
int | is_DiscontinuousGalerkin | ||
) |
Evaluate the mass matrix v = Mu.
vd | Pointer to the first entry of the output vector. |
ud | Pointer to the first entry of the input vector. |
Elements | List of elements of a finite element triangulation (CSR-format, row major order). |
ElementLabels | List of element Labels. |
Verts | List of vertices (row major order). |
VertexLabels | List of vertex Labels. |
K_Omega | Number of rows and columns in M. Example: If you use continuous Galerkin basis functions and want to solve a scalar problem K_Omega = J. |
J | Number of elements in the triangulation. |
nP | Number of quadrature points in the outer integral |
P | (nPx, d) Pointer to quadrature points. |
dx | (nPx,) Pointer to quadrature weights. |
dim | Dimension of the domain Omega (2 or 3). |
outdim | Dimension of the kernel |
is_DiscontinuousGalerkin | switch for discontinuous Galerkin ansatz spaces |
Definition at line 193 of file Cassemble.cpp.
void par_forcing | ( | MeshType & | mesh, |
QuadratureType & | quadRule, | ||
ConfigurationType & | conf | ||
) |
Multi threaded assembly of the forcing term. Forcing functions can be defined in model see model_f() for more information. The resulting vector is written to the disk.
mesh | The mesh is of type MeshType and contains all information about the finite element descretization. See MeshStruct for more information. |
quadRule | Quadrature rules for inner, and outer elements as well as for the singular kernels. |
conf | General configuration, namely kernel, and forcing functions, as well as integration method. |
Definition at line 831 of file Cassemble.cpp.
void par_system | ( | MeshType & | mesh, |
QuadratureType & | quadRule, | ||
ConfigurationType & | conf | ||
) |
Multi threaded assembly of a nonlocal operator using a finite element approach. Kernel functions can be defined in model see model_kernel() for more information.
This function assembles the stiffness matrix corresponding to the operator
It loops over all elements aT with element Label != 0 and adds up their contribution to the global stiffness matrix. The integration starts with the domain aT x aT and then proceeds with the neighboring elements of aT in the mesh until the interaction domain is exceeded. For each pair aT, bT an integration routine (options are e.g. integrate_retriangulate(), integrate_baryCenter(), integrate_baryCenterRT()) is called. If the total contribution is 0 the elements are considered as non-interacting. The integrals, and the interaction sets again depend on the truncation routines, which are called inside integrate(). This approach allows to define interaction sets directly in the truncation routines. The code then automatically find the interacting elements bT without traversing all of them. This approach requires setting up the dual graph of the mesh which is stored mesh.neighbors. The resulting stiffness matrix is written to the disk in sparse format.
mesh | The mesh is of type MeshType and contains all information about the finite element descretization. See MeshStruct for more information. |
quadRule | Quadrature rules for inner, and outer elements as well as for the singular kernels. |
conf | General configuration, namely kernel, and forcing functions, as well as integration method. |
Definition at line 466 of file Cassemble.cpp.
map<string, void (*)(const double * x, const MeshType &mesh, double * forcing_out)> lookup_f |
Options for focing terms
Definition at line 23 of file Cassemble.cpp.
map<string, int (*)(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)> lookup_integrate |
Options for integration methods
Definition at line 63 of file Cassemble.cpp.
map<string, void (*)(const double * x, const ElementType &aT, const double * y, const ElementType &bT, const MeshType &mesh, double * kernel_val)> lookup_kernel |
Options for kernels
Definition at line 42 of file Cassemble.cpp.