G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsElTimeIntegrator.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 : stiffAssembler(stiffAssembler_),
31 massAssembler(massAssembler_)
32{
33 initialized = false;
35 m_ddof = stiffAssembler.allFixedDofs();
36 numIters = 0;
37 hasSavedState = false;
41}
42
43template <class T>
45{
46 gsOptionList opt = Base::defaultOptions();
47 opt.addInt("Scheme","Time integration scheme",time_integration::implicit_linear);
48 opt.addReal("Beta","Parameter beta for the time integration scheme, see Wriggers, Nonlinear FEM, p.213 ",0.25);
49 opt.addReal("Gamma","Parameter gamma for the time integration scheme, see Wriggers, Nonlinear FEM, p.213 ",0.5);
50 opt.addInt("Verbosity","Amount of information printed to the terminal: none, some, all",solver_verbosity::none);
51 return opt;
52}
53
54template <class T>
56{
57 stiffAssembler.assemble(solVector,m_ddof);
58 massAssembler.assemble();
59 gsSparseSolver<>::SimplicialLDLT solver(massAssembler.matrix());
60 accVector = solver.solve(stiffAssembler.rhs().middleRows(0,massAssembler.numDofs()));
61 initialized = true;
62}
63
64template <class T>
66{
67 if (!initialized)
68 initialize();
69
70 tStep = timeStep;
71 if (m_options.getInt("Scheme") == time_integration::implicit_linear)
72 newSolVector = implicitLinear();
73 if (m_options.getInt("Scheme") == time_integration::implicit_nonlinear)
74 newSolVector = implicitNonlinear();
75 oldVelVector = velVector;
76 dispVectorDiff = (newSolVector - solVector).middleRows(0,massAssembler.numDofs());
77 velVector = alpha4()*dispVectorDiff + alpha5()*oldVelVector + alpha6()*accVector;
78 accVector = alpha1()*dispVectorDiff - alpha2()*oldVelVector - alpha3()*accVector;
79 solVector = newSolVector;
80}
81
82template <class T>
84{
85 if (massAssembler.numDofs() == stiffAssembler.numDofs())
86 { // displacement formulation
87 m_system.matrix() = alpha1()*massAssembler.matrix() + stiffAssembler.matrix();
88 m_system.matrix().makeCompressed();
89 m_system.rhs() = massAssembler.matrix()*(alpha1()*solVector.middleRows(0,massAssembler.numDofs())
90 + alpha2()*velVector + alpha3()*accVector) + stiffAssembler.rhs();
91 }
92 else
93 { // displacement-pressure formulation
94 m_system.matrix() = stiffAssembler.matrix();
95 tempMassBlock = alpha1()*massAssembler.matrix();
96 tempMassBlock.conservativeResize(stiffAssembler.numDofs(),massAssembler.numDofs());
97 m_system.matrix().leftCols(massAssembler.numDofs()) += tempMassBlock;
98 m_system.matrix().makeCompressed();
99 m_system.rhs() = stiffAssembler.rhs();
100 m_system.rhs().middleRows(0,massAssembler.numDofs()) +=
101 massAssembler.matrix()*(alpha1()*solVector.middleRows(0,massAssembler.numDofs())
102 + alpha2()*velVector + alpha3()*accVector);
103 }
104
105#ifdef GISMO_WITH_PARDISO
106 gsSparseSolver<>::PardisoLDLT solver(m_system.matrix());
107 return solver.solve(m_system.rhs());
108#else
109 gsSparseSolver<>::SimplicialLDLT solver(m_system.matrix());
110 return solver.solve(m_system.rhs());
111#endif
112 numIters = 1;
113 return gsMatrix<T>();
114}
115
116template <class T>
118{
119 gsIterative<T> solver(*this,solVector);
120 solver.options().setInt("Verbosity",m_options.getInt("Verbosity"));
121 solver.options().setInt("Solver",linear_solver::LDLT);
122 solver.solve();
123 numIters = solver.numberIterations();
124 return solver.solution();
125}
126
127template <class T>
128int gsElTimeIntegrator<T>::numDofs() const { return stiffAssembler.numDofs(); }
129
130template <class T>
132 const std::vector<gsMatrix<T> > & fixedDoFs)
133{
134 stiffAssembler.assemble(solutionVector,fixedDoFs);
135 if (massAssembler.numDofs() == stiffAssembler.numDofs())
136 { // displacement formulation
137 m_system.matrix() = alpha1()*massAssembler.matrix() + stiffAssembler.matrix();
138 m_system.matrix().makeCompressed();
139 m_system.rhs() = stiffAssembler.rhs() +
140 massAssembler.matrix()*(alpha1()*(solVector-solutionVector) + alpha2()*velVector + alpha3()*accVector);
141 }
142 else
143 { // displacement-pressure formulation
144 m_system.matrix() = stiffAssembler.matrix();
145 tempMassBlock = alpha1()*massAssembler.matrix();
146 tempMassBlock.conservativeResize(stiffAssembler.numDofs(),massAssembler.numDofs());
147 m_system.matrix().leftCols(massAssembler.numDofs()) += tempMassBlock;
148 m_system.matrix().makeCompressed();
149 m_system.rhs() = stiffAssembler.rhs();
150 m_system.rhs().middleRows(0,massAssembler.numDofs()) +=
151 massAssembler.matrix()*(alpha1()*(solVector-solutionVector).middleRows(0,massAssembler.numDofs())
152 + alpha2()*velVector + alpha3()*accVector);
153 }
154
155 return true;
156}
157
158template <class T>
160{
161 stiffAssembler.constructSolution(solVector,m_ddof,displacement);
162}
163
164template <class T>
166{
167 GISMO_ENSURE(stiffAssembler.numDofs() > massAssembler.numDofs(),
168 "This is a displacement-only formulation. Can't construct pressure");
169 stiffAssembler.constructSolution(solVector,m_ddof,displacement,pressure);
170}
171
172template <class T>
174{
175 if (!initialized)
176 initialize();
177 solVecSaved = solVector;
178 velVecSaved = velVector;
179 accVecSaved = accVector;
180 ddofsSaved = m_ddof;
181 hasSavedState = true;
182}
183
184template <class T>
186{
187 GISMO_ENSURE(hasSavedState,"No state saved!");
188 solVector = solVecSaved;
189 velVector = velVecSaved;
190 accVector = accVecSaved;
191 m_ddof = ddofsSaved;
192}
193
194
195} // namespace ends
gsOptionList m_options
Options.
Definition gsAssembler.h:285
std::vector< gsMatrix< T > > m_ddof
Definition gsAssembler.h:295
Time integation for equations of dynamic elasticity with implicit schemes.
Definition gsElTimeIntegrator.h:29
bool hasSavedState
saved state
Definition gsElTimeIntegrator.h:135
gsMassAssembler< T > & massAssembler
assembler object that generates the mass matrix
Definition gsElTimeIntegrator.h:118
gsElasticityAssembler< T > & stiffAssembler
assembler object that generates the static system
Definition gsElTimeIntegrator.h:116
gsMatrix< T > accVector
vector of acceleration DoFs
Definition gsElTimeIntegrator.h:128
gsElTimeIntegrator(gsElasticityAssembler< T > &stiffAssembler_, gsMassAssembler< T > &massAssembler_)
Definition gsElTimeIntegrator.hpp:28
virtual void assemble()
assemble the linear system for the nonlinear solver
Definition gsBaseAssembler.h:40
index_t numIters
number of iterations Newton's method took to converge at the last time step
Definition gsElTimeIntegrator.h:133
void constructSolution(gsMultiPatch< T > &displacement) const
construct displacement using the stiffness assembler
Definition gsElTimeIntegrator.hpp:159
virtual bool assemble(const gsMatrix< T > &solutionVector, const std::vector< gsMatrix< T > > &fixedDoFs)
Definition gsElTimeIntegrator.hpp:131
void makeTimeStep(T timeStep)
make a time step according to a chosen scheme
Definition gsElTimeIntegrator.hpp:65
gsMatrix< T > implicitLinear()
time integraton schemes
Definition gsElTimeIntegrator.hpp:83
void recoverState()
recover solver state from saved state
Definition gsElTimeIntegrator.hpp:185
gsMatrix< T > solVector
vector of displacement DoFs ( + possibly pressure)
Definition gsElTimeIntegrator.h:124
void saveState()
save solver state
Definition gsElTimeIntegrator.hpp:173
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsElTimeIntegrator.hpp:44
bool initialized
initialization flag
Definition gsElTimeIntegrator.h:120
virtual int numDofs() const
return the number of free degrees of freedom
Definition gsElTimeIntegrator.hpp:128
gsMatrix< T > velVector
vector of velocity DoFs
Definition gsElTimeIntegrator.h:126
Assembles the stiffness matrix and the right-hand side vector for linear and nonlinear elasticity for...
Definition gsElasticityAssembler.h:32
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
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 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
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
Provides time integration for dynamical elasticity.
Provides linear and nonlinear elasticity systems for 2D plain strain and 3D continua.
A class providing an iterative solver for nonlinear problems.
Provides mass matrix for elasticity systems in 2D plain strain and 3D continua.
The G+Smo namespace, containing all definitions for the library.
@ LDLT
LU decomposition: direct, no matrix requirements, robust but a bit slow, Eigen and Pardiso available.
Definition gsBaseUtils.h:70
@ 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