G+Smo
24.08.0
Geometry + Simulation Modules
|
This module is responsible for generating matrices that discretize a PDE.
Namespaces | |
gismo::expr | |
This namespace contains expressions used for FE computations. | |
Classes | |
class | gsAssembler< T > |
The assembler class provides generic routines for volume and boundary integrals that are used for for matrix and right-hand side generation. More... | |
class | gsCDRAssembler< T > |
Implementation of an (multiple righ-hand side) Poisson solver. More... | |
class | gsCPPInterface< T > |
Provides a mapping between the corresponding sides of two patches sharing an interface, by means of a closest point projection. More... | |
class | gsGaussRule< T > |
Class that represents the (tensor) Gauss-Legendre quadrature rule. More... | |
class | gsGenericAssembler< T > |
Assembles mass and stiffness matrices and right-hand sides on a given domain. More... | |
class | gsHeatEquation< T > |
Constructs the assembler for the discretized isogeometric heat equation. More... | |
class | gsLobattoRule< T > |
Class that represents the (tensor) Gauss-Lobatto quadrature rule. More... | |
class | gsNewtonCotesRule< T > |
Class that represents the (tensor) Newton-Cotes quadrature rule. More... | |
class | gsOverIntegrateRule< T > |
Class that defines a mixed quadrature rule with different rules for the interior and the boundaries. More... | |
class | gsPatchRule< T > |
Class that represents the (tensor) patch quadrature rule. More... | |
class | gsPoissonAssembler< T > |
Implementation of an (multiple right-hand side) Poisson assembler. More... | |
class | gsQuadRule< T > |
Class representing a reference quadrature rule. More... | |
class | gsRemapInterface< T > |
Provides a mapping between the corresponding sides of two patches sharing an interface. More... | |
class | gsSparseSystem< T, symm > |
A sparse linear system indexed by sets of degrees of freedom. More... | |
class | gsVisitorBiharmonic< T > |
Visitor for the biharmonic equation. More... | |
class | gsVisitorCDR< T > |
Visitor for the convection-diffusion-reaction equation. More... | |
class | gsVisitorDg< T > |
Visitor for adding the interface conditions for the interior penalty methods of the Poisson problem. More... | |
class | gsVisitorGradGrad< T > |
The visitor computes element grad-grad integrals. More... | |
class | gsVisitorMass< T > |
The visitor computes element mass integrals. More... | |
class | gsVisitorMoments< T, paramCoef > |
Visitor for the moment vector of a function. More... | |
class | gsVisitorNeumann< T > |
Implementation of a Neumann BC for elliptic assemblers. More... | |
class | gsVisitorNeumannBiharmonic< T > |
Visitor for Neumann boundary condition for the biharmonic equation. More... | |
class | gsVisitorNitsche< T > |
Visitor for adding the terms associated to weak (Nitsche-type) imposition of the Dirichlet boundary conditions. More... | |
class | gsVisitorNitscheBiharmonic< T > |
Visitor for the weak imposition of the first-type dirichlet boundary condition. More... | |
class | gsVisitorPoisson< T, paramCoef > |
Visitor for the Poisson equation. More... | |
Functions | |
template<bool _coarsen, bool _admissible> | |
void | _markElements (const std::vector< T > &elError, const index_t refCriterion, const std::vector< gsHBoxCheck< 2, T > * > &predicates, HBoxContainer &elMarked) const |
Marks elements/cells for refinement. More... | |
template<class T > | |
void | gsMarkElementsForRef (const std::vector< T > &elError, int refCriterion, T refParameter, std::vector< bool > &elMarked) |
Marks elements/cells for refinement. More... | |
gsPoissonAssembler (gsMultiPatch< T > const &patches, gsMultiBasis< T > const &basis, gsBoundaryConditions< T > const &bconditions, const gsFunction< T > &rhs, dirichlet::strategy dirStrategy=dirichlet::elimination, iFace::strategy intStrategy=iFace::glue) | |
Constructor of the assembler object. More... | |
template<class T > | |
void | gsRefineMarkedElements (gsMultiBasis< T > &basis, const std::vector< bool > &elMarked, index_t refExtension=0) |
Refine a gsMultiBasis, based on a vector of element-markings. More... | |
template<class T > | |
void | gsUnrefineMarkedElements (gsMultiBasis< T > &basis, const std::vector< bool > &elMarked, index_t refExtension=0) |
Unrefine a gsMultiBasis, based on a vector of element-markings. More... | |
|
private |
Marks elements/cells for refinement.
Let the global error/error estimate \(\eta\) be a sum of element/cell-wise local contributions:
\[ \eta = \sum_{K} \eta_k \quad \mathrm{or} \quad \eta^2 = \sum_K \eta_K^2 \]
Computes a threshold \(\Theta\) and marks all elements \(K\) for refinement, for which
\[ \eta_K \geq \Theta \]
holds. Three criteria for computing \(\Theta\) are currently (26.Nov.2014) implemented:
Let \(\rho\) denote the m_input parameter parameter.
refCriterion = 1 = treshold, GARU-criterion (greatest appearing eRror utilization):
Threshold computed based on the largest of all appearing local errors:
\[ \Theta = \rho \cdot \max_K \{ \eta_K \} \]
The actual number of marked elements can vary in each refinement step, depending on the distribution of the error.
refCriterion = 2 = cellPercentage, PUCA-criterion (percentile-utilizing cutoff ascertainment):
In each step, a certain percentage of all elements are marked.
\[ \Theta = (1-\rho)\cdot 100\ \textrm{-percentile of}\ \{ \eta_K \}_K \]
For example, if \(\rho = 0.8\), those 20% of all elements which have the largest local errors are marked for refinement.
refCriterion = 3 = errorFraction, BULK-criterion ("Doerfler-marking"):
The threshold is chosen in such a manner that the local errors on the marked cells sum up to a certain fraction of the global error:
\[ \sum_{ K:\ \eta_K \geq \Theta } \eta_K \geq (1-\rho) \cdot \eta \]
elError | std::vector of local errors on some elements. | |
refCriterion | selects the criterion (see above) for marking elements. | |
parameter | parameter \( \rho \) for refinement criterion (see above). \(\rho = 0\) corresponds to global refinement, \( \rho=1\) corresponds to (almost) no refinement. | |
[out] | elMarked | std::vector of Booleans indicating whether the corresponding element is marked or not. |
void gismo::gsMarkElementsForRef | ( | const std::vector< T > & | elError, |
int | refCriterion, | ||
T | refParameter, | ||
std::vector< bool > & | elMarked | ||
) |
Marks elements/cells for refinement.
Let the global error/error estimate \(\eta\) be a sum of element/cell-wise local contributions:
\[ \eta = \sum_{K} \eta_k \quad \mathrm{or} \quad \eta^2 = \sum_K \eta_K^2 \]
Computes a threshold \(\Theta\) and marks all elements \(K\) for refinement, for which
\[ \eta_K \geq \Theta \]
holds. Three criteria for computing \(\Theta\) are currently (26.Nov.2014) implemented:
Let \(\rho\) denote the input parameter refParameter.
refCriterion = 1 = treshold, GARU-criterion (greatest appearing eRror utilization):
Threshold computed based on the largest of all appearing local errors:
\[ \Theta = \rho \cdot \max_K \{ \eta_K \} \]
The actual number of marked elements can vary in each refinement step, depending on the distribution of the error.
refCriterion = 2 = cellPercentage, PUCA-criterion (percentile-utilizing cutoff ascertainment):
In each step, a certain percentage of all elements are marked.
\[ \Theta = (1-\rho)\cdot 100\ \textrm{-percentile of}\ \{ \eta_K \}_K \]
For example, if \(\rho = 0.8\), those 20% of all elements which have the largest local errors are marked for refinement.
refCriterion = 3 = errorFraction, BULK-criterion ("Doerfler-marking"):
The threshold is chosen in such a manner that the local errors on the marked cells sum up to a certain fraction of the global error:
\[ \sum_{ K:\ \eta_K \geq \Theta } \eta_K \geq (1-\rho) \cdot \eta \]
elError | std::vector of local errors on some elements. | |
refCriterion | selects the criterion (see above) for marking elements. | |
refParameter | parameter \( \rho \) for refinement criterion (see above). \(\rho = 0\) corresponds to global refinement, \( \rho=1\) corresponds to (almost) no refinement. | |
[out] | elMarked | std::vector of Booleans indicating whether the corresponding element is marked or not. |
|
inline |
Constructor of the assembler object.
[in] | patches | is a gsMultiPatch object describing the geometry. |
[in] | basis | a multi-basis that contains patch-wise bases |
[in] | bconditions | is a gsBoundaryConditions object that holds all boundary conditions. |
[in] | rhs | is the right-hand side of the Poisson equation, \(\mathbf{f}\). |
[in] | dirStrategy | option for the treatment of Dirichlet boundary |
[in] | intStrategy | option for the treatment of patch interfaces |
void gismo::gsRefineMarkedElements | ( | gsMultiBasis< T > & | basis, |
const std::vector< bool > & | elMarked, | ||
index_t | refExtension = 0 |
||
) |
Refine a gsMultiBasis, based on a vector of element-markings.
Given the vector of element-markings (see gsRefineMarkedElements()), the corresponding element in the mesh underlying basis is refined.
It is possible to extend the refinement to the neighbouring elements of the marked area by setting the parameter refExtension.
This parameter is given as number of cells at the level of the marked element before refinement.
basis | gsMultiBasis to be refined adaptively. |
elMarked | std::vector of Booleans indicating for each element of the mesh underlying basis, whether it should be refined or not. |
refExtension | Specifies how large the refinement extension should be. Given as number of cells at the level before refinement. |
void gismo::gsUnrefineMarkedElements | ( | gsMultiBasis< T > & | basis, |
const std::vector< bool > & | elMarked, | ||
index_t | refExtension = 0 |
||
) |
Unrefine a gsMultiBasis, based on a vector of element-markings.
Given the vector of element-markings (see gsRefineMarkedElements()), the corresponding element in the mesh underlying basis is refined.
It is possible to extend the refinement to the neighbouring elements of the marked area by setting the parameter refExtension.
This parameter is given as number of cells at the level of the marked element before refinement.
basis | gsMultiBasis to be refined adaptively. |
elMarked | std::vector of Booleans indicating for each element of the mesh underlying basis, whether it should be refined or not. |
refExtension | Specifies how large the refinement extension should be. Given as number of cells at the level before refinement. |