G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsElTimeIntegrator< T > Class Template Reference

Detailed Description

template<class T>
class gismo::gsElTimeIntegrator< T >

Time integation for equations of dynamic elasticity with implicit schemes.

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

Public Member Functions

const std::vector< gsMatrix< T > > & allFixedDofs () const
 Returns all the Dirichlet values (if applicable)
 
virtual void assemble ()
 assemble the linear system for the nonlinear solver
 
virtual void assemble (bool)
 assemble the linear system for the nonlinear solver
 
virtual bool assemble (const gsMatrix< T > &solutionVector, const std::vector< gsMatrix< T > > &fixedDDoFs)=0
 assemble the linear system for the nonlinear solver
 
virtual bool assemble (const gsMatrix< T > &solutionVector, const std::vector< gsMatrix< T > > &fixedDoFs)
 
virtual void assemble (const gsMultiPatch< T > &)
 assemble the linear system for the nonlinear solver
 
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.
 
virtual void constructSolution (const gsMatrix< T > &, const std::vector< gsMatrix< T > > &, gsMultiPatch< T > &) const
 construct displacement and pressure (if applicable) using the stiffness assembler
 
virtual void constructSolution (const gsMatrix< T > &, gsMultiPatch< T > &, const gsVector< index_t > &) const
 construct displacement and pressure (if applicable) using the stiffness assembler
 
virtual void constructSolution (const gsMatrix< T > &, gsMultiPatch< T > &, short_t=0) const
 construct displacement and pressure (if applicable) using the stiffness assembler
 
virtual void constructSolution (const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDDofs, gsMultiPatch< T > &result, const gsVector< index_t > &unknowns) const
 construct displacement and pressure (if applicable) using the stiffness assembler
 
void constructSolution (gsMultiPatch< T > &displacement) const
 construct displacement using the stiffness assembler
 
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.
 
virtual void eliminateFixedDofs ()
 Eliminates new Dirichelt degrees of fredom.
 
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)
 
virtual void getFixedDofs (index_t patch, boxSide side, gsMatrix< T > &ddofs) const
 
 gsElTimeIntegrator (gsElasticityAssembler< T > &stiffAssembler_, gsMassAssembler< T > &massAssembler_)
 
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.
 
void makeTimeStep (T timeStep)
 make a time step according to a chosen scheme
 
gsBaseAssembler< T > & mAssembler ()
 assemblers' accessors
 
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 numberIterations () const
 number of iterations Newton's method required at the last time step
 
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.
 
virtual int numDofs () const
 return the number of free degrees of freedom
 
virtual index_t numFixedDofs () const
 get the size of the Dirichlet vector for elimination
 
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 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 recoverState ()
 recover solver state from saved state
 
virtual void refresh ()
 Creates the mappers and setups the sparse system. to be implemented in derived classes, see scalarProblemGalerkinRefresh() for a possible implementation.
 
const gsMatrix< T > & rhs () const
 Returns the left-hand side vector(s) ( multiple right hand sides possible )
 
void saveState ()
 save solver state
 
void setDisplacementVector (const gsMatrix< T > &displacementVector)
 set intial conditions
 
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
 
virtual void setFixedDofs (const std::vector< gsMatrix< T > > &ddofs)
 set all fixed degrees of freedom
 
virtual void setFixedDofs (index_t patch, boxSide side, const gsMatrix< T > &ddofs, bool oneUnk=false)
 Set Dirichet degrees of freedom on a given side of a given patch from a given matrix.
 
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 gsMatrix< T > & solutionVector () const
 returns complete solution vector (displacement + possibly pressure)
 
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.
 
const gsMatrix< T > & velocityVector () const
 returns vector of displacement DoFs
 

Static Public Member Functions

static gsOptionList defaultOptions ()
 Returns the list of default options for assembly.
 

Protected Member Functions

alpha1 ()
 time integration scheme coefficients
 
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.
 
gsMatrix< T > implicitLinear ()
 time integraton schemes
 
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

gsMatrix< T > accVector
 vector of acceleration DoFs
 
