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

Detailed Description

template<class T, int MatOrder>
class gismo::gsINSAssembler< T, MatOrder >

A base class for incompressible Navier-Stokes assemblers.

Template Parameters
Treal number type
MatOrdersparse matrix storage order (ColMajor/RowMajor)
+ Inheritance diagram for gsINSAssembler< T, MatOrder >:
+ Collaboration diagram for gsINSAssembler< T, MatOrder >:

Public Member Functions

computeFlowRate (index_t patch, boxSide side, gsMatrix< T > solution) const
 Compute flow rate through a side of a given patch.
 
virtual gsField< T > constructSolution (const gsMatrix< T > &solVector, index_t unk) const
 Construct solution from computed solution vector for unknown unk.
 
virtual void fillStokesSystem (gsSparseMatrix< T, MatOrder > &stokesMat, gsMatrix< T > &stokesRhs)
 Fill the matrix and right-hand side for the Stokes problem.
 
gsAssemblerOptions getAssemblerOptions () const
 Returns the assembler options.
 
std::vector< gsMultiBasis< T > > & getBases ()
 Returns a reference to the discretization bases.
 
const gsBoundaryConditions< T > & getBCs () const
 Returns a const reference to the boundary conditions.
 
gsSparseMatrix< T, MatOrder > getBlockPU () const
 Returns the pressure-velocity block of the linear system.
 
virtual gsSparseMatrix< T, MatOrder > getBlockPUcomp (index_t i) const
 Returns part of pressure-velocity block for i-th velocity component.
 
const gsSparseMatrix< T, MatOrder > & getBlockUP () const
 Returns the velocity-pressure block of the linear system.
 
virtual gsSparseMatrix< T, MatOrder > getBlockUPcomp (index_t i) const
 Returns the part of velocity-pressure block for i-th velocity component.
 
virtual gsSparseMatrix< T, MatOrder > getBlockUU (bool linPartOnly=false)
 Returns the velocity-velocity block of the linear system.
 
virtual gsSparseMatrix< T, MatOrder > getBlockUUcompDiag (index_t i=0)
 Returns the diagonal block of velocity-velocity block for i-th component.
 
const std::vector< gsMatrix< T > > & getDirichletDofs () const
 Returns a const reference to the vectors of coefficients at the Dirichlet boundaries.
 
const std::vector< gsDofMapper > & getMappers () const
 Returns a const reference to the DOF mappers.
 
virtual gsSparseMatrix< T, MatOrder > & getMassMatrix (index_t unkID)
 Returns the mass matrix for unknown with index unk. There is also a const version.
 
const gsMultiPatch< T > & getPatches () const
 Returns a const reference to the multipatch representing the computational domain.
 
index_t getPdofs () const
 Returns the number of pressure DOFs.
 
index_t getPshift () const
 Returns the DOF shift of pressure (i.e. the total number of velocity DOFs).
 
const gsFunction< T > * getRhsFcn () const
 Returns a pointer to the right-hand-side function.
 
gsMatrix< T > getRhsP () const
 Returns the pressure part of the right-hand side.
 
virtual gsMatrix< T > getRhsU () const
 ///
 
virtual gsMatrix< T > getRhsUcomp (index_t i) const
 ///
 
const gsMatrix< T > & getSolution () const
 Returns a const reference to the current computed solution.
 
short_t getTarDim () const
 Returns the target dimension.
 
index_t getUdofs () const
 Returns the number of velocity DOFs (one velocity component).
 
getViscosity () const
 Returns the viscosity value.
 
virtual void initialize ()
 Initialize the assembler.
 
bool isInitialized ()
 Returns true if the assembler has been initialized.
 
virtual void markDofsAsEliminatedZeros (const std::vector< gsMatrix< index_t > > &boundaryDofs, const index_t unk)
 Eliminate given DOFs as homogeneous Dirichlet boundary.
 
virtual const gsSparseMatrix< T, MatOrder > & matrix () const
 Returns the assembled matrix.
 
virtual const gsSparseMatrix< T, MatOrder > & matrix (index_t unk) const
 Returns the assembled matrix for unknown with index unk (e.g. from two-equation turbulence models).
 
index_t numDofs () const
 Returns the number of degrees of freedom (DOFs).
 
index_t numDofsUnk (index_t i)
 Returns the number of DOFs for the i-th unknown.
 
