30 : stiffAssembler(stiffAssembler_),
31 massAssembler(massAssembler_)
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);
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()));
72 newSolVector = implicitLinear();
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;
85 if (massAssembler.numDofs() == stiffAssembler.numDofs())
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();
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);
105 #ifdef GISMO_WITH_PARDISO
106 gsSparseSolver<>::PardisoLDLT solver(m_system.matrix());
107 return solver.solve(m_system.rhs());
109 gsSparseSolver<>::SimplicialLDLT solver(m_system.matrix());
110 return solver.solve(m_system.rhs());
120 solver.options().setInt(
"Verbosity",m_options.getInt(
"Verbosity"));
123 numIters = solver.numberIterations();
124 return solver.solution();
131 stiffAssembler.assemble(solutionVector,fixedDoFs);
132 if (massAssembler.numDofs() == stiffAssembler.numDofs())
134 m_system.matrix() = alpha1()*massAssembler.matrix() + stiffAssembler.matrix();
135 m_system.matrix().makeCompressed();
136 m_system.rhs() = stiffAssembler.rhs() +
137 massAssembler.matrix()*(alpha1()*(solVector-solutionVector) + alpha2()*velVector + alpha3()*accVector);
141 m_system.matrix() = stiffAssembler.matrix();
142 tempMassBlock = alpha1()*massAssembler.matrix();
143 tempMassBlock.conservativeResize(stiffAssembler.numDofs(),massAssembler.numDofs());
144 m_system.matrix().leftCols(massAssembler.numDofs()) += tempMassBlock;
145 m_system.matrix().makeCompressed();
146 m_system.rhs() = stiffAssembler.rhs();
147 m_system.rhs().middleRows(0,massAssembler.numDofs()) +=
148 massAssembler.matrix()*(alpha1()*(solVector-solutionVector).middleRows(0,massAssembler.numDofs())
149 + alpha2()*velVector + alpha3()*accVector);
158 stiffAssembler.constructSolution(solVector,m_ddof,displacement);
164 GISMO_ENSURE(stiffAssembler.numDofs() > massAssembler.numDofs(),
165 "This is a displacement-only formulation. Can't construct pressure");
166 stiffAssembler.constructSolution(solVector,m_ddof,displacement,pressure);
174 solVecSaved = solVector;
175 velVecSaved = velVector;
176 accVecSaved = accVector;
178 hasSavedState =
true;
185 solVector = solVecSaved;
186 velVector = velVecSaved;
187 accVector = accVecSaved;
bool initialized
initialization flag
Definition: gsElTimeIntegrator.h:121
gsElasticityAssembler< T > & stiffAssembler
assembler object that generates the static system
Definition: gsElTimeIntegrator.h:117
implicit scheme with linear problem (theta-scheme)
Definition: gsBaseUtils.h:60
gsOptionList m_options
Options.
Definition: gsAssembler.h:285
Provides mass matrix for elasticity systems in 2D plain strain and 3D continua.
void recoverState()
recover solver state from saved state
Definition: gsElTimeIntegrator.hpp:182
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition: gsElTimeIntegrator.hpp:44
void makeTimeStep(T timeStep)
make a time step according to a chosen scheme
Definition: gsElTimeIntegrator.hpp:65
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition: gsMatrix.h:38
Assembles the stiffness matrix and the right-hand side vector for linear and nonlinear elasticity for...
Definition: gsElasticityAssembler.h:31
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
Assembles the mass matrix and right-hand side vector for linear and nonlinear elasticity for 2D plain...
Definition: gsElTimeIntegrator.h:26
Provides time integration for dynamical elasticity.
gsMassAssembler< T > & massAssembler
assembler object that generates the mass matrix
Definition: gsElTimeIntegrator.h:119
Provides linear and nonlinear elasticity systems for 2D plain strain and 3D continua.
virtual bool assemble(const gsMatrix< T > &solutionVector, const std::vector< gsMatrix< T > > &fixedDoFs)
assemble the linear system for the nonlinear solver
Definition: gsElTimeIntegrator.hpp:128
bool hasSavedState
saved state
Definition: gsElTimeIntegrator.h:136
gsElTimeIntegrator(gsElasticityAssembler< T > &stiffAssembler_, gsMassAssembler< T > &massAssembler_)
Definition: gsElTimeIntegrator.hpp:28
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
A general iterative solver for nonlinear problems. An equation to solve is specified by an assembler ...
Definition: gsALE.h:28
explicit scheme with lumped mass matrix
Definition: gsBaseUtils.h:59
std::vector< gsMatrix< T > > m_ddof
Definition: gsAssembler.h:295
void saveState()
save solver state
Definition: gsElTimeIntegrator.hpp:170
gsMatrix< T > implicitLinear()
time integraton schemes
Definition: gsElTimeIntegrator.hpp:83
index_t numIters
number of iterations Newton's method took to converge at the last time step
Definition: gsElTimeIntegrator.h:134
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 constructSolution(gsMultiPatch< T > &displacement) const
construct displacement using the stiffness assembler
Definition: gsElTimeIntegrator.hpp:156
gsMatrix< T > solVector
vector of displacement DoFs ( + possibly pressure)
Definition: gsElTimeIntegrator.h:125
gsMatrix< T > accVector
vector of acceleration DoFs
Definition: gsElTimeIntegrator.h:129
LU decomposition: direct, no matrix requirements, robust but a bit slow, Eigen and Pardiso available...
Definition: gsBaseUtils.h:70
Time integation for equations of dynamic elasticity with implicit schemes.
Definition: gsElTimeIntegrator.h:31
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
gsMatrix< T > velVector
vector of velocity DoFs
Definition: gsElTimeIntegrator.h:127
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition: gsAssembler.hpp:51
A class providing an iterative solver for nonlinear problems.