bool hasSavedState
 saved state
 
bool initialized
 initialization flag
 
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.
 
gsMassAssembler< T > & massAssembler
 assembler object that generates the mass matrix
 
gsMatrix< T > newSolVector
 temporary objects for memory efficiency
 
index_t numIters
 number of iterations Newton's method took to converge at the last time step
 
gsMatrix< T > solVector
 vector of displacement DoFs ( + possibly pressure)
 
gsElasticityAssembler< T > & stiffAssembler
 assembler object that generates the static system
 
tStep
 time step length
 
gsMatrix< T > velVector
 vector of velocity DoFs
 

Constructor & Destructor Documentation

◆ gsElTimeIntegrator()

template<class T >
gsElTimeIntegrator ( gsElasticityAssembler< T > &  stiffAssembler_,
gsMassAssembler< T > &  massAssembler_ 
)

constructor method. requires a gsElasticityAssembler for construction of the static linear system and a gsMassAssembler for the mass matrix

Member Function Documentation

◆ apply()

template<class T >
template<class ElementVisitor >
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.

◆ assemble()

template<class T >
bool assemble ( const gsMatrix< T > &  solutionVector,
const std::vector< gsMatrix< T > > &  fixedDDoFs 
)
virtual

Assembles the tangential linear system for Newton's method given the current solution in the form of free and fixed/Dirichelt degrees of freedom. Checks if the current solution is valid (Newton's solver can exit safely if invalid).

Implements gsBaseAssembler< T >.

◆ computeDirichletDofs()

template<class T >
void computeDirichletDofs ( short_t  unk = 0)
inherited

Triggers computation of the Dirichlet dofs.

Parameters
[in]unkthe considered unknown

◆ computeDirichletDofsIntpl()

template<class T >
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

◆ computeDirichletDofsL2Proj()

template<class T >
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

◆ fixedDofs()

template<class T >
const gsMatrix< T > & fixedDofs ( short_t  unk = 0) const
inlineinherited

Returns the Dirichlet values for a unknown (if applicable)

Parameters
[in]unkthe considered unknown

◆ getFixedDofs()

template<class T >
void getFixedDofs ( index_t  patch,
boxSide  side,
gsMatrix< T > &  ddofs 
) const
virtualinherited

get fixed degrees of freedom corresponding to a given part of the bdry. each column of the resulting matrix correspond to one variable/component of the vector-valued vairable

◆ homogenizeFixedDofs()

template<class T >
void homogenizeFixedDofs ( short_t  unk = 0)
inlineinherited

Sets any Dirichlet values to homogeneous (if applicable)

Parameters
[in]unkthe considered unknown

◆ penalizeDirichletDofs()

template<class T >
void penalizeDirichletDofs ( short_t  unk = 0)
inherited

Enforce Dirichlet boundary conditions by diagonal penalization

Parameters
[in]unkthe considered unknown

◆ setFixedDofs() [1/2]

template<class T >
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

◆ setFixedDofs() [2/2]

template<class T >
void setFixedDofs ( index_t  patch,
boxSide  side,
const gsMatrix< T > &  ddofs,
bool  oneUnk = false 
)
virtualinherited

Set Dirichet degrees of freedom on a given side of a given patch from a given matrix.

A usual source of degrees of freedom is another geometry where it is known that the corresponding bases match. The function is not safe in that it assumes a correct numbering of DoFs in a given matrix. To allocate space for these DoFs in the assembler, add an empty/zero Dirichlet boundary condition to gsBoundaryCondtions container that is passed to the assembler constructor.

Reimplemented in gsNsTimeIntegrator< T >.

◆ setFixedDofVector()

template<class T >
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

◆ updateSolution()

template<class T >
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.

◆ velocityVector()

template<class T >
const gsMatrix< T > & velocityVector ( ) const
inline

returns vector of displacement DoFs

returns vector of velocity DoFs

Member Data Documentation

◆ m_bases

template<class T >
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)

◆ m_ddof

template<class T >
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().

◆ m_pde_ptr

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

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