G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsNsTimeIntegrator.hpp
Go to the documentation of this file.
1
16#pragma once
17
19
23
24namespace gismo
25{
26
27template <class T>
29 gsMassAssembler<T> & massAssembler_,
30 gsMultiPatch<T> * ALEvelocity,
31 gsBoundaryInterface * ALEpatches)
32 : stiffAssembler(stiffAssembler_),
33 massAssembler(massAssembler_),
34 velocityALE(ALEvelocity),
35 interface(ALEpatches)
36{
37 initialized = false;
39 m_ddof = stiffAssembler.allFixedDofs();
40 numIters = 0;
41 hasSavedState = false;
42}
43
44template <class T>
46{
47 gsOptionList opt = Base::defaultOptions();
48 opt.addInt("Scheme","Time integration scheme",time_integration::implicit_nonlinear);
49 opt.addReal("Theta","Time integration parametrer: 0 - explicit Euler, 1 - implicit Euler, 0.5 - Crank-Nicolson",0.5);
50 opt.addInt("Verbosity","Amount of information printed to the terminal: none, some, all",solver_verbosity::none);
51 opt.addReal("AbsTol","Absolute tolerance for the convergence cretiria",1e-10);
52 opt.addReal("RelTol","Relative tolerance for the stopping criteria",1e-7);
53 opt.addSwitch("ALE","ALE deformation is applied to the flow domain",false);
54 return opt;
55}
56
57template <class T>
59{
60 GISMO_ENSURE(solVector.rows() == stiffAssembler.numDofs(),"No initial conditions provided!");
61 stiffAssembler.assemble(solVector,m_ddof);
62 massAssembler.setFixedDofs(m_ddof);
63 massAssembler.assemble(true);
64 // IMEX stuff
65 oldSolVector = solVector;
66 oldTimeStep = 1.;
67
68 initialized = true;
69}
70
71template <class T>
73{
74 if (!initialized)
75 initialize();
76
77 tStep = timeStep;
78 if (m_options.getInt("Scheme") == time_integration::implicit_nonlinear)
79 implicitNonlinear();
80 if (m_options.getInt("Scheme") == time_integration::implicit_linear)
81 implicitLinear();
82}
83
84template <class T>
86{
87 T theta = m_options.getReal("Theta");
88 index_t numDofsVel = massAssembler.numDofs();
89 stiffAssembler.options().setInt("Assembly",ns_assembly::ossen);
90
91 // rhs = M*u_n - dt*(1-theta)*A(u_n)*u_n + dt*(1-theta)*F_n + dt*theta*F_n+1
92 // rhs: dt*(1-theta)*F_n
93 m_system.rhs() = tStep*(1-theta)*stiffAssembler.rhs();
94 // rhs: -dt*(1-theta)*A(u_n)*u_n
95 m_system.rhs().middleRows(0,numDofsVel) -= tStep*(1-theta)*stiffAssembler.matrix().block(0,0,numDofsVel,numDofsVel) *
96 solVector.middleRows(0,numDofsVel);
97 // rhs: M*u_n
98 m_system.rhs().middleRows(0,numDofsVel) += massAssembler.matrix() *
99 solVector.middleRows(0,numDofsVel);
100 // rhs: M_FD*u_DDOFS_n
101 m_system.rhs().middleRows(0,numDofsVel) -= massAssembler.rhs();
102
103 gsMultiPatch<T> velocity, pressure;
104 stiffAssembler.constructSolution(solVector + tStep/oldTimeStep*(solVector-oldSolVector),
105 stiffAssembler.allFixedDofs(),velocity,pressure);
106 if (m_options.getSwitch("ALE"))
107 for (size_t p = 0; p < interface->patches.size(); ++p)
108 velocity.patch(interface->patches[p].second).coefs() -=
109 velocityALE->patch(interface->patches[p].first).coefs();
110 stiffAssembler.assemble(velocity,pressure);
111
112 massAssembler.setFixedDofs(stiffAssembler.allFixedDofs());
113 if (m_options.getSwitch("ALE"))
114 massAssembler.assemble();
115 else
116 massAssembler.eliminateFixedDofs();
117
118 // rhs: dt*theta*F_n+1
119 m_system.rhs() += tStep*theta*stiffAssembler.rhs();
120 // rhs: -M_FD*u_DDOFS_n+1
121 m_system.rhs().middleRows(0,numDofsVel) += massAssembler.rhs();
122 // matrix = M + dt*theta*A(u_exp)
123 m_system.matrix() = tStep*stiffAssembler.matrix();
124 // we need to modify the (0,0,numDofsVel,numDofsVel) block of the matrix.
125 // unfortunately, eigen provides only a read-only block interfaces.
126 // the following is an ugly way to overcome this
127 gsSparseMatrix<T> tempVelocityBlock = m_system.matrix().block(0,0,numDofsVel,numDofsVel);
128 tempVelocityBlock *= (theta-1);
129 tempVelocityBlock += massAssembler.matrix();
130 tempVelocityBlock.conservativeResize(stiffAssembler.numDofs(),numDofsVel);
131 m_system.matrix().leftCols(numDofsVel) += tempVelocityBlock;
132 m_system.matrix().makeCompressed();
133
134 oldSolVector = solVector;
135 oldTimeStep = tStep;
136 m_ddof = stiffAssembler.allFixedDofs();
137 numIters = 1;
138
139#ifdef GISMO_WITH_PARDISO
140 gsSparseSolver<>::PardisoLU solver(m_system.matrix());
141 solVector = solver.solve(m_system.rhs());
142#else
143 gsSparseSolver<>::LU solver(m_system.matrix());
144 solVector = solver.solve(m_system.rhs());
145#endif
146}
147
148template <class T>
150{
151 stiffAssembler.options().setInt("Assembly",ns_assembly::newton_next);
152 T theta = m_options.getReal("Theta");
153 index_t numDofsVel = massAssembler.numDofs();
154
155 constRHS = tStep*(1-theta)*stiffAssembler.rhs();
156 constRHS.middleRows(0,numDofsVel).noalias() -= tStep*(1-theta)*stiffAssembler.matrix().block(0,0,numDofsVel,numDofsVel)*solVector.middleRows(0,numDofsVel);
157 constRHS.middleRows(0,numDofsVel).noalias() += massAssembler.matrix()*solVector.middleRows(0,numDofsVel);
158 constRHS.middleRows(0,numDofsVel).noalias() -= massAssembler.rhs();
159 massAssembler.setFixedDofs(stiffAssembler.allFixedDofs());
160 if (m_options.getSwitch("ALE"))
161 massAssembler.assemble();
162 else
163 massAssembler.eliminateFixedDofs();
164 constRHS.middleRows(0,numDofsVel).noalias() += massAssembler.rhs();
165
166 gsIterative<T> solver(*this,solVector,m_ddof);
167 solver.options().setInt("Verbosity",m_options.getInt("Verbosity"));
168 solver.options().setInt("Solver",linear_solver::LU);
169 solver.options().setInt("IterType",iteration_type::next);
170 solver.options().setReal("AbsTol",m_options.getReal("AbsTol"));
171 solver.options().setReal("RelTol",m_options.getReal("RelTol"));
172 solver.solve();
173
174 solVector = solver.solution();
175 m_ddof = stiffAssembler.allFixedDofs();
176 numIters = solver.numberIterations();
177}
178
179template <class T>
181 const std::vector<gsMatrix<T> > & fixedDoFs)
182{
183 T theta = m_options.getReal("Theta");
184
185 gsMultiPatch<T> velocity, pressure;
186 stiffAssembler.constructSolution(solutionVector,fixedDoFs,velocity,pressure);
187 if (m_options.getSwitch("ALE"))
188 for (size_t p = 0; p < interface->patches.size(); ++p)
189 velocity.patch(interface->patches[p].second).coefs() -=
190 velocityALE->patch(interface->patches[p].first).coefs();
191 stiffAssembler.assemble(velocity,pressure);
192
193 m_system.matrix() = tStep*stiffAssembler.matrix();
194 index_t numDofsVel = massAssembler.numDofs();
195 gsSparseMatrix<T> tempVelocityBlock = m_system.matrix().block(0,0,numDofsVel,numDofsVel);
196 tempVelocityBlock *= (theta-1);
197 tempVelocityBlock += massAssembler.matrix();
198 tempVelocityBlock.conservativeResize(stiffAssembler.numDofs(),numDofsVel);
199 m_system.matrix().leftCols(numDofsVel) += tempVelocityBlock;
200 m_system.matrix().makeCompressed();
201
202 m_system.rhs() = tStep*theta*stiffAssembler.rhs() + constRHS;
203 return true;
204}
205
206template <class T>
208{
209 stiffAssembler.constructSolution(solVector,m_ddof,velocity,pressure);
210}
211
212template <class T>
214
215template <class T>
216gsBaseAssembler<T> & gsNsTimeIntegrator<T>::assembler() { return stiffAssembler; }
217
218template <class T>
220{
221 if (!initialized)
222 initialize();
223
224 velVecSaved = solVector;
225 oldVecSaved = oldSolVector;
226 massRhsSaved = massAssembler.rhs();
227 stiffRhsSaved = stiffAssembler.rhs();
228 stiffMatrixSaved = stiffAssembler.matrix();
229 ddofsSaved = m_ddof;
230 hasSavedState = true;
231}
232
233template <class T>
235{
236 GISMO_ENSURE(hasSavedState,"No state saved!");
237 solVector = velVecSaved;
238 oldSolVector = oldVecSaved;
239 massAssembler.setRHS(massRhsSaved);
240 stiffAssembler.setMatrix(stiffMatrixSaved);
241 stiffAssembler.setRHS(stiffRhsSaved);
242 m_ddof = ddofsSaved;
243}
244
245} // namespace ends
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
A general iterative solver for nonlinear problems. An equation to solve is specified by an assembler ...
Definition gsIterative.h:43
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
gsGeometry< T > & patch(size_t i) const
Return the i-th patch.
Definition gsMultiPatch.h:292
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
void implicitLinear()
time integraton schemes
Definition gsNsTimeIntegrator.hpp:85
virtual void assemble()
assemble the linear system for the nonlinear solver
Definition gsBaseAssembler.h:40
gsNsTimeIntegrator(gsNsAssembler< T > &stiffAssembler_, gsMassAssembler< T > &massAssembler_, gsMultiPatch< T > *ALEvelocity=nullptr, gsBoundaryInterface *interfaceALE2Flow=nullptr)
Definition gsNsTimeIntegrator.hpp:28
gsNsAssembler< T > & stiffAssembler
assembler object that generates the static system
Definition gsNsTimeIntegrator.h:106
void makeTimeStep(T timeStep)
make a time step according to a chosen scheme
Definition gsNsTimeIntegrator.hpp:72
void recoverState()
recover solver state from the previously saved state
Definition gsNsTimeIntegrator.hpp:234
void saveState()
save solver state
Definition gsNsTimeIntegrator.hpp:219
gsBaseAssembler< T > & mAssembler()
assemblers' accessors
Definition gsNsTimeIntegrator.hpp:213
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsNsTimeIntegrator.hpp:45
bool initialized
initialization flag
Definition gsNsTimeIntegrator.h:110
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
void setInt(const std::string &label, const index_t &value)
Sets an existing option label to be equal to value.
Definition gsOptionList.cpp:158
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:211
void addSwitch(const std::string &label, const std::string &desc, const bool &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:235
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
#define index_t
Definition gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
A class providing an iterative solver for nonlinear problems.
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...
Provides time integration for incompressible Navier-Stokes equations.
The G+Smo namespace, containing all definitions for the library.
@ next
each iteration yields an update
Definition gsBaseUtils.h:94
@ newton_next
newton's method, 2nd order, yields updates to the solution
Definition gsBaseUtils.h:48
@ implicit_nonlinear
implicit scheme with linear problem (theta-scheme)
Definition gsBaseUtils.h:60
@ implicit_linear
explicit scheme with lumped mass matrix
Definition gsBaseUtils.h:59