G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsNsTimeIntegrator.h
Go to the documentation of this file.
1
15#pragma once
16
21
22namespace gismo
23{
24
27template <class T>
29{
30public:
34 gsNsTimeIntegrator(gsNsAssembler<T> & stiffAssembler_,
35 gsMassAssembler<T> & massAssembler_,
36 gsMultiPatch<T> * ALEvelocity = nullptr,
37 gsBoundaryInterface * interfaceALE2Flow = nullptr);
38
41
44 {
45 GISMO_ENSURE(solutionVector.rows() == stiffAssembler.numDofs(),"Wrong size of the solution vector: " + util::to_string(solutionVector.rows()) +
46 ". Must be: " + util::to_string(stiffAssembler.numDofs()));
48 initialized = false;
49 }
52 virtual void setFixedDofs(const std::vector<gsMatrix<T> > & ddofs)
53 {
54 Base::setFixedDofs(ddofs);
55 initialized = false;
56 }
57
59 void makeTimeStep(T timeStep);
60
62 using Base::assemble;
63 virtual bool assemble(const gsMatrix<T> & solutionVector,
64 const std::vector<gsMatrix<T> > & fixedDoFs);
65
67 virtual int numDofs() const { return stiffAssembler.numDofs(); }
68
71 {
72 GISMO_ENSURE(solVector.rows() == stiffAssembler.numDofs(),
73 "No initial conditions provided!");
74 return solVector;
75 }
76
78 void saveState();
79
81 void recoverState();
82
84 index_t numberIterations() const { return numIters;}
85
88 void constructSolution(gsMultiPatch<T> & velocity, gsMultiPatch<T> & pressure) const;
89
92 gsBaseAssembler<T> & assembler();
93
95 const gsBoundaryInterface & aleInterface() const {return *interface;}
96
97protected:
98 void initialize();
99
101 void implicitLinear();
102 void implicitNonlinear();
103
104protected:
115 using Base::m_system;
116 using Base::m_options;
117 using Base::m_ddof;
118
121 gsMatrix<T> oldSolVector;
122
125 index_t numIters;
126
130 gsBoundaryInterface * interface;
131
134 gsMatrix<T> velVecSaved;
135 gsMatrix<T> oldVecSaved;
136 gsMatrix<T> massRhsSaved;
137 gsMatrix<T> stiffRhsSaved;
138 gsSparseMatrix<T> stiffMatrixSaved;
139 std::vector<gsMatrix<T> > ddofsSaved;
140};
141
142}
143
144#ifndef GISMO_BUILD_LIB
145#include GISMO_HPP_HEADER(gsNsTimeIntegrator.hpp)
146#endif
gsSparseSystem< T > m_system
Global sparse linear system.
Definition gsAssembler.h:290
gsOptionList m_options
Options.
Definition gsAssembler.h:285
std::vector< gsMatrix< T > > m_ddof
Definition gsAssembler.h:295
Extends the gsAssembler class by adding functionality necessary for a general nonlinear solver....
Definition gsBaseAssembler.h:27
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition gsBaseAssembler.h:40
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.
Definition gsBaseAssembler.hpp:61
virtual void constructSolution(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDDofs, gsMultiPatch< T > &result, const gsVector< index_t > &unknowns) const
Constructs solution as a gsMultiPatch object from the solution vector and fixed DoFs.
Definition gsBaseAssembler.hpp:23
Assembles the mass matrix and right-hand side vector for linear and nonlinear elasticity for 2D plain...
Definition gsMassAssembler.h:31
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
TODO: write.
Definition gsNsAssembler.h:28
Time integation for incompressible Navier-Stokes equations.
Definition gsNsTimeIntegrator.h:29
bool hasSavedState
saved state
Definition gsNsTimeIntegrator.h:133
gsMultiPatch< T > * velocityALE
ALE velocity.
Definition gsNsTimeIntegrator.h:128
void implicitLinear()
time integraton schemes
Definition gsNsTimeIntegrator.hpp:85
gsMassAssembler< T > & massAssembler
assembler object that generates the mass matrix
Definition gsNsTimeIntegrator.h:108
virtual void assemble()
assemble the linear system for the nonlinear solver
Definition gsBaseAssembler.h:40
index_t numberIterations() const
number of iterations Newton's method required at the last time step; always 1 for IMEX
Definition gsNsTimeIntegrator.h:84
void setSolutionVector(const gsMatrix< T > &solutionVector)
set intial conditions
Definition gsNsTimeIntegrator.h:43
const gsMatrix< T > & solutionVector() const
returns solution vector
Definition gsNsTimeIntegrator.h:70
gsNsAssembler< T > & stiffAssembler
assembler object that generates the static system
Definition gsNsTimeIntegrator.h:106
virtual void setFixedDofs(const std::vector< gsMatrix< T > > &ddofs)
set all fixed degrees of freedom
Definition gsNsTimeIntegrator.h:52
void makeTimeStep(T timeStep)
make a time step according to a chosen scheme
Definition gsNsTimeIntegrator.hpp:72
const gsBoundaryInterface & aleInterface() const
get mapping between the flow domain patches and the ALE mapping patches (if only some patches of the ...
Definition gsNsTimeIntegrator.h:95
void recoverState()
recover solver state from the previously saved state
Definition gsNsTimeIntegrator.hpp:234
T tStep
time step length
Definition gsNsTimeIntegrator.h:112
gsMatrix< T > solVector
solution vector
Definition gsNsTimeIntegrator.h:114
void saveState()
save solver state
Definition gsNsTimeIntegrator.hpp:219
gsBaseAssembler< T > & mAssembler()
assemblers' accessors
Definition gsNsTimeIntegrator.hpp:213
T oldTimeStep
IMEX stuff.
Definition gsNsTimeIntegrator.h:120
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsNsTimeIntegrator.hpp:45
bool initialized
initialization flag
Definition gsNsTimeIntegrator.h:110
virtual int numDofs() const
returns number of degrees of freedom
Definition gsNsTimeIntegrator.h:67
gsMatrix< T > constRHS
Newton stuff.
Definition gsNsTimeIntegrator.h:124
gsBoundaryInterface * interface
mapping between the geometry patches and the ALE velocity patches
Definition gsNsTimeIntegrator.h:130
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
std::string to_string(const C &value)
Converts value to string, assuming "operator<<" defined on C.
Definition gsUtils.h:56
Base class for assemblers of gsElasticity.
Provides several simple utility and naming classes.
#define index_t
Definition gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
Provides mass matrix for elasticity systems in 2D plain strain and 3D continua.
Provides matrix and rhs assebmly for stationary and transient incompressible Stokes and Navier-Stokes...
The G+Smo namespace, containing all definitions for the library.