gsOptionList options () const
 Returns the flow solver option list.
 
virtual const gsMatrix< T > & rhs () const
 Returns the assembled right-hand side.
 
virtual const gsMatrix< T > & rhs (index_t unk) const
 Returns the assembled right-hand side for unknown with index unk (e.g. from two-equation turbulence models).
 
virtual void update (const gsMatrix< T > &solVector, bool updateSol=true)
 

Protected Member Functions

void assembleBlock (gsFlowVisitor< T, MatOrder > &visitor, index_t testBasisID, gsSparseMatrix< T, MatOrder > &block, gsMatrix< T > &blockRhs)
 Assemble a matrix block.
 
virtual void assembleLinearPart ()
 Assemble the linear part of the problem.
 
virtual void assembleNonlinearPart ()
 Assemble the linear part of the problem.
 
void assembleRhs (gsFlowVisitor< T, MatOrder > &visitor, index_t testBasisID, gsMatrix< T > &rhs)
 Assemble the right-hand side.
 
void computeDirichletDofs (const index_t unk, const index_t basisID, gsMatrix< T > &ddofVector)
 Compute the coefficients of the basis functions at the Dirichlet boundaries.
 
void computeDirichletDofsIntpl (const index_t unk, const gsDofMapper &mapper, const gsMultiBasis< T > &mbasis, gsMatrix< T > &ddofVector)
 Compute the coefficients of the basis functions at the Dirichlet boundaries using interpolation.
 
void computeDirichletDofsL2Proj (const index_t unk, const gsDofMapper &mapper, const gsMultiBasis< T > &mbasis, gsMatrix< T > &ddofVector)
 Compute the coefficients of the basis functions at the Dirichlet boundaries using L2-projection.
 
virtual void fillBaseSystem ()
 Fill the linear part of the global matrix and right-hand side.
 
void fillGlobalMat_PP (gsSparseMatrix< T, MatOrder > &globalMat, const gsSparseMatrix< T, MatOrder > &sourceMat)
 Fill the pressure-pressure block into the global saddle-point matrix.
 
void fillGlobalMat_PU (gsSparseMatrix< T, MatOrder > &globalMat, const gsSparseMatrix< T, MatOrder > &sourceMat)
 Fill the pressure-velocity block into the global saddle-point matrix.
 
void fillGlobalMat_UP (gsSparseMatrix< T, MatOrder > &globalMat, const gsSparseMatrix< T, MatOrder > &sourceMat)
 Fill the velocity-pressure block into the global saddle-point matrix.
 
void fillGlobalMat_UU (gsSparseMatrix< T, MatOrder > &globalMat, const gsSparseMatrix< T, MatOrder > &sourceMat)
 Fill the velocity-velocity block into the global saddle-point matrix.
 
virtual void fillSystem ()
 Add the nonlinear part to the given matrix and right-hand side.
 
void initMembers ()
 Initialize the class members.
 
virtual void updateAssembly ()
 Assemble all that needs to be updated in each nonlinear iteration.
 
virtual void updateCurrentSolField (const gsMatrix< T > &solVector, bool updateSol)
 Update the current solution field stored in the assembler (used as convection coefficient).
 
virtual void updateDofMappers ()
 Update the DOF mappers in all visitors (when DOF numbers change, e.g. after markDofsAsEliminatedZeros()).
 
virtual void updateSizes ()
 Update sizes of members (when DOF numbers change, e.g. after markDofsAsEliminatedZeros()).
 

Member Function Documentation

◆ assembleBlock()

template<class T , int MatOrder>
void assembleBlock ( gsFlowVisitor< T, MatOrder > &  visitor,
index_t  testBasisID,
gsSparseMatrix< T, MatOrder > &  block,
gsMatrix< T > &  blockRhs 
)
protectedinherited

Assemble a matrix block.

Parameters
[in]visitorvisitor for the required block
[in]testBasisIDID of the test basis
[out]blockthe resulting matrix block
[out]blockRhsright-hand side for the matrix block (arising from eliminated Dirichlet DOFs)

◆ assembleRhs()

template<class T , int MatOrder>
void assembleRhs ( gsFlowVisitor< T, MatOrder > &  visitor,
index_t  testBasisID,
gsMatrix< T > &  rhs 
)
protectedinherited

Assemble the right-hand side.

Parameters
[in]visitorvisitor for the right-hand side
[in]testBasisIDID of the test basis
[out]rhsthe resulting right-hand side vector

