G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsNsTimeIntegrator.h
Go to the documentation of this file.
1 
15 #pragma once
16 
19 
20 namespace gismo
21 {
22 
23 template <class T>
24 class gsNsAssembler;
25 template <class T>
26 class gsMassAssembler;
27 
30 template <class T>
32 {
33 public:
34  typedef gsBaseAssembler<T> Base;
37  gsNsTimeIntegrator(gsNsAssembler<T> & stiffAssembler_,
38  gsMassAssembler<T> & massAssembler_,
39  gsMultiPatch<T> * ALEvelocity = nullptr,
40  gsBoundaryInterface * interfaceALE2Flow = nullptr);
41 
44 
47  {
48  GISMO_ENSURE(solutionVector.rows() == stiffAssembler.numDofs(),"Wrong size of the solution vector: " + util::to_string(solutionVector.rows()) +
49  ". Must be: " + util::to_string(stiffAssembler.numDofs()));
51  initialized = false;
52  }
54  virtual void setFixedDofs(const std::vector<gsMatrix<T> > & ddofs)
55  {
56  Base::setFixedDofs(ddofs);
57  initialized = false;
58  }
59 
61  void makeTimeStep(T timeStep);
62 
64  virtual bool assemble(const gsMatrix<T> & solutionVector,
65  const std::vector<gsMatrix<T> > & fixedDoFs);
66 
68  virtual int numDofs() const { return stiffAssembler.numDofs(); }
69 
71  const gsMatrix<T> & solutionVector() const
72  {
73  GISMO_ENSURE(solVector.rows() == stiffAssembler.numDofs(),
74  "No initial conditions provided!");
75  return solVector;
76  }
77 
79  void saveState();
80 
82  void recoverState();
83 
85  index_t numberIterations() const { return numIters;}
86 
88  void constructSolution(gsMultiPatch<T> & velocity, gsMultiPatch<T> & pressure) const;
89 
92  gsBaseAssembler<T> & assembler();
93 
95  const gsBoundaryInterface & aleInterface() const {return *interface;}
96 
97 protected:
98  void initialize();
99 
101  void implicitLinear();
102  void implicitNonlinear();
103 
104 protected:
112  T tStep;
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
const gsMatrix< T > & solutionVector() const
returns solution vector
Definition: gsNsTimeIntegrator.h:71
void recoverState()
recover solver state from the previously saved state
Definition: gsNsTimeIntegrator.hpp:234
gsOptionList m_options
Options.
Definition: gsAssembler.h:285
std::string to_string(const C &value)
Converts value to string, assuming &quot;operator&lt;&lt;&quot; defined on C.
Definition: gsUtils.h:56
gsMassAssembler< T > & massAssembler
assembler object that generates the mass matrix
Definition: gsNsTimeIntegrator.h:108
void setSolutionVector(const gsMatrix< T > &solutionVector)
set intial conditions
Definition: gsNsTimeIntegrator.h:46
T oldTimeStep
IMEX stuff.
Definition: gsNsTimeIntegrator.h:120
T tStep
time step length
Definition: gsNsTimeIntegrator.h:112
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 saveState()
save solver state
Definition: gsNsTimeIntegrator.hpp:219
bool hasSavedState
saved state
Definition: gsNsTimeIntegrator.h:133
#define index_t
Definition: gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
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
gsBoundaryInterface * interface
mapping between the geometry patches and the ALE velocity patches
Definition: gsNsTimeIntegrator.h:130
gsNsAssembler< T > & stiffAssembler
assembler object that generates the static system
Definition: gsNsTimeIntegrator.h:106
void implicitLinear()
time integraton schemes
Definition: gsNsTimeIntegrator.hpp:85
Provides several simple utility and naming classes.
index_t numberIterations() const
number of iterations Newton&#39;s method required at the last time step; always 1 for IMEX ...
Definition: gsNsTimeIntegrator.h:85
virtual int numDofs() const
returns number of degrees of freedom
Definition: gsNsTimeIntegrator.h:68
Assembles the mass matrix and right-hand side vector for linear and nonlinear elasticity for 2D plain...
Definition: gsElTimeIntegrator.h:26
bool initialized
initialization flag
Definition: gsNsTimeIntegrator.h:110
gsMultiPatch< T > * velocityALE
ALE velocity.
Definition: gsNsTimeIntegrator.h:128
Base class for assemblers of gsElasticity.
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
Extends the gsAssembler class by adding functionality necessary for a general nonlinear solver...
Definition: gsALE.h:26
virtual void setFixedDofs(const std::vector< gsMatrix< T > > &ddofs)
set all fixed degrees of freedom
Definition: gsNsTimeIntegrator.h:54
gsSparseSystem< T > m_system
Global sparse linear system.
Definition: gsAssembler.h:290
std::vector< gsMatrix< T > > m_ddof
Definition: gsAssembler.h:295
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition: gsNsTimeIntegrator.hpp:45
void constructSolution(gsMultiPatch< T > &velocity, gsMultiPatch< T > &pressure) const
construct the solution using the stiffness matrix assembler
Definition: gsNsTimeIntegrator.hpp:207
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
gsMatrix< T > solVector
solution vector
Definition: gsNsTimeIntegrator.h:114
TODO: write.
Definition: gsNsAssembler.h:27
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition: gsAssembler.hpp:51
gsBaseAssembler< T > & mAssembler()
assemblers&#39; accessors
Definition: gsNsTimeIntegrator.hpp:213
Time integation for incompressible Navier-Stokes equations.
Definition: gsNsTimeIntegrator.h:31
gsMatrix< T > constRHS
Newton stuff.
Definition: gsNsTimeIntegrator.h:124
void makeTimeStep(T timeStep)
make a time step according to a chosen scheme
Definition: gsNsTimeIntegrator.hpp:72
gsNsTimeIntegrator(gsNsAssembler< T > &stiffAssembler_, gsMassAssembler< T > &massAssembler_, gsMultiPatch< T > *ALEvelocity=nullptr, gsBoundaryInterface *interfaceALE2Flow=nullptr)
Definition: gsNsTimeIntegrator.hpp:28