G+Smo
25.01.0
Geometry + Simulation Modules
|
Assembles mass and stiffness matrices and right-hand sides on a given domain.
Public Member Functions | |
const std::vector< gsMatrix< T > > & | allFixedDofs () const |
Returns all the Dirichlet values (if applicable) | |
virtual void | assemble () |
Main assemble routine, to be implemented in derived classes. | |
virtual void | assemble (const gsMultiPatch< T > &curSolution) |
Main non-linear assemble routine with input from current solution. | |
const gsSparseMatrix< T > & | assembleDG (const boundaryInterface &iFace, bool refresh=true) |
Assemble dG interface terms. | |
const gsSparseMatrix< T > & | assembleMass (const index_t patchIndex=-1, bool refresh=true) |
Mass assembly routine. | |
const gsSparseMatrix< T > & | assembleMass2 () |
Mass assembly routine. | |
const gsMatrix< T > & | assembleMoments (const gsFunction< T > &func, index_t patchIndex=-1, bool refresh=true) |
Moments assembly routine. | |
const gsSparseMatrix< T > & | assembleNeumann (const boundary_condition< T > &bc, bool refresh=true) |
Assemble Neumann boundary terms. | |
const gsSparseMatrix< T > & | assembleNitsche (const boundary_condition< T > &bc, bool refresh=true) |
Assemble Nitsche terms for weakly imposing Dirichlet conditions. | |
const gsSparseMatrix< T > & | assembleStiffness (const index_t patchIndex=-1, const bool refresh=true) |
Stiffness assembly routine. | |
bool | check () |
checks for consistency and legal values of the stored members. | |
virtual gsAssembler * | clone () const |
Clone this Assembler, making a deep copy. | |
void | computeDirichletDofs (short_t unk=0) |
Triggers computation of the Dirichlet dofs. | |
virtual void | constructSolution (const gsMatrix< T > &solVector, gsMultiPatch< T > &result, const gsVector< index_t > &unknowns) const |
Construct solution from computed solution vector for a set of unknowns. The result is a vectorfield, where each component is given the corresponding entry of. | |
virtual void | constructSolution (const gsMatrix< T > &solVector, gsMultiPatch< T > &result, short_t unk=0) const |
Construct solution from computed solution vector for a single unknows. | |
virtual gsAssembler * | create () const |
Create an empty Assembler of the derived type and return a pointer to it. Call the initialize functions to set the members. | |
void | finalize () |
finishes the assembling of the system matrix, i.e. calls its .makeCompressed() method. | |
const gsMatrix< T > & | fixedDofs (short_t unk=0) const |
Returns the Dirichlet values for a unknown (if applicable) | |
gsSparseMatrix< T >::fullView | fullMatrix () |
Returns an expression of the "full" assembled sparse matrix. | |
const gsSparseMatrix< T >::constFullView | fullMatrix () const |
Returns an expression of the "full" assembled sparse matrix. | |
gsGenericAssembler (const gsMultiPatch< T > &patches, const gsMultiBasis< T > &bases, const gsOptionList &opt=defaultOptions(), const gsBoundaryConditions< T > *bc=NULL) | |
Constructor. | |
void | homogenizeFixedDofs (short_t unk=0) |
Sets any Dirichlet values to homogeneous (if applicable) | |
void | initialize (const gsPde< T > &pde, const gsBasisRefs< T > &basis, const gsOptionList &opt=defaultOptions()) |
Intitialize function for, sets data fields using the pde, a vector of bases and assembler options. | |
void | initialize (const gsPde< T > &pde, const gsMultiBasis< T > &bases, const gsOptionList &opt=defaultOptions()) |
Intitialize function for, sets data fields using the pde, a multi-basis and assembler options. | |
void | initialize (const gsPde< T > &pde, const gsStdVectorRef< gsMultiBasis< T > > &bases, const gsOptionList &opt=defaultOptions()) |
Intitialize function for, sets data fields using the pde, a vector of multi-basis and assembler options. | |
void | initialize (typename gsPde< T >::Ptr pde, const gsStdVectorRef< gsMultiBasis< T > > &bases, const gsOptionList &opt=defaultOptions()) |
Intitialize function for, sets data fields using the pde, a vector of multi-basis and assembler options. | |
const gsSparseMatrix< T > & | matrix () const |
Returns the left-hand global matrix. | |
gsMultiBasis< T > & | multiBasis (index_t k=0) |
Return the multi-basis. Note: if the basis is altered, then refresh() should be called. | |
const gsMultiBasis< T > & | multiBasis (index_t k=0) const |
Return the multi-basis. | |
index_t | numColNz () const |
Provides an estimation of the number of non-zero matrix entries per column. This value can be used for sparse matrix memory allocation. | |
index_t | numDofs () const |
Returns the number of (free) degrees of freedom. | |
size_t | numMultiBasis () const |
Returns the number of multi-bases. | |
const gsMultiPatch< T > & | patches () const |
Return the multipatch. | |
const gsPde< T > & | pde () const |
Return the Pde. | |
void | penalizeDirichletDofs (short_t unk=0) |
T | penalty (index_t k) const |
Penalty constant for patch k, used for Nitsche and / Discontinuous Galerkin methods. | |
template<class ElementVisitor > | |
void | push () |
Iterates over all elements of the domain and applies the ElementVisitor. | |
template<class BElementVisitor > | |
void | push (const bcContainer &BCs) |
Iterates over all elements of the boundaries BCs and applies the BElementVisitor. | |
template<class BElementVisitor > | |
void | push (const BElementVisitor &visitor, const boundary_condition< T > &BC) |
Applies the BElementVisitor to the boundary condition BC. | |
template<class ElementVisitor > | |
void | push (const ElementVisitor &visitor) |
Iterates over all elements of the domain and applies the ElementVisitor. | |
template<class InterfaceVisitor > | |
void | pushInterface () |
Iterates over all elements of interfaces and applies the InterfaceVisitor. | |
void | refresh () |
Refreshes the sparse system (deletes assembled matrices and vectors) | |
void | refresh (gsDofMapper mapper) |
Refreshes the sparse system based on given dof mapper. | |
const gsMatrix< T > & | rhs () const |
Returns the left-hand side vector(s) ( multiple right hand sides possible ) | |
void | setFixedDofs (const gsMatrix< T > &coefMatrix, short_t unk=0, size_t patch=0) |
the user can manually set the dirichlet Dofs for a given patch and unknown, based on the Basis coefficients | |
void | setFixedDofVector (gsMatrix< T > vals, short_t unk=0) |
the user can manually set the dirichlet Dofs for a given patch and unknown. | |
void | setSparseSystem (gsSparseSystem< T > &sys) |
Swaps the actual sparse system with the given one. | |
const gsSparseSystem< T > & | system () const |
Returns the left-hand global matrix. | |
virtual void | updateSolution (const gsMatrix< T > &solVector, gsMultiPatch< T > &result, T theta=(T)(1)) const |
Update solution by adding the computed solution vector to the current solution specified by. | |
Static Public Member Functions | |
static gsOptionList | defaultOptions () |
Returns the list of default options for assembly. | |
Protected Member Functions | |
template<class ElementVisitor > | |
void | apply (ElementVisitor &visitor, size_t patchIndex=0, boxSide side=boundary::none) |
Generic assembly routine for volume or boundary integrals. | |
template<class InterfaceVisitor > | |
void | apply (InterfaceVisitor &visitor, const boundaryInterface &bi) |
Generic assembly routine for patch-interface integrals. | |
void | computeDirichletDofsIntpl (const gsDofMapper &mapper, const gsMultiBasis< T > &mbasis, const short_t unk_=0) |
calculates the values of the eliminated dofs based on Interpolation. | |
void | computeDirichletDofsL2Proj (const gsDofMapper &mapper, const gsMultiBasis< T > &mbasis, const short_t unk_=0) |
calculates the values of the eliminated dofs based on L2 Projection. | |
void | scalarProblemGalerkinRefresh () |
A prototype of the refresh function for a "standard" scalar problem. Creats one mapper based on the set options and initializes the sparse system (without allocating memory. | |
Private Attributes | |
std::vector< gsMultiBasis< T > > | m_bases |
std::vector< gsMatrix< T > > | m_ddof |
gsOptionList | m_options |
Options. | |
memory::shared_ptr< gsPde< T > > | m_pde_ptr |
gsSparseSystem< T > | m_system |
Global sparse linear system. | |
|
inline |
Constructor.
patches | The geometry as gsMultiPatch |
bases | The bases as gsMultiBasis |
opt | The assembler options |
bc | The boundary conditons (used for setting up gsDofMapper) |
|
protectedinherited |
Generic assembly routine for volume or boundary integrals.
[in] | visitor | The visitor for the boundary or volume integral |
[in] | patchIndex | The considered patch |
[in] | side | The considered boundary side, only necessary for boundary integrals. |
const gsSparseMatrix< T > & assembleDG | ( | const boundaryInterface & | iFace, |
bool | refresh = true |
||
) |
Assemble dG interface terms.
iFace | The interface for which assembling is done |
refresh | Calles member refresh See gsVisiorDg for possible options |
const gsSparseMatrix< T > & assembleMass | ( | const index_t | patchIndex = -1 , |
bool | refresh = true |
||
) |
Mass assembly routine.
patchIndex | If non-negative, only assemble matrix for specified patch |
refresh | Calles member refresh |
|
inline |
Mass assembly routine.
This routine uses the expression assembler and always overwrites already assembled matrices.
const gsMatrix< T > & assembleMoments | ( | const gsFunction< T > & | func, |
index_t | patchIndex = -1 , |
||
bool | refresh = true |
||
) |
Moments assembly routine.
func | Right-hand-side (source) function |
patchIndex | If non-negative, only assemble matrix for specified patch |
refresh | Calles member refresh |
const gsSparseMatrix< T > & assembleNeumann | ( | const boundary_condition< T > & | bc, |
bool | refresh = true |
||
) |
Assemble Neumann boundary terms.
bc | The boundary condition for which assembling is done |
refresh | Calles member refresh |
const gsSparseMatrix< T > & assembleNitsche | ( | const boundary_condition< T > & | bc, |
bool | refresh = true |
||
) |
Assemble Nitsche terms for weakly imposing Dirichlet conditions.
bc | The boundary condition for which assembling is done |
refresh | Calles member refresh |
See gsVisiorNitsche for possible options
const gsSparseMatrix< T > & assembleStiffness | ( | const index_t | patchIndex = -1 , |
const bool | refresh = true |
||
) |
Stiffness assembly routine.
patchIndex | If non-negative, only assemble matrix for specified patch |
refresh | Calles member refresh |
|
inherited |
Triggers computation of the Dirichlet dofs.
[in] | unk | the considered unknown |
|
protectedinherited |
calculates the values of the eliminated dofs based on Interpolation.
[in] | mapper | the dofMapper for the considered unknown |
[in] | mbasis | the multipabasis for the considered unknown |
[in] | unk_ | the considered unknown |
|
protectedinherited |
calculates the values of the eliminated dofs based on L2 Projection.
[in] | mapper | the dofMapper for the considered unknown |
[in] | mbasis | the multipabasis for the considered unknown |
[in] | unk_ | the considered unknown |
|
virtualinherited |
Construct solution from computed solution vector for a set of unknowns. The result is a vectorfield, where each component is given the corresponding entry of.
[in] | solVector | the solution vector obtained from the linear system |
[out] | result | the solution seen as vectorfield in form of a gsMultiBasis |
[in] | unknowns | the considered vector of unknowns |
Reimplemented in gsBaseAssembler< T >, gsElasticityAssembler< T >, gsElTimeIntegrator< T >, gsNsAssembler< T >, gsNsTimeIntegrator< T >, and gsBiharmonicAssembler< T, bhVisitor >.
|
virtualinherited |
Construct solution from computed solution vector for a single unknows.
[in] | solVector | the solution vector obtained from the linear system |
[out] | result | the solution in form of a gsMultiBasis |
[in] | unk | the considered unknown |
Reimplemented in gsBaseAssembler< T >, gsElasticityAssembler< T >, gsElTimeIntegrator< T >, gsNsAssembler< T >, gsNsTimeIntegrator< T >, and gsBiharmonicAssembler< T, bhVisitor >.
Returns the Dirichlet values for a unknown (if applicable)
[in] | unk | the considered unknown |
|
inline |
Returns an expression of the "full" assembled sparse matrix.
Note that matrix() might return a lower diagonal matrix, if we exploit possible symmetry during assembly (check: m_matrix.symmetry() == true )
|
inline |
Returns an expression of the "full" assembled sparse matrix.
Note that matrix() might return a lower diagonal matrix, if we exploit possible symmetry during assembly (check: m_matrix.symmetry() == true )
|
inlineinherited |
Sets any Dirichlet values to homogeneous (if applicable)
[in] | unk | the considered unknown |
|
inherited |
Enforce Dirichlet boundary conditions by diagonal penalization
[in] | unk | the considered unknown |
|
inherited |
the user can manually set the dirichlet Dofs for a given patch and unknown, based on the Basis coefficients
[in] | coefMatrix | the coefficients of the function |
[in] | unk | the consideren unknown |
[in] | patch | the patch index |
the user can manually set the dirichlet Dofs for a given patch and unknown.
[in] | vals | the values of the eliminated dofs. |
[in] | unk | the considered unknown |
|
virtualinherited |
Update solution by adding the computed solution vector to the current solution specified by.
[in] | solVector | the solution vector obtained from the linear system |
[out] | result | the solution in form of a gsMultiBasis, |
[in] | theta | damping factor for update, theta = 1 corresponds to a full step. |
|
private |
The discretization bases corresponding to the patches and to the number of solution fields that are to be computed. One entry of the vector corresponds to one or more unknown (in case of a system of PDEs)
|
private |
Fixed DoF values (if applicable, for instance eliminated Dirichlet DoFs). One entry of the vector corresponds to exactly one unknown, i.e. the size of m_ddof must fit m_system.colBlocks().
|
private |
The PDE: contains multi-patch domain, boundary conditions and coeffcient functions