◆ computeDirichletDofs()

template<class T , int MatOrder>
void computeDirichletDofs ( const index_t  unk,
const index_t  basisID,
gsMatrix< T > &  ddofVector 
)
protectedinherited

Compute the coefficients of the basis functions at the Dirichlet boundaries.

Parameters
[in]unkthe considered unknown (0 - velocity, 1 - pressure)
[in]basisIDthe index of the basis corresponding to unk (same as unk in this case)
[in]ddofVectorreference to the vector where computed coefficients will be stored

◆ computeDirichletDofsIntpl()

template<class T , int MatOrder>
void computeDirichletDofsIntpl ( const index_t  unk,
const gsDofMapper mapper,
const gsMultiBasis< T > &  mbasis,
gsMatrix< T > &  ddofVector 
)
protectedinherited

Compute the coefficients of the basis functions at the Dirichlet boundaries using interpolation.

Parameters
[in]unkthe considered unknown (0 - velocity, 1 - pressure)
[in]mapperreference to the DOF mapper for unk
[in]mbasisreference to the basis corresponding to unk
[out]ddofVectorreference to the vector where computed coefficients will be stored

◆ computeDirichletDofsL2Proj()

template<class T , int MatOrder>
void computeDirichletDofsL2Proj ( const index_t  unk,
const gsDofMapper mapper,
const gsMultiBasis< T > &  mbasis,
gsMatrix< T > &  ddofVector 
)
protectedinherited

Compute the coefficients of the basis functions at the Dirichlet boundaries using L2-projection.

Parameters
[in]unkthe considered unknown (0 - velocity, 1 - pressure)
[in]mapperreference to the DOF mapper for unk
[in]mbasisreference to the basis corresponding to unk
[out]ddofVectorreference to the vector where computed coefficients will be stored

◆ computeFlowRate()

template<class T , int MatOrder>
T computeFlowRate ( index_t  patch,
boxSide  side,
gsMatrix< T >  solution 
) const

Compute flow rate through a side of a given patch.

Parameters
[in]patchthe given patch ID
[in]sidethe given patch side
[in]solutionsolution vector to compute the flow rate from

◆ constructSolution()

template<class T , int MatOrder>
gsField< T > constructSolution ( const gsMatrix< T > &  solVector,
index_t  unk 
) const
virtual

Construct solution from computed solution vector for unknown unk.

Parameters
[in]solVectorthe solution vector obtained from the linear system
[out]resultthe resulting solution as a gsMultiPatch object
[in]unkthe considered unknown

Reimplemented from gsFlowAssemblerBase< T, MatOrder >.

◆ fillGlobalMat_PP()

template<class T , int MatOrder>
void fillGlobalMat_PP ( gsSparseMatrix< T, MatOrder > &  globalMat,
const gsSparseMatrix< T, MatOrder > &  sourceMat 
)
protected

Fill the pressure-pressure block into the global saddle-point matrix.

Parameters
globalMat[out]global saddle-point matrix
sourceMat[in]pressure-pressure block

◆ fillGlobalMat_PU()

template<class T , int MatOrder>
void fillGlobalMat_PU ( gsSparseMatrix< T, MatOrder > &  globalMat,
const gsSparseMatrix< T, MatOrder > &  sourceMat 
)
protected

Fill the pressure-velocity block into the global saddle-point matrix.

Parameters
globalMat[out]global saddle-point matrix
sourceMat[in]pressure-velocitye block

◆ fillGlobalMat_UP()

template<class T , int MatOrder>
void fillGlobalMat_UP ( gsSparseMatrix< T, MatOrder > &  globalMat,
const gsSparseMatrix< T, MatOrder > &  sourceMat 
)
protected

Fill the velocity-pressure block into the global saddle-point matrix.

Parameters
globalMat[out]global saddle-point matrix
sourceMat[in]velocity-pressure block

◆ fillGlobalMat_UU()

template<class T , int MatOrder>
void fillGlobalMat_UU ( gsSparseMatrix< T, MatOrder > &  globalMat,
const gsSparseMatrix< T, MatOrder > &  sourceMat 
)
protected

Fill the velocity-velocity block into the global saddle-point matrix.

Parameters
globalMat[out]global saddle-point matrix
sourceMat[in]velocity-velocity block (either for one velocity component or the whole block for all components)

◆ getBases()

