G+Smo  23.12.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsGenericAssembler< T > Class Template Reference

Detailed Description

template<class T>
class gismo::gsGenericAssembler< T >

Assembles mass and stiffness matrices and right-hand sides on a given domain.

+ Inheritance diagram for gsGenericAssembler< T >:
+ Collaboration diagram for gsGenericAssembler< T >:

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. More...
 
const gsSparseMatrix< T > & assembleMass (const index_t patchIndex=-1, bool refresh=true)
 Mass assembly routine. More...
 
const gsSparseMatrix< T > & assembleMass2 ()
 Mass assembly routine. More...
 
const gsMatrix< T > & assembleMoments (const gsFunction< T > &func, index_t patchIndex=-1, bool refresh=true)
 Moments assembly routine. More...
 
const gsSparseMatrix< T > & assembleNeumann (const boundary_condition< T > &bc, bool refresh=true)
 Assemble Neumann boundary terms. More...
 
const gsSparseMatrix< T > & assembleNitsche (const boundary_condition< T > &bc, bool refresh=true)
 Assemble Nitsche terms for weakly imposing Dirichlet conditions. More...
 
const gsSparseMatrix< T > & assembleStiffness (const index_t patchIndex=-1, const bool refresh=true)
 Stiffness assembly routine. More...
 
bool check ()
 checks for consistency and legal values of the stored members.
 
virtual gsAssemblerclone () const
 Clone this Assembler, making a deep copy.
 
void computeDirichletDofs (short_t unk=0)
 Triggers computation of the Dirichlet dofs. More...
 
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. More...
 
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. More...
 
virtual gsAssemblercreate () 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) More...
 
gsSparseMatrix< T >::fullView fullMatrix ()
 Returns an expression of the "full" assembled sparse matrix. More...
 
const gsSparseMatrix< T >
::constFullView 
fullMatrix () const
 Returns an expression of the "full" assembled sparse matrix. More...
 
 gsGenericAssembler (const gsMultiPatch< T > &patches, const gsMultiBasis< T > &bases, const gsOptionList &opt=defaultOptions(), const gsBoundaryConditions< T > *bc=NULL)
 Constructor. More...
 
void homogenizeFixedDofs (short_t unk=0)
 Sets any Dirichlet values to homogeneous (if applicable) More...
 
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.
 
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 gsBasisRefs< T > &basis, const gsOptionList &opt=defaultOptions())
 Intitialize function for, sets data fields using the pde, a vector of bases and assembler options.
 
const gsSparseMatrix< T > & matrix () const
 Returns the left-hand global matrix.
 
const gsMultiBasis< T > & multiBasis (index_t k=0) const
 Return the multi-basis.
 
gsMultiBasis< T > & multiBasis (index_t k=0)
 Return the multi-basis. Note: if the basis is altered, then refresh() should be called.
 
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)
 
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 ElementVisitor >
void push (const ElementVisitor &visitor)
 Iterates over all elements of the domain and applies the ElementVisitor.
 
template<class BElementVisitor >
void push (const BElementVisitor &visitor, const boundary_condition< T > &BC)
 Applies the BElementVisitor to the boundary condition BC.
 
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 More...
 
void setFixedDofVector (gsMatrix< T > vals, short_t unk=0)
 the user can manually set the dirichlet Dofs for a given patch and unknown. More...
 
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. More...
 

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. More...
 
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. More...
 
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. More...
 
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.
 

Protected 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.
 

Constructor & Destructor Documentation

gsGenericAssembler ( const gsMultiPatch< T > &  patches,
const gsMultiBasis< T > &  bases,
const gsOptionList opt = defaultOptions(),
const gsBoundaryConditions< T > *  bc = NULL 
)
inline

Constructor.

Parameters
patchesThe geometry as gsMultiPatch
basesThe bases as gsMultiBasis
optThe assembler options
bcThe boundary conditons (used for setting up gsDofMapper)

Member Function Documentation

void apply ( ElementVisitor &  visitor,
size_t  patchIndex = 0,
boxSide  side = boundary::none 
)
protectedinherited

Generic assembly routine for volume or boundary integrals.

Parameters
[in]visitorThe visitor for the boundary or volume integral
[in]patchIndexThe considered patch
[in]sideThe considered boundary side, only necessary for boundary integrals.
const gsSparseMatrix< T > & assembleDG ( const boundaryInterface iFace,
bool  refresh = true 
)

Assemble dG interface terms.

Parameters
iFaceThe interface for which assembling is done
refreshCalles member refresh See gsVisiorDg for possible options
const gsSparseMatrix< T > & assembleMass ( const index_t  patchIndex = -1,
bool  refresh = true 
)

Mass assembly routine.

Parameters
patchIndexIf non-negative, only assemble matrix for specified patch
refreshCalles member refresh
const gsSparseMatrix<T>& assembleMass2 ( )
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.

