32template <
class T, 
bool _NL>
 
   58                    const Damping_t     & Damping,
 
   59                    const Stiffness_t   & Stiffness,
 
   63    Base(Mass,Damping,Stiffness,Force)
 
 
   69                    const Damping_t     & Damping,
 
   70                    const Stiffness_t   & Stiffness,
 
   71                    const TForce_t      & TForce
 
   74    Base(Mass,Damping,Stiffness,TForce)
 
 
   80                    const Damping_t     & Damping,
 
   81                    const Jacobian_t    & Jacobian,
 
   82                    const Residual_t    & Residual
 
   85    Base(Mass,Damping,Jacobian,Residual)
 
 
   91                    const Damping_t     & Damping,
 
   92                    const Jacobian_t    & Jacobian,
 
   93                    const TResidual_t   & TResidual
 
   96    Base(Mass,Damping,Jacobian,TResidual)
 
 
  102                    const Damping_t     & Damping,
 
  103                    const TJacobian_t   & TJacobian,
 
  104                    const TResidual_t   & TResidual
 
  107    Base(Mass,Damping,TJacobian,TResidual)
 
 
  112                    const TMass_t       & TMass,
 
  113                    const TDamping_t    & TDamping,
 
  114                    const TJacobian_t   & TJacobian,
 
  115                    const TResidual_t   & TResidual
 
  118    Base(TMass,TDamping,TJacobian,TResidual)
 
 
  126    void _initOutput() 
const;
 
  127    void _stepOutput(
const index_t it, 
const T resnorm, 
const T updatenorm) 
const;
 
  140    using Base::m_solver;
 
  144    using Base::m_options;
 
  148    template <
bool _nonlinear>
 
  149    typename std::enable_if<(_nonlinear==
false), 
gsStatus>::type 
 
  152    template <
bool _nonlinear>
 
  153    typename std::enable_if<(_nonlinear==
true), 
gsStatus>::type 
 
 
  159#ifndef GISMO_BUILD_LIB 
  160#include GISMO_HPP_HEADER(gsDynamicImplicitEuler.hpp) 
Performs the arc length method to solve a nonlinear system of equations.
Definition gsDynamicBase.h:35
virtual void _computeResidual(const gsVector< T > &U, const T time, gsVector< T > &R) const
Compute the residual.
Definition gsDynamicBase.h:268
virtual void _computeJacobian(const gsVector< T > &U, const T time, gsSparseMatrix< T > &K) const
Compute the Jacobian matrix.
Definition gsDynamicBase.h:304
index_t m_numIterations
Number of iterations performed.
Definition gsDynamicBase.h:348
virtual void _computeMass(const T time, gsSparseMatrix< T > &M) const
Compute the mass matrix.
Definition gsDynamicBase.h:275
virtual void _computeMassInverse(const gsSparseMatrix< T > &M, gsSparseMatrix< T > &Minv) const
Compute the mass matrix.
Definition gsDynamicBase.h:282
virtual void _computeDamping(const gsVector< T > &U, const T time, gsSparseMatrix< T > &C) const
Compute the damping matrix.
Definition gsDynamicBase.h:297
virtual void _computeForce(const T time, gsVector< T > &F) const
Compute the residual.
Definition gsDynamicBase.h:261
Performs the arc length method to solve a nonlinear system of equations.
Definition gsDynamicImplicitEuler.h:34
gsDynamicImplicitEuler(const Mass_t &Mass, const Damping_t &Damping, const Stiffness_t &Stiffness, const TForce_t &TForce)
Constructor.
Definition gsDynamicImplicitEuler.h:67
gsStatus _step(const T t, const T dt, gsVector< T > &U, gsVector< T > &V, gsVector< T > &A) const override
Initialize the ALM.
Definition gsDynamicImplicitEuler.hpp:167
gsDynamicImplicitEuler(const TMass_t &TMass, const TDamping_t &TDamping, const TJacobian_t &TJacobian, const TResidual_t &TResidual)
Constructor.
Definition gsDynamicImplicitEuler.h:111
gsDynamicImplicitEuler(const Mass_t &Mass, const Damping_t &Damping, const Stiffness_t &Stiffness, const Force_t &Force)
Constructor.
Definition gsDynamicImplicitEuler.h:56
gsDynamicImplicitEuler(const Mass_t &Mass, const Damping_t &Damping, const Jacobian_t &Jacobian, const Residual_t &Residual)
Constructor.
Definition gsDynamicImplicitEuler.h:78
gsDynamicImplicitEuler(const Mass_t &Mass, const Damping_t &Damping, const TJacobian_t &TJacobian, const TResidual_t &TResidual)
Constructor.
Definition gsDynamicImplicitEuler.h:100
gsDynamicImplicitEuler(const Mass_t &Mass, const Damping_t &Damping, const Jacobian_t &Jacobian, const TResidual_t &TResidual)
Constructor.
Definition gsDynamicImplicitEuler.h:89
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
Base class to perform time integration of second-order structural dynamics systems.
This is the main header file that collects wrappers of Eigen for linear algebra.
Provides a list of labeled parameters/options that can be set and accessed easily.
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
std::function< bool(gsVector< T > const &, const T, gsSparseMatrix< T > &) > TJacobian_t
Jacobian.
Definition gsStructuralAnalysisTypes.h:88
std::function< bool(gsVector< T > &)> Force_t
Force.
Definition gsStructuralAnalysisTypes.h:63
std::function< bool(gsVector< T > const &, const T, gsVector< T > &)> TResidual_t
Time-dependent Residual Fint(t)-Fext(t)
Definition gsStructuralAnalysisTypes.h:72
std::function< bool(gsVector< T > const &, gsSparseMatrix< T > &) > Jacobian_t
Jacobian.
Definition gsStructuralAnalysisTypes.h:86
std::function< bool(const T, gsSparseMatrix< T > &) > TMass_t
Time-dependent mass matrix.
Definition gsStructuralAnalysisTypes.h:77
std::function< bool(gsSparseMatrix< T > &) > Mass_t
Mass matrix.
Definition gsStructuralAnalysisTypes.h:75
std::function< bool(const T, gsVector< T > &)> TForce_t
Time-dependent force.
Definition gsStructuralAnalysisTypes.h:65
std::function< bool(gsVector< T > const &, const T, gsSparseMatrix< T > &) > TDamping_t
Time-dependent Damping matrix.
Definition gsStructuralAnalysisTypes.h:81
std::function< bool(gsSparseMatrix< T > &) > Stiffness_t
Stiffness matrix.
Definition gsStructuralAnalysisTypes.h:84
std::function< bool(gsVector< T > const &, gsSparseMatrix< T > &) > Damping_t
Damping matrix.
Definition gsStructuralAnalysisTypes.h:79
std::function< bool(gsVector< T > const &, gsVector< T > &)> Residual_t
Residual, Fint-Fext.
Definition gsStructuralAnalysisTypes.h:68