template<class T , int MatOrder>
std::vector< gsMultiBasis< T > > & getBases ( )
inline

Returns a reference to the discretization bases.

In the case of velocity and pressure, the velocity basis is stored first, the pressure basis is second.

There is also a const version returning a const reference.

◆ getBlockUU()

template<class T , int MatOrder>
gsSparseMatrix< T, MatOrder > getBlockUU ( bool  linPartOnly = false)
virtual

Returns the velocity-velocity block of the linear system.

Parameters
[in]linPartOnlyif true, returns only the linear part of the velocity-velocity block

◆ getDirichletDofs()

template<class T , int MatOrder>
const std::vector< gsMatrix< T > > & getDirichletDofs ( ) const
inlineinherited

Returns a const reference to the vectors of coefficients at the Dirichlet boundaries.

In the case of velocity and pressure, the vector of velocity coefficients is stored first, the vector of pressure coefficients is second.

◆ getMappers()

template<class T , int MatOrder>
const std::vector< gsDofMapper > & getMappers ( ) const
inlineinherited

Returns a const reference to the DOF mappers.

In the case of velocity and pressure, the mapper for velocity is stored first, the mapper for pressure is second.

◆ getMassMatrix()

template<class T , int MatOrder>
virtual gsSparseMatrix< T, MatOrder > & getMassMatrix ( index_t  unkID)
inlinevirtual

Returns the mass matrix for unknown with index unk. There is also a const version.

Parameters
[in]unkIDindex of the unknown (0 - velocity, 1 - pressure)

Reimplemented from gsFlowAssemblerBase< T, MatOrder >.

◆ getRhsU()

template<class T , int MatOrder>
gsMatrix< T > getRhsU ( ) const
virtual

///

Returns the velocity part of the right-hand side.

Reimplemented in gsINSAssemblerUnsteady< T, MatOrder >.

◆ getRhsUcomp()

template<class T , int MatOrder>
virtual gsMatrix< T > getRhsUcomp ( index_t  i) const
inlinevirtual

///

Returns part of the right-hand side for i-th velocity component.

◆ markDofsAsEliminatedZeros()

template<class T , int MatOrder>
void markDofsAsEliminatedZeros ( const std::vector< gsMatrix< index_t > > &  boundaryDofs,
const index_t  unk 
)
virtual

Eliminate given DOFs as homogeneous Dirichlet boundary.

Parameters
[in]boundaryDofsindices of the given boundary DOFs
[in]unkthe considered unknown

Reimplemented from gsFlowAssemblerBase< T, MatOrder >.

◆ matrix()

template<class T , int MatOrder>
virtual const gsSparseMatrix< T, MatOrder > & matrix ( index_t  unk) const
inlinevirtualinherited

Returns the assembled matrix for unknown with index unk (e.g. from two-equation turbulence models).

Parameters
[in]unkindex of the unknown

◆ numDofsUnk()

template<class T , int MatOrder>
index_t numDofsUnk ( index_t  i)

Returns the number of DOFs for the i-th unknown.

Parameters
[in]i0 - velocity, 1 - pressure

◆ rhs()

template<class T , int MatOrder>
virtual const gsMatrix< T > & rhs ( index_t  unk) const
inlinevirtualinherited

Returns the assembled right-hand side for unknown with index unk (e.g. from two-equation turbulence models).

Parameters
[in]unkindex of the unknown

◆ update()

template<class T , int MatOrder>
void update ( const gsMatrix< T > &  solVector,
bool  updateSol = true 
)
virtual
Parameters
solVector
[in]updateSoltrue - save solVector into m_solution (false is used in the inner Picard iteration for unsteady problem)

Reimplemented from gsFlowAssemblerBase< T, MatOrder >.

Reimplemented in gsINSAssemblerUnsteady< T, MatOrder >.

◆ updateCurrentSolField()

template<class T , int MatOrder>
void updateCurrentSolField ( const gsMatrix< T > &  solVector,
bool  updateSol 
)
protectedvirtual

Update the current solution field stored in the assembler (used as convection coefficient).

Parameters
[in]solVectornew solution vector
[in]updateSoltrue - save solVector into m_solution (false is used in the inner nonlinear iteration for unsteady problem)

Reimplemented from gsFlowAssemblerBase< T, MatOrder >.

Reimplemented in gsINSAssemblerUnsteady< T, MatOrder >.