Parameters
funcRight-hand-side (source) function
patchIndexIf non-negative, only assemble matrix for specified patch
refreshCalles member refresh
const gsSparseMatrix< T > & assembleNeumann ( const boundary_condition< T > &  bc,
bool  refresh = true 
)

Assemble Neumann boundary terms.

Parameters
bcThe boundary condition for which assembling is done
refreshCalles member refresh
const gsSparseMatrix< T > & assembleNitsche ( const boundary_condition< T > &  bc,
bool  refresh = true 
)

Assemble Nitsche terms for weakly imposing Dirichlet conditions.

Parameters
bcThe boundary condition for which assembling is done
refreshCalles member refresh

See gsVisiorNitsche for possible options

const gsSparseMatrix< T > & assembleStiffness ( const index_t  patchIndex = -1,
const bool  refresh = true 
)

Stiffness assembly routine.

Parameters
patchIndexIf non-negative, only assemble matrix for specified patch
refreshCalles member refresh
void computeDirichletDofs ( short_t  unk = 0)
inherited

Triggers computation of the Dirichlet dofs.

Parameters
[in]unkthe considered unknown
void computeDirichletDofsIntpl ( const gsDofMapper mapper,
const gsMultiBasis< T > &  mbasis,
const short_t  unk_ = 0 
)
protectedinherited

calculates the values of the eliminated dofs based on Interpolation.

Parameters
[in]mapperthe dofMapper for the considered unknown
[in]mbasisthe multipabasis for the considered unknown
[in]unk_the considered unknown
void computeDirichletDofsL2Proj ( const gsDofMapper mapper,
const gsMultiBasis< T > &  mbasis,
const short_t  unk_ = 0 
)
protectedinherited

calculates the values of the eliminated dofs based on L2 Projection.

Parameters
[in]mapperthe dofMapper for the considered unknown
[in]mbasisthe multipabasis for the considered unknown
[in]unk_the considered unknown
void constructSolution ( const gsMatrix< T > &  solVector,
gsMultiPatch< T > &  result,
short_t  unk = 0 
) const
virtualinherited

Construct solution from computed solution vector for a single unknows.

Parameters
[in]solVectorthe solution vector obtained from the linear system
[out]resultthe solution in form of a gsMultiBasis
[in]unkthe considered unknown
void constructSolution ( const gsMatrix< T > &  solVector,
gsMultiPatch< T > &  result,
const gsVector< index_t > &  unknowns 
) const
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.

unknowns. This method assumes that the specified unknowns have the
same basis.
Parameters
[in]solVectorthe solution vector obtained from the linear system
[out]resultthe solution seen as vectorfield in form of a gsMultiBasis
[in]unknownsthe considered vector of unknowns
const gsMatrix<T>& fixedDofs ( short_t  unk = 0) const
inlineinherited

Returns the Dirichlet values for a unknown (if applicable)

Parameters
[in]unkthe considered unknown
gsSparseMatrix<T>::fullView fullMatrix ( )
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 )

const gsSparseMatrix<T>::constFullView fullMatrix ( ) const
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 )

void homogenizeFixedDofs ( short_t  unk = 0)
inlineinherited

Sets any Dirichlet values to homogeneous (if applicable)

Parameters
[in]unkthe considered unknown
void penalizeDirichletDofs ( short_t  unk = 0)
inherited

Enforce Dirichlet boundary conditions by diagonal penalization

Parameters
[in]unkthe considered unknown
void setFixedDofs ( const gsMatrix< T > &  coefMatrix,
short_t  unk = 0,
size_t  patch = 0 
)
inherited

the user can manually set the dirichlet Dofs for a given patch and unknown, based on the Basis coefficients

Parameters
[in]coefMatrixthe coefficients of the function
[in]unkthe consideren unknown
[in]patchthe patch index
void setFixedDofVector ( gsMatrix< T >  vals,
short_t  unk = 0 
)
inherited

the user can manually set the dirichlet Dofs for a given patch and unknown.

Parameters
[in]valsthe values of the eliminated dofs.
[in]unkthe considered unknown
void updateSolution ( const gsMatrix< T > &  solVector,
gsMultiPatch< T > &  result,
theta = (T)(1) 
) const
virtualinherited

Update solution by adding the computed solution vector to the current solution specified by.

result. This method assumes that all
unknowns have the same basis.
Parameters
[in]solVectorthe solution vector obtained from the linear system
[out]resultthe solution in form of a gsMultiBasis,
solVector is added to the
coefficients of result.
Parameters
[in]thetadamping factor for update, theta = 1 corresponds to a full step.

Member Data Documentation

std::vector< gsMultiBasis<T> > m_bases
protectedinherited

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)

std::vector<gsMatrix<T> > m_ddof
protectedinherited

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().

memory::shared_ptr<gsPde<T> > m_pde_ptr
protectedinherited

The PDE: contains multi-patch domain, boundary conditions and coeffcient functions