FrictionQPotFEM 0.23.3
Loading...
Searching...
No Matches
FrictionQPotFEM::Generic2d::System Class Reference

System with elastic elements and plastic elements (GMatElastoPlasticQPot). More...

#include <FrictionQPotFEM/Generic2d.h>

Inheritance diagram for FrictionQPotFEM::Generic2d::System:
Collaboration diagram for FrictionQPotFEM::Generic2d::System:

Public Member Functions

template<class C , class E , class L , class M , class Y >
 System (const C &coor, const E &conn, const E &dofs, const L &iip, const L &elastic_elem, const M &elastic_K, const M &elastic_G, const L &plastic_elem, const M &plastic_K, const M &plastic_G, const Y &plastic_epsy, double dt, double rho, double alpha, double eta)
 Define the geometry, including boundary conditions and element sets. More...
 
virtual size_t N () const
 Return the linear system size (in number of blocks). More...
 
virtual std::string type () const
 Return the type of system. More...
 
double rho () const
 Mass density. More...
 
void setRho (double rho)
 Overwrite the mass density to a homogeneous quantity. More...
 
template<class T >
void setMassMatrix (const T &val_elem)
 Overwrite mass matrix, based on certain density that is uniform per element. More...
 
void setEta (double eta)
 Overwrite the value of the damping at the interface. More...
 
double eta () const
 Get the damping at the interface. More...
 
void setAlpha (double alpha)
 Overwrite background damping density (proportional to the velocity), To use a non-homogeneous density use setDampingMatrix(). More...
 
double alpha () const
 Background damping density. More...
 
template<class T >
void setDampingMatrix (const T &val_elem)
 Overwrite damping matrix, based on certain density (taken uniform per element). More...
 
bool isHomogeneousElastic () const
 Check if elasticity is homogeneous. More...
 
double dt () const
 Get time step. More...
 
void setDt (double dt)
 Overwrite time step. More...
 
template<class T >
void setU (const T &u)
 Overwrite nodal displacements. More...
 
template<class T >
void setV (const T &v)
 Overwrite nodal velocities. More...
 
template<class T >
void setA (const T &a)
 Overwrite nodal accelerations. More...
 
virtual void updated_u ()
 Evaluate relevant forces when m_u is updated. More...
 
void updated_v ()
 Evaluate relevant forces when m_v is updated. More...
 
template<class T >
void setFext (const T &fext)
 Overwrite external force. More...
 
const auto & elastic_elem () const
 List of elastic elements. More...
 
const auto & plastic_elem () const
 List of plastic elements. More...
 
const auto & conn () const
 Connectivity. More...
 
const auto & coor () const
 Nodal coordinates. More...
 
const auto & dofs () const
 DOFs per node. More...
 
const auto & u () const
 Nodal displacements. More...
 
const auto & v () const
 Nodal velocities. More...
 
const auto & a () const
 Nodal accelerations. More...
 
auto & mass () const
 Mass matrix, see System::m_M. More...
 
auto & damping () const
 Damping matrix, see setAlpha() and setDampingMatrix(). More...
 
const auto & fext ()
 External force. More...
 
void applyShearStress (double sigxy)
 Modify the external force such that a shear stress is applied to the top boundary. More...
 
const auto & fint ()
 Internal force. More...
 
const auto & fmaterial () const
 Force deriving from elasticity. More...
 
const auto & fdamp () const
 Force deriving from damping. More...
 
double residual () const
 Norm of the relative residual force (the external force at the top/bottom boundaries is used for normalisation). More...
 
void quench ()
 Set nodal velocities and accelerations equal to zero. More...
 
double t () const
 Current time. More...
 
void setT (double arg)
 Overwrite time. More...
 
size_t inc () const
 The increment, see setInc(). More...
 
virtual void setInc (size_t arg)
 Set increment. More...
 
const auto & dV () const
 Integration point volume (of each integration point) More...
 
const GooseFEM::VectorPartitionedvector () const
 GooseFEM vector definition. More...
 
const GooseFEM::Element::Quad4::Quadraturequad () const
 GooseFEM quadrature definition. More...
 
GMatElastoPlasticQPot::Cartesian2d::Elastic< 2 > & elastic ()
 Elastic material mode, used for all elastic elements. More...
 
GMatElastoPlasticQPot::Cartesian2d::Cusp< 2 > & plastic ()
 Elastic material mode, used for all elastic elements. More...
 
array_type::tensor< double, 2 > K () const
 Bulk modulus per integration point. More...
 
array_type::tensor< double, 2 > G () const
 Shear modulus per integration point. More...
 
const array_type::tensor< double, 4 > & Sig ()
 Stress tensor of each integration point. More...
 
const array_type::tensor< double, 4 > & Eps ()
 Strain tensor of each integration point. More...
 
array_type::tensor< double, 4 > Epsdot () const
 Strain-rate tensor (the symmetric gradient of the nodal velocities) of each integration point. More...
 
array_type::tensor< double, 4 > Epsddot () const
 Symmetric gradient of the nodal accelerations of each integration point. More...
 
GooseFEM::MatrixPartitioned stiffness () const
 Stiffness tensor of each integration point. More...
 
const array_type::tensor< double, 4 > & plastic_Epsdot ()
 Strain-rate tensor of integration points of plastic elements only, see System::plastic. More...
 
void refresh ()
 Recompute all forces. More...
 
template<class T >
double eventDriven_setDeltaU (const T &delta_u, bool autoscale=true)
 Set purely elastic and linear response to specific boundary conditions. More...
 
const auto & eventDriven_deltaU () const
 Get displacement perturbation used for event driven code, see eventDriven_setDeltaU(). More...
 
virtual double eventDrivenStep (double deps, bool kick, int direction=1, bool yield_element=false, bool iterative=false)
 Add event driven step for the boundary conditions that correspond to the displacement perturbation set in eventDriven_setDeltaU(). More...
 
void timeStep ()
 Make a time-step: apply velocity-Verlet integration. More...
 
long timeSteps (size_t n, size_t nmargin=0)
 Make a number of time steps, see timeStep(). More...
 
size_t timeStepsUntilEvent (size_t nmargin=0, double tol=1e-5, size_t niter_tol=20, size_t max_iter=1e7)
 Perform a series of time-steps until the next plastic event, or equilibrium. More...
 
std::tuple< array_type::tensor< double, 1 >, array_type::tensor< double, 2 >, array_type::tensor< double, 2 > > minimise_highfrequency (const array_type::tensor< size_t, 1 > &nodes, const array_type::tensor< size_t, 1 > &top, size_t t_step=1, size_t nmargin=0, double tol=1e-5, size_t niter_tol=20, size_t max_iter=1e7)
 Minimise energy while doing a high-frequency measurement of: More...
 
template<class T >
long flowSteps (size_t n, const T &v, size_t nmargin=0)
 Make a number of steps with the following protocol. More...
 
long minimise (size_t nmargin=0, double tol=1e-5, size_t niter_tol=20, size_t max_iter=1e7, bool time_activity=false, bool max_iter_is_error=true)
 Minimise energy: run System::timeStep until a mechanical equilibrium has been reached. More...
 
size_t quasistaticActivityFirst () const
 Increment with the first plastic event. More...
 
size_t quasistaticActivityLast () const
 Increment with the last plastic event. More...
 
template<class T >
size_t minimise_truncate (const T &idx_n, size_t A_truncate=0, size_t S_truncate=0, double tol=1e-5, size_t niter_tol=20, size_t max_iter=1e7)
 Like Generic2d::minimise(), but stops when a certain number of blocks has yielded at least once. More...
 
size_t minimise_truncate (size_t A_truncate=0, size_t S_truncate=0, double tol=1e-5, size_t niter_tol=20, size_t max_iter=1e7)
 Like Generic2d::minimise(), but stops when a certain number of blocks has yielded at least once. More...
 
array_type::tensor< double, 2 > affineSimpleShear (double delta_gamma) const
 Get the displacement field that corresponds to an affine simple shear of a certain strain. More...
 
array_type::tensor< double, 2 > affineSimpleShearCentered (double delta_gamma) const
 Get the displacement field that corresponds to an affine simple shear of a certain strain. More...
 

Protected Member Functions

void computeEpsSig ()
 Compute m_Sig and m_Eps using m_u. More...
 
void computeForceMaterial ()
 Update m_fmaterial based on the current displacement field m_u. More...
 
virtual void computeInternalExternalResidualForce ()
 Compute: More...
 

Protected Attributes

array_type::tensor< double, 2 > m_coor
 Nodal coordinates, see coor(). More...
 
size_t m_N
 Number of plastic elements, alias of m_nelem_plas. More...
 
size_t m_nelem
 Number of elements. More...
 
size_t m_nelem_elas
 Number of elastic elements. More...
 
size_t m_nelem_plas
 Number of plastic elements. More...
 
size_t m_nne
 Number of nodes per element. More...
 
size_t m_nnode
 Number of nodes. More...
 
size_t m_ndim
 Number of spatial dimensions. More...
 
size_t m_nip
 Number of integration points. More...
 
array_type::tensor< size_t, 1 > m_elem_elas
 Elastic elements. More...
 
array_type::tensor< size_t, 1 > m_elem_plas
 Plastic elements. More...
 
GooseFEM::Element::Quad4::Quadrature m_quad
 Quadrature for all elements. More...
 
GooseFEM::Element::Quad4::Quadrature m_quad_elas
 m_quad for elastic elements only. More...
 
GooseFEM::Element::Quad4::Quadrature m_quad_plas
 m_quad for plastic elements only. More...
 
GooseFEM::VectorPartitioned m_vector
 Convert vectors between 'nodevec', 'elemvec', .... More...
 
GooseFEM::VectorPartitioned m_vector_elas
 m_vector for elastic elements only. More...
 
GooseFEM::VectorPartitioned m_vector_plas
 m_vector for plastic elements only. More...
 
GooseFEM::MatrixDiagonalPartitioned m_M
 Mass matrix (diagonal). More...
 
GooseFEM::MatrixDiagonal m_D
 Damping matrix (diagonal). More...
 
GMatElastoPlasticQPot::Cartesian2d::Elastic< 2 > m_elas
 Material for elastic el. More...
 
GMatElastoPlasticQPot::Cartesian2d::Cusp< 2 > m_plas
 Material for plastic el. More...
 
array_type::tensor< double, 2 > m_u
 Nodal displacements. More...
 
array_type::tensor< double, 2 > m_v
 Nodal velocities. More...
 
array_type::tensor< double, 2 > m_a
 Nodal accelerations. More...
 
array_type::tensor< double, 2 > m_v_n
 Nodal velocities last time-step. More...
 
array_type::tensor< double, 2 > m_a_n
 Nodal accelerations last time-step. More...
 
array_type::tensor< double, 3 > m_ue
 Element vector (used for displacements). More...
 
array_type::tensor< double, 3 > m_fe
 Element vector (used for forces). More...
 
array_type::tensor< double, 3 > m_ue_elas
 m_ue for elastic element only More...
 
array_type::tensor< double, 3 > m_fe_elas
 m_fe for elastic element only More...
 
array_type::tensor< double, 3 > m_ue_plas
 m_ue for plastic element only More...
 
array_type::tensor< double, 3 > m_fe_plas
 m_fe for plastic element only More...
 
array_type::tensor< double, 2 > m_fmaterial
 Nodal force, deriving from elasticity. More...
 
array_type::tensor< double, 2 > m_felas
 Nodal force from elasticity of elastic elements. More...
 
array_type::tensor< double, 2 > m_fplas
 Nodal force from elasticity of plastic elements. More...
 
array_type::tensor< double, 2 > m_fdamp
 Nodal force from damping. More...
 
array_type::tensor< double, 2 > m_fvisco
 Nodal force from damping at the interface. More...
 
array_type::tensor< double, 2 > m_ftmp
 Temporary for internal use. More...
 
array_type::tensor< double, 2 > m_fint
 Nodal force: total internal force. More...
 
array_type::tensor< double, 2 > m_fext
 Nodal force: total external force (reaction force) More...
 
array_type::tensor< double, 2 > m_fres
 Nodal force: residual force. More...
 
array_type::tensor< double, 4 > m_Eps
 Quad-point tensor: strain. More...
 
array_type::tensor< double, 4 > m_Sig
 Quad-point tensor: stress. More...
 
array_type::tensor< double, 4 > m_Epsdot_plas
 Quad-point tensor: strain-rate for plastic el. More...
 
GooseFEM::Matrix m_K_elas
 Stiffness matrix for elastic elements only. More...
 
size_t m_qs_inc_first = 0
 First increment with plastic activity during minimisation. More...
 
size_t m_qs_inc_last = 0
 Last increment with plastic activity during minimisation. More...
 
size_t m_inc
 Current increment (current time = m_dt * m_inc). More...
 
double m_dt
 Time-step. More...
 
double m_eta
 Damping at the interface. More...
 
double m_rho
 Mass density (non-zero only if homogeneous). More...
 
double m_alpha
 Background damping density (non-zero only if homogeneous). More...
 
bool m_set_D
 Use m_D only if it is non-zero. More...
 
bool m_set_visco
 Use m_eta only if it is non-zero. More...
 
bool m_full_outdated
 Keep track of the need to recompute fields on full geometry. More...
 
array_type::tensor< double, 2 > m_pert_u
 See eventDriven_setDeltaU() More...
 
array_type::tensor< double, 4 > m_pert_Epsd_plastic
 Strain deviator for m_pert_u. More...
 

Detailed Description

System with elastic elements and plastic elements (GMatElastoPlasticQPot).

For efficiency, the nodal forces for the elastic elements are evaluated using the tangent. This means that getting stresses and strains in those elements is not for free. Therefore, there are separate methods to get stresses and strains only in the plastic elements (as they are readily available as they are needed for the force computation).

Definition at line 183 of file Generic2d.h.

Constructor & Destructor Documentation

◆ ~System()

virtual FrictionQPotFEM::Generic2d::System::~System ( )
inlinevirtual

Definition at line 189 of file Generic2d.h.

◆ System()

template<class C , class E , class L , class M , class Y >
FrictionQPotFEM::Generic2d::System::System ( const C &  coor,
const E &  conn,
const E &  dofs,
const L &  iip,
const L &  elastic_elem,
const M &  elastic_K,
const M &  elastic_G,
const L &  plastic_elem,
const M &  plastic_K,
const M &  plastic_G,
const Y &  plastic_epsy,
double  dt,
double  rho,
double  alpha,
double  eta 
)
inline

Define the geometry, including boundary conditions and element sets.

Template Parameters
CType of nodal coordinates, e.g. array_type::tensor<double, 2>
EType of connectivity and DOFs, e.g. array_type::tensor<size_t, 2>
LType of node/element lists, e.g. array_type::tensor<size_t, 1>
Parameters
coorNodal coordinates.
connConnectivity.
dofsDOFs per node.
iipDOFs whose displacement is fixed.
elastic_elemElastic elements.
elastic_KBulk modulus per quad. point of each elastic element, see setElastic().
elastic_GShear modulus per quad. point of each elastic element, see setElastic().
plastic_elemPlastic elements.
plastic_KBulk modulus per quad. point of each plastic element, see Plastic().
plastic_GShear modulus per quad. point of each plastic element, see Plastic().
plastic_epsyYield strain per quad. point of each plastic element, see Plastic().
dtTime step, set setDt().
rhoMass density, see setMassMatrix().
alphaBackground damping density, see setDampingMatrix().
etaDamping at the interface (homogeneous), see setEta().

Definition at line 215 of file Generic2d.h.

Member Function Documentation

◆ a()

const auto & FrictionQPotFEM::Generic2d::System::a ( ) const
inline

Nodal accelerations.

Returns
Nodal accelerations.

Definition at line 748 of file Generic2d.h.

◆ affineSimpleShear()

array_type::tensor< double, 2 > FrictionQPotFEM::Generic2d::System::affineSimpleShear ( double  delta_gamma) const
inline

Get the displacement field that corresponds to an affine simple shear of a certain strain.

The displacement of the bottom boundary is zero, while it is maximal for the top boundary.

Parameters
delta_gammaStrain to add (the shear component of deformation gradient is twice that).
Returns
Nodal displacements.

Definition at line 2041 of file Generic2d.h.

◆ affineSimpleShearCentered()

array_type::tensor< double, 2 > FrictionQPotFEM::Generic2d::System::affineSimpleShearCentered ( double  delta_gamma) const
inline

Get the displacement field that corresponds to an affine simple shear of a certain strain.

Similar to affineSimpleShear() with the difference that the displacement is zero exactly in the middle, while the displacement at the bottom and the top boundary is maximal (with a negative displacement for the bottom boundary).

Parameters
delta_gammaStrain to add (the shear component of deformation gradient is twice that).
Returns
Nodal displacements.

Definition at line 2061 of file Generic2d.h.

◆ alpha()

double FrictionQPotFEM::Generic2d::System::alpha ( ) const
inline

Background damping density.

The output is non-zero only if the density is homogeneous. I.e. a zero value does not mean that there is no damping.

Returns
Float

Definition at line 507 of file Generic2d.h.

◆ applyShearStress()

void FrictionQPotFEM::Generic2d::System::applyShearStress ( double  sigxy)
inline

Modify the external force such that a shear stress is applied to the top boundary.

Parameters
sigxyThe shear stress to be applied.

Definition at line 789 of file Generic2d.h.

◆ computeEpsSig()

void FrictionQPotFEM::Generic2d::System::computeEpsSig ( )
inlineprotected

Compute m_Sig and m_Eps using m_u.

Definition at line 1134 of file Generic2d.h.

◆ computeForceMaterial()

void FrictionQPotFEM::Generic2d::System::computeForceMaterial ( )
inlineprotected

Update m_fmaterial based on the current displacement field m_u.

Definition at line 1179 of file Generic2d.h.

◆ computeInternalExternalResidualForce()

virtual void FrictionQPotFEM::Generic2d::System::computeInternalExternalResidualForce ( )
inlineprotectedvirtual

Compute:

  • m_fint = m_fmaterial + m_fdamp
  • m_fext[iip] = m_fint[iip]
  • m_fres = m_fext - m_fint

Internal rule: all relevant forces are expected to be updated before this function is called.

Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System, and FrictionQPotFEM::UniformSingleLayerThermal2d::System.

Definition at line 1203 of file Generic2d.h.

◆ conn()

const auto & FrictionQPotFEM::Generic2d::System::conn ( ) const
inline

Connectivity.

Shape: [System::m_nelem, System::m_nne].

Returns
Connectivity.

Definition at line 698 of file Generic2d.h.

◆ coor()

const auto & FrictionQPotFEM::Generic2d::System::coor ( ) const
inline

Nodal coordinates.

Shape: [System::m_nnode, System::m_ndim].

Returns
Nodal coordinates.

Definition at line 708 of file Generic2d.h.

◆ damping()

auto & FrictionQPotFEM::Generic2d::System::damping ( ) const
inline

Damping matrix, see setAlpha() and setDampingMatrix().

Note that there can be second source of damping, see setEta().

Returns
Damping matrix (diagonal).

Definition at line 769 of file Generic2d.h.

◆ dofs()

const auto & FrictionQPotFEM::Generic2d::System::dofs ( ) const
inline

DOFs per node.

Returns
DOFs per node.

Definition at line 718 of file Generic2d.h.

◆ dt()

double FrictionQPotFEM::Generic2d::System::dt ( ) const
inline

Get time step.

Definition at line 566 of file Generic2d.h.

◆ dV()

const auto & FrictionQPotFEM::Generic2d::System::dV ( ) const
inline

Integration point volume (of each integration point)

Returns
Integration point volume (System::m_quad::dV).

Definition at line 928 of file Generic2d.h.

◆ elastic()

GMatElastoPlasticQPot::Cartesian2d::Elastic< 2 > & FrictionQPotFEM::Generic2d::System::elastic ( )
inline

Elastic material mode, used for all elastic elements.

Returns
GMatElastoPlasticQPot::Cartesian2d::Elastic<2>&

Definition at line 959 of file Generic2d.h.

◆ elastic_elem()

const auto & FrictionQPotFEM::Generic2d::System::elastic_elem ( ) const
inline

List of elastic elements.

Shape: [System::m_nelem_elas].

Returns
List of element numbers.

Definition at line 678 of file Generic2d.h.

◆ Eps()

const array_type::tensor< double, 4 > & FrictionQPotFEM::Generic2d::System::Eps ( )
inline

Strain tensor of each integration point.

Returns
Integration point tensor (pointer). Shape: [nelem, nip, 2, 2].

Definition at line 1054 of file Generic2d.h.

◆ Epsddot()

array_type::tensor< double, 4 > FrictionQPotFEM::Generic2d::System::Epsddot ( ) const
inline

Symmetric gradient of the nodal accelerations of each integration point.

Returns
Integration point tensor (copy). Shape: [nelem, nip, 2, 2].

Definition at line 1076 of file Generic2d.h.

◆ Epsdot()

array_type::tensor< double, 4 > FrictionQPotFEM::Generic2d::System::Epsdot ( ) const
inline

Strain-rate tensor (the symmetric gradient of the nodal velocities) of each integration point.

Returns
Integration point tensor (copy). Shape: [nelem, nip, 2, 2].

Definition at line 1066 of file Generic2d.h.

◆ eta()

double FrictionQPotFEM::Generic2d::System::eta ( ) const
inline

Get the damping at the interface.

Returns
Float

Definition at line 478 of file Generic2d.h.

◆ eventDriven_deltaU()

const auto & FrictionQPotFEM::Generic2d::System::eventDriven_deltaU ( ) const
inline

Get displacement perturbation used for event driven code, see eventDriven_setDeltaU().

Returns
Nodal displacements.

Definition at line 1261 of file Generic2d.h.

◆ eventDriven_setDeltaU()

template<class T >
double FrictionQPotFEM::Generic2d::System::eventDriven_setDeltaU ( const T &  delta_u,
bool  autoscale = true 
)
inline

Set purely elastic and linear response to specific boundary conditions.

Since this response is linear and elastic it can be scaled freely to transverse a fully elastic interval at once, without running any dynamics, and run an event driven code using eventDrivenStep().

Parameters
delta_uNodal displacement field.
autoscaleScale such that the perturbation is small enough.
Returns
Value with which the input perturbation is scaled, see also eventDriven_deltaU().

Definition at line 1234 of file Generic2d.h.

◆ eventDrivenStep()

virtual double FrictionQPotFEM::Generic2d::System::eventDrivenStep ( double  deps,
bool  kick,
int  direction = 1,
bool  yield_element = false,
bool  iterative = false 
)
inlinevirtual

Add event driven step for the boundary conditions that correspond to the displacement perturbation set in eventDriven_setDeltaU().

Todo:
The current implementation is suited for unloading, but only until the first yield strain. An improved implementation is needed to deal with the symmetry of the equivalent strain. A FRICTIONQPOTFEM_WIP assertion is made.
Parameters
depsSize of the local stain kick to apply.
kickIf false, increment displacements to deps / 2 of yielding again. If true, increment displacements by a affine simple shear increment deps.
directionIf +1: apply shear to the right. If -1 applying shear to the left.
yield_elementIf true and kick == true then the element closest to yielding is selected (as normal), but of that element the displacement update is the maximum of the element, such that all integration points of the element are forced to yield.
iterativeIf true the step is iteratively searched. This is more costly by recommended, if the perturbation is non-analytical (and contains numerical errors).
Returns
Factor with which the displacement perturbation, see eventDriven_deltaU(), is scaled.

Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System, and FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System.

Definition at line 1300 of file Generic2d.h.

◆ fdamp()

const auto & FrictionQPotFEM::Generic2d::System::fdamp ( ) const
inline

Force deriving from damping.

This force is the sum of damping at the interface plus that of background damping (or just one of both if just one is specified).

Returns
Nodal force. Shape [nnode, ndim] .

Definition at line 854 of file Generic2d.h.

◆ fext()

const auto & FrictionQPotFEM::Generic2d::System::fext ( )
inline

External force.

Returns
Nodal force. Shape [nnode, ndim] .

Definition at line 779 of file Generic2d.h.

◆ fint()

const auto & FrictionQPotFEM::Generic2d::System::fint ( )
inline

Internal force.

Returns
Nodal force. Shape [nnode, ndim] .

Definition at line 831 of file Generic2d.h.

◆ flowSteps()

template<class T >
long FrictionQPotFEM::Generic2d::System::flowSteps ( size_t  n,
const T &  v,
size_t  nmargin = 0 
)
inline

Make a number of steps with the following protocol.

  1. Add a displacement \( \underline{v} \Delta t \) to each of the nodes.
  2. Make a timeStep().
Parameters
nNumber of steps to make.
vNodal velocity [nnode, ndim].s
nmarginNumber of potentials to leave as margin.
  • nmargin > 0: stop if the yield-index of any box is nmargin from the end, and return a negative number.
  • nmargin == 0: no bounds-check is performed (the last potential is assumed infinitely elastic to the right).
Returns
Number of time-steps made (negative if failure).

Definition at line 1803 of file Generic2d.h.

◆ fmaterial()

const auto & FrictionQPotFEM::Generic2d::System::fmaterial ( ) const
inline

Force deriving from elasticity.

Returns
Nodal force. Shape [nnode, ndim] .

Definition at line 842 of file Generic2d.h.

◆ G()

array_type::tensor< double, 2 > FrictionQPotFEM::Generic2d::System::G ( ) const
inline

Shear modulus per integration point.

Convenience function: assembles output from elastic() and plastic().

Returns
Integration point scalar (copy). Shape: [nelem, nip].

Definition at line 1012 of file Generic2d.h.

◆ inc()

size_t FrictionQPotFEM::Generic2d::System::inc ( ) const
inline

The increment, see setInc().

Returns
size_t.

Definition at line 909 of file Generic2d.h.

◆ isHomogeneousElastic()

bool FrictionQPotFEM::Generic2d::System::isHomogeneousElastic ( ) const
inline

Check if elasticity is homogeneous.

Returns
true is elasticity is homogeneous (false otherwise).

Definition at line 554 of file Generic2d.h.

◆ K()

array_type::tensor< double, 2 > FrictionQPotFEM::Generic2d::System::K ( ) const
inline

Bulk modulus per integration point.

Convenience function: assembles output from elastic() and plastic().

Returns
Integration point scalar (copy). Shape: [nelem, nip].

Definition at line 980 of file Generic2d.h.

◆ mass()

auto & FrictionQPotFEM::Generic2d::System::mass ( ) const
inline

Mass matrix, see System::m_M.

Returns
Mass matrix (diagonal, partitioned).

Definition at line 758 of file Generic2d.h.

◆ minimise()

long FrictionQPotFEM::Generic2d::System::minimise ( size_t  nmargin = 0,
double  tol = 1e-5,
size_t  niter_tol = 20,
size_t  max_iter = 1e7,
bool  time_activity = false,
bool  max_iter_is_error = true 
)
inline

Minimise energy: run System::timeStep until a mechanical equilibrium has been reached.

Can also be used to run n time-steps (or < n is equilibrium is reached before). In that case, call with minimise(max_iter=n, max_iter_is_error=False)

Parameters
nmarginNumber of potentials to leave as margin.
  • nmargin > 0: stop if the yield-index of any box is nmargin from the end and return a negative number.
  • nmargin == 0: no bounds-check is performed (the last potential is assumed infinitely elastic to the right).
tolRelative force tolerance for equilibrium. See System::residual for definition.
niter_tolEnforce the residual check for niter_tol consecutive increments.
max_iterMaximum number of time-steps. Throws std::runtime_error otherwise.
time_activityIf true plastic activity is timed. After this function you can find:
max_iter_is_errorIf true an error is thrown when the maximum number of time-steps is reached. If false the function simply returns max_iter.
Returns
  • 0: if stopped when the residual is reached (and number of steps < max_iter).
  • max_iter: if no residual was reached, and max_iter_is_error = false.
  • Negative number: if stopped because of a yield-index margin.

Definition at line 1868 of file Generic2d.h.

◆ minimise_highfrequency()

std::tuple< array_type::tensor< double, 1 >, array_type::tensor< double, 2 >, array_type::tensor< double, 2 > > FrictionQPotFEM::Generic2d::System::minimise_highfrequency ( const array_type::tensor< size_t, 1 > &  nodes,
const array_type::tensor< size_t, 1 > &  top,
size_t  t_step = 1,
size_t  nmargin = 0,
double  tol = 1e-5,
size_t  niter_tol = 20,
size_t  max_iter = 1e7 
)
inline

Minimise energy while doing a high-frequency measurement of:

  • The average force at the top boundary.
  • The displacement, velocity, and acceleration of selected nodes.
Parameters
nodesList of nodes to sample.
topNodes of the top boundary. f_x = sum(fext[top, 0]) / (top.size - 1) (accounting for periodicity).
t_stepSample every t_step time steps.
nmarginNumber of potentials to have as margin initially.
tolRelative force tolerance for equilibrium. See System::residual for definition.
niter_tolEnforce the residual check for niter_tol consecutive increments.
max_iterMaximum number of iterations. Throws std::runtime_error otherwise.
Returns
A tuple (f_x, u_x, u_y). For the latter three the rows correspond to the selected nodes.

Definition at line 1702 of file Generic2d.h.

◆ minimise_truncate() [1/2]

template<class T >
size_t FrictionQPotFEM::Generic2d::System::minimise_truncate ( const T &  idx_n,
size_t  A_truncate = 0,
size_t  S_truncate = 0,
double  tol = 1e-5,
size_t  niter_tol = 20,
size_t  max_iter = 1e7 
)
inline

Like Generic2d::minimise(), but stops when a certain number of blocks has yielded at least once.

In that case the function returns zero (in all other cases it returns a positive number).

Note
A_truncate and S_truncate are defined on the first integration point.
Parameters
A_truncateTruncate if A_truncate blocks have yielded at least once.
S_truncateTruncate if the number of times blocks yielded is equal to S_truncate. Warning This option is reserved for future use, but for the moment does nothing.
tolRelative force tolerance for equilibrium. See System::residual for definition.
niter_tolEnforce the residual check for niter_tol consecutive increments.
max_iterMaximum number of time-steps. Throws std::runtime_error otherwise.

This overload can be used to specify a reference state when manually triggering. If triggering is done before calling this function, already one block yielded, making A_truncate and S_truncate inaccurate.

Parameters
idx_nReference potential index of the first integration point.

Definition at line 1962 of file Generic2d.h.

◆ minimise_truncate() [2/2]

size_t FrictionQPotFEM::Generic2d::System::minimise_truncate ( size_t  A_truncate = 0,
size_t  S_truncate = 0,
double  tol = 1e-5,
size_t  niter_tol = 20,
size_t  max_iter = 1e7 
)
inline

Like Generic2d::minimise(), but stops when a certain number of blocks has yielded at least once.

In that case the function returns zero (in all other cases it returns a positive number).

Note
A_truncate and S_truncate are defined on the first integration point.
Parameters
A_truncateTruncate if A_truncate blocks have yielded at least once.
S_truncateTruncate if the number of times blocks yielded is equal to S_truncate. Warning This option is reserved for future use, but for the moment does nothing.
tolRelative force tolerance for equilibrium. See System::residual for definition.
niter_tolEnforce the residual check for niter_tol consecutive increments.
max_iterMaximum number of time-steps. Throws std::runtime_error otherwise.

Definition at line 2023 of file Generic2d.h.

◆ N()

virtual size_t FrictionQPotFEM::Generic2d::System::N ( ) const
inlinevirtual

Return the linear system size (in number of blocks).

Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System.

Definition at line 385 of file Generic2d.h.

◆ plastic()

GMatElastoPlasticQPot::Cartesian2d::Cusp< 2 > & FrictionQPotFEM::Generic2d::System::plastic ( )
inline

Elastic material mode, used for all elastic elements.

Returns
GMatElastoPlasticQPot::Cartesian2d::Cusp<2>&

Definition at line 969 of file Generic2d.h.

◆ plastic_elem()

const auto & FrictionQPotFEM::Generic2d::System::plastic_elem ( ) const
inline

List of plastic elements.

Shape: [System::m_nelem_plas].

Returns
List of element numbers.

Definition at line 688 of file Generic2d.h.

◆ plastic_Epsdot()

const array_type::tensor< double, 4 > & FrictionQPotFEM::Generic2d::System::plastic_Epsdot ( )
inline

Strain-rate tensor of integration points of plastic elements only, see System::plastic.

Returns
Integration point tensor (pointer). Shape: [plastic_elem().size(), nip, 2, 2].

Definition at line 1119 of file Generic2d.h.

◆ quad()

const GooseFEM::Element::Quad4::Quadrature & FrictionQPotFEM::Generic2d::System::quad ( ) const
inline

GooseFEM quadrature definition.

Takes case of interpolation, integration, and differentiation.

Returns
GooseFEM::Element::Quad4::Quadrature (System::m_quad)

Definition at line 949 of file Generic2d.h.

◆ quasistaticActivityFirst()

size_t FrictionQPotFEM::Generic2d::System::quasistaticActivityFirst ( ) const
inline

Increment with the first plastic event.

This value is only relevant if time_activity = true was used in the last call of minimise().

Returns
Increment.

Definition at line 1936 of file Generic2d.h.

◆ quasistaticActivityLast()

size_t FrictionQPotFEM::Generic2d::System::quasistaticActivityLast ( ) const
inline

Increment with the last plastic event.

This value is only relevant if time_activity = true was used in the last call of minimise().

Returns
Increment.

Definition at line 1947 of file Generic2d.h.

◆ quench()

void FrictionQPotFEM::Generic2d::System::quench ( )
inline

Set nodal velocities and accelerations equal to zero.

Call this function after an energy minimisation (taken care of in System::minimise).

Definition at line 880 of file Generic2d.h.

◆ refresh()

void FrictionQPotFEM::Generic2d::System::refresh ( )
inline

Recompute all forces.

Under normal circumstances, this function is called automatically when needed.

Definition at line 1215 of file Generic2d.h.

◆ residual()

double FrictionQPotFEM::Generic2d::System::residual ( ) const
inline

Norm of the relative residual force (the external force at the top/bottom boundaries is used for normalisation).

Returns
Relative residual.

Definition at line 865 of file Generic2d.h.

◆ rho()

double FrictionQPotFEM::Generic2d::System::rho ( ) const
inline

Mass density.

The output is non-zero only if the density is homogeneous. I.e. a zero value does not mean that there is no mass.

Returns
Float

Definition at line 406 of file Generic2d.h.

◆ setA()

template<class T >
void FrictionQPotFEM::Generic2d::System::setA ( const T &  a)
inline

Overwrite nodal accelerations.

Parameters
a[nnode, ndim].

Definition at line 615 of file Generic2d.h.

◆ setAlpha()

void FrictionQPotFEM::Generic2d::System::setAlpha ( double  alpha)
inline

Overwrite background damping density (proportional to the velocity), To use a non-homogeneous density use setDampingMatrix().

Parameters
alphaDamping parameter.

Definition at line 489 of file Generic2d.h.

◆ setDampingMatrix()

template<class T >
void FrictionQPotFEM::Generic2d::System::setDampingMatrix ( const T &  val_elem)
inline

Overwrite damping matrix, based on certain density (taken uniform per element).

This function allows for more heterogeneity than the constructor. Note that you can specify either setEta() or setDampingMatrix() or both. To use a homogeneous system, use setAlpha().

Parameters
val_elemDamping per element.

Definition at line 521 of file Generic2d.h.

◆ setDt()

void FrictionQPotFEM::Generic2d::System::setDt ( double  dt)
inline

Overwrite time step.

Using for example in System::timeStep and System::minimise.

Definition at line 574 of file Generic2d.h.

◆ setEta()

void FrictionQPotFEM::Generic2d::System::setEta ( double  eta)
inline

Overwrite the value of the damping at the interface.

Note that you can specify either setEta() or setDampingMatrix() or both.

Parameters
etaDamping parameter

Definition at line 461 of file Generic2d.h.

◆ setFext()

template<class T >
void FrictionQPotFEM::Generic2d::System::setFext ( const T &  fext)
inline

Overwrite external force.

Note: the external force on the DOFs whose displacement are prescribed are response forces computed during timeStep(). Internally on the system of unknown DOFs is solved, so any change to the response forces is ignored.

Parameters
fext[nnode, ndim].

Definition at line 667 of file Generic2d.h.

◆ setInc()

virtual void FrictionQPotFEM::Generic2d::System::setInc ( size_t  arg)
inlinevirtual

Set increment.

Parameters
argsize_t.

Reimplemented in FrictionQPotFEM::UniformSingleLayerThermal2d::System.

Definition at line 918 of file Generic2d.h.

◆ setMassMatrix()

template<class T >
void FrictionQPotFEM::Generic2d::System::setMassMatrix ( const T &  val_elem)
inline

Overwrite mass matrix, based on certain density that is uniform per element.

This function allows for more heterogeneity than the constructor. To use a homogeneous system, use setRho().

Template Parameters
Te.g. array_type::tensor<double, 1>.
Parameters
val_elemDensity per element.

Definition at line 431 of file Generic2d.h.

◆ setRho()

void FrictionQPotFEM::Generic2d::System::setRho ( double  rho)
inline

Overwrite the mass density to a homogeneous quantity.

To use a non-homogeneous density use setMassMatrix().

Parameters
rhoMass density.

Definition at line 417 of file Generic2d.h.

◆ setT()

void FrictionQPotFEM::Generic2d::System::setT ( double  arg)
inline

Overwrite time.

Definition at line 899 of file Generic2d.h.

◆ setU()

template<class T >
void FrictionQPotFEM::Generic2d::System::setU ( const T &  u)
inline

Overwrite nodal displacements.

Internally, this updates the relevant forces using updated_u().

Parameters
u[nnode, ndim].

Definition at line 586 of file Generic2d.h.

◆ setV()

template<class T >
void FrictionQPotFEM::Generic2d::System::setV ( const T &  v)
inline

Overwrite nodal velocities.

Internally, this updates the relevant forces using updated_v().

Parameters
v[nnode, ndim].

Definition at line 601 of file Generic2d.h.

◆ Sig()

const array_type::tensor< double, 4 > & FrictionQPotFEM::Generic2d::System::Sig ( )
inline

Stress tensor of each integration point.

Returns
Integration point tensor (pointer). Shape: [nelem, nip, 2, 2].

Definition at line 1043 of file Generic2d.h.

◆ stiffness()

GooseFEM::MatrixPartitioned FrictionQPotFEM::Generic2d::System::stiffness ( ) const
inline

Stiffness tensor of each integration point.

Convenience function: assembles output from elastic() and plastic().

Returns
Integration point tensor (copy). Shape: [nelem, nip, 2, 2, 2, 2].

Definition at line 1087 of file Generic2d.h.

◆ t()

double FrictionQPotFEM::Generic2d::System::t ( ) const
inline

Current time.

Definition at line 891 of file Generic2d.h.

◆ timeStep()

void FrictionQPotFEM::Generic2d::System::timeStep ( )
inline

Make a time-step: apply velocity-Verlet integration.

Forces are computed where needed using: updated_u(), updated_v(), and computeInternalExternalResidualForce().

Definition at line 1541 of file Generic2d.h.

◆ timeSteps()

long FrictionQPotFEM::Generic2d::System::timeSteps ( size_t  n,
size_t  nmargin = 0 
)
inline

Make a number of time steps, see timeStep().

Parameters
nNumber of steps to make.
nmarginNumber of potentials to leave as margin.
  • nmargin > 0: stop if the yield-index of any box is nmargin from the end and return a negative number.
  • nmargin == 0: no bounds-check is performed (the last potential is assumed infinitely elastic to the right).
Returns
  • Number of steps: == n.
  • Negative number: if stopped because of a yield-index margin.

Definition at line 1601 of file Generic2d.h.

◆ timeStepsUntilEvent()

size_t FrictionQPotFEM::Generic2d::System::timeStepsUntilEvent ( size_t  nmargin = 0,
double  tol = 1e-5,
size_t  niter_tol = 20,
size_t  max_iter = 1e7 
)
inline

Perform a series of time-steps until the next plastic event, or equilibrium.

Parameters
nmarginNumber of potentials to have as margin initially.
tolRelative force tolerance for equilibrium. See System::residual for definition.
niter_tolEnforce the residual check for niter_tol consecutive increments.
max_iterMaximum number of iterations. Throws std::runtime_error otherwise.
Returns
  • Number of steps.
  • 0 if there was no plastic activity and the residual was reached.

Definition at line 1638 of file Generic2d.h.

◆ type()

virtual std::string FrictionQPotFEM::Generic2d::System::type ( ) const
inlinevirtual

◆ u()

const auto & FrictionQPotFEM::Generic2d::System::u ( ) const
inline

Nodal displacements.

Returns
Nodal displacements.

Definition at line 728 of file Generic2d.h.

◆ updated_u()

virtual void FrictionQPotFEM::Generic2d::System::updated_u ( )
inlinevirtual

Evaluate relevant forces when m_u is updated.

Under normal circumstances, this function is called automatically when needed.

Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System, and FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System.

Definition at line 626 of file Generic2d.h.

◆ updated_v()

void FrictionQPotFEM::Generic2d::System::updated_v ( )
inline

Evaluate relevant forces when m_v is updated.

Under normal circumstances, this function is called automatically when needed.

Definition at line 635 of file Generic2d.h.

◆ v()

const auto & FrictionQPotFEM::Generic2d::System::v ( ) const
inline

Nodal velocities.

Returns
Nodal velocities.

Definition at line 738 of file Generic2d.h.

◆ vector()

const GooseFEM::VectorPartitioned & FrictionQPotFEM::Generic2d::System::vector ( ) const
inline

GooseFEM vector definition.

Takes care of bookkeeping.

Returns
GooseFEM::VectorPartitioned (System::m_vector)

Definition at line 938 of file Generic2d.h.

Member Data Documentation

◆ m_a

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_a
protected

Nodal accelerations.

Definition at line 2114 of file Generic2d.h.

◆ m_a_n

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_a_n
protected

Nodal accelerations last time-step.

Definition at line 2116 of file Generic2d.h.

◆ m_alpha

double FrictionQPotFEM::Generic2d::System::m_alpha
protected

Background damping density (non-zero only if homogeneous).

Definition at line 2142 of file Generic2d.h.

◆ m_coor

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_coor
protected

Nodal coordinates, see coor().

Definition at line 2076 of file Generic2d.h.

◆ m_D

GooseFEM::MatrixDiagonal FrictionQPotFEM::Generic2d::System::m_D
protected

Damping matrix (diagonal).

Definition at line 2094 of file Generic2d.h.

◆ m_dt

double FrictionQPotFEM::Generic2d::System::m_dt
protected

Time-step.

Definition at line 2139 of file Generic2d.h.

◆ m_elas

GMatElastoPlasticQPot::Cartesian2d::Elastic<2> FrictionQPotFEM::Generic2d::System::m_elas
protected

Material for elastic el.

Definition at line 2095 of file Generic2d.h.

◆ m_elem_elas

array_type::tensor<size_t, 1> FrictionQPotFEM::Generic2d::System::m_elem_elas
protected

Elastic elements.

Definition at line 2085 of file Generic2d.h.

◆ m_elem_plas

array_type::tensor<size_t, 1> FrictionQPotFEM::Generic2d::System::m_elem_plas
protected

Plastic elements.

Definition at line 2086 of file Generic2d.h.

◆ m_Eps

array_type::tensor<double, 4> FrictionQPotFEM::Generic2d::System::m_Eps
protected

Quad-point tensor: strain.

Definition at line 2132 of file Generic2d.h.

◆ m_Epsdot_plas

array_type::tensor<double, 4> FrictionQPotFEM::Generic2d::System::m_Epsdot_plas
protected

Quad-point tensor: strain-rate for plastic el.

Definition at line 2134 of file Generic2d.h.

◆ m_eta

double FrictionQPotFEM::Generic2d::System::m_eta
protected

Damping at the interface.

Definition at line 2140 of file Generic2d.h.

◆ m_fdamp

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fdamp
protected

Nodal force from damping.

Definition at line 2126 of file Generic2d.h.

◆ m_fe

array_type::tensor<double, 3> FrictionQPotFEM::Generic2d::System::m_fe
protected

Element vector (used for forces).

Definition at line 2118 of file Generic2d.h.

◆ m_fe_elas

array_type::tensor<double, 3> FrictionQPotFEM::Generic2d::System::m_fe_elas
protected

m_fe for elastic element only

Definition at line 2120 of file Generic2d.h.

◆ m_fe_plas

array_type::tensor<double, 3> FrictionQPotFEM::Generic2d::System::m_fe_plas
protected

m_fe for plastic element only

Definition at line 2122 of file Generic2d.h.

◆ m_felas

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_felas
protected

Nodal force from elasticity of elastic elements.

Definition at line 2124 of file Generic2d.h.

◆ m_fext

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fext
protected

Nodal force: total external force (reaction force)

Definition at line 2130 of file Generic2d.h.

◆ m_fint

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fint
protected

Nodal force: total internal force.

Definition at line 2129 of file Generic2d.h.

◆ m_fmaterial

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fmaterial
protected

Nodal force, deriving from elasticity.

Definition at line 2123 of file Generic2d.h.

◆ m_fplas

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fplas
protected

Nodal force from elasticity of plastic elements.

Definition at line 2125 of file Generic2d.h.

◆ m_fres

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fres
protected

Nodal force: residual force.

Definition at line 2131 of file Generic2d.h.

◆ m_ftmp

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_ftmp
protected

Temporary for internal use.

Definition at line 2128 of file Generic2d.h.

◆ m_full_outdated

bool FrictionQPotFEM::Generic2d::System::m_full_outdated
protected

Keep track of the need to recompute fields on full geometry.

Definition at line 2145 of file Generic2d.h.

◆ m_fvisco

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_fvisco
protected

Nodal force from damping at the interface.

Definition at line 2127 of file Generic2d.h.

◆ m_inc

size_t FrictionQPotFEM::Generic2d::System::m_inc
protected

Current increment (current time = m_dt * m_inc).

Definition at line 2138 of file Generic2d.h.

◆ m_K_elas

GooseFEM::Matrix FrictionQPotFEM::Generic2d::System::m_K_elas
protected

Stiffness matrix for elastic elements only.

Definition at line 2135 of file Generic2d.h.

◆ m_M

GooseFEM::MatrixDiagonalPartitioned FrictionQPotFEM::Generic2d::System::m_M
protected

Mass matrix (diagonal).

Definition at line 2093 of file Generic2d.h.

◆ m_N

size_t FrictionQPotFEM::Generic2d::System::m_N
protected

Number of plastic elements, alias of m_nelem_plas.

Definition at line 2077 of file Generic2d.h.

◆ m_ndim

size_t FrictionQPotFEM::Generic2d::System::m_ndim
protected

Number of spatial dimensions.

Definition at line 2083 of file Generic2d.h.

◆ m_nelem

size_t FrictionQPotFEM::Generic2d::System::m_nelem
protected

Number of elements.

Definition at line 2078 of file Generic2d.h.

◆ m_nelem_elas

size_t FrictionQPotFEM::Generic2d::System::m_nelem_elas
protected

Number of elastic elements.

Definition at line 2079 of file Generic2d.h.

◆ m_nelem_plas

size_t FrictionQPotFEM::Generic2d::System::m_nelem_plas
protected

Number of plastic elements.

Definition at line 2080 of file Generic2d.h.

◆ m_nip

size_t FrictionQPotFEM::Generic2d::System::m_nip
protected

Number of integration points.

Definition at line 2084 of file Generic2d.h.

◆ m_nne

size_t FrictionQPotFEM::Generic2d::System::m_nne
protected

Number of nodes per element.

Definition at line 2081 of file Generic2d.h.

◆ m_nnode

size_t FrictionQPotFEM::Generic2d::System::m_nnode
protected

Number of nodes.

Definition at line 2082 of file Generic2d.h.

◆ m_pert_Epsd_plastic

array_type::tensor<double, 4> FrictionQPotFEM::Generic2d::System::m_pert_Epsd_plastic
protected

Strain deviator for m_pert_u.

Definition at line 2147 of file Generic2d.h.

◆ m_pert_u

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_pert_u
protected

See eventDriven_setDeltaU()

Definition at line 2146 of file Generic2d.h.

◆ m_plas

GMatElastoPlasticQPot::Cartesian2d::Cusp<2> FrictionQPotFEM::Generic2d::System::m_plas
protected

Material for plastic el.

Definition at line 2096 of file Generic2d.h.

◆ m_qs_inc_first

size_t FrictionQPotFEM::Generic2d::System::m_qs_inc_first = 0
protected

First increment with plastic activity during minimisation.

Definition at line 2136 of file Generic2d.h.

◆ m_qs_inc_last

size_t FrictionQPotFEM::Generic2d::System::m_qs_inc_last = 0
protected

Last increment with plastic activity during minimisation.

Definition at line 2137 of file Generic2d.h.

◆ m_quad

GooseFEM::Element::Quad4::Quadrature FrictionQPotFEM::Generic2d::System::m_quad
protected

Quadrature for all elements.

Definition at line 2087 of file Generic2d.h.

◆ m_quad_elas

GooseFEM::Element::Quad4::Quadrature FrictionQPotFEM::Generic2d::System::m_quad_elas
protected

m_quad for elastic elements only.

Definition at line 2088 of file Generic2d.h.

◆ m_quad_plas

GooseFEM::Element::Quad4::Quadrature FrictionQPotFEM::Generic2d::System::m_quad_plas
protected

m_quad for plastic elements only.

Definition at line 2089 of file Generic2d.h.

◆ m_rho

double FrictionQPotFEM::Generic2d::System::m_rho
protected

Mass density (non-zero only if homogeneous).

Definition at line 2141 of file Generic2d.h.

◆ m_set_D

bool FrictionQPotFEM::Generic2d::System::m_set_D
protected

Use m_D only if it is non-zero.

Definition at line 2143 of file Generic2d.h.

◆ m_set_visco

bool FrictionQPotFEM::Generic2d::System::m_set_visco
protected

Use m_eta only if it is non-zero.

Definition at line 2144 of file Generic2d.h.

◆ m_Sig

array_type::tensor<double, 4> FrictionQPotFEM::Generic2d::System::m_Sig
protected

Quad-point tensor: stress.

Definition at line 2133 of file Generic2d.h.

◆ m_u

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_u
protected

Nodal displacements.

Warning
To make sure that the right forces are computed at the right time, always call updated_u() after manually updating m_u (setU() automatically take care of this).

Definition at line 2104 of file Generic2d.h.

◆ m_ue

array_type::tensor<double, 3> FrictionQPotFEM::Generic2d::System::m_ue
protected

Element vector (used for displacements).

Definition at line 2117 of file Generic2d.h.

◆ m_ue_elas

array_type::tensor<double, 3> FrictionQPotFEM::Generic2d::System::m_ue_elas
protected

m_ue for elastic element only

Definition at line 2119 of file Generic2d.h.

◆ m_ue_plas

array_type::tensor<double, 3> FrictionQPotFEM::Generic2d::System::m_ue_plas
protected

m_ue for plastic element only

Definition at line 2121 of file Generic2d.h.

◆ m_v

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_v
protected

Nodal velocities.

Warning
To make sure that the right forces are computed at the right time, always call updated_v() after manually updating m_v (setV() automatically take care of this).

Definition at line 2112 of file Generic2d.h.

◆ m_v_n

array_type::tensor<double, 2> FrictionQPotFEM::Generic2d::System::m_v_n
protected

Nodal velocities last time-step.

Definition at line 2115 of file Generic2d.h.

◆ m_vector

GooseFEM::VectorPartitioned FrictionQPotFEM::Generic2d::System::m_vector
protected

Convert vectors between 'nodevec', 'elemvec', ....

Definition at line 2090 of file Generic2d.h.

◆ m_vector_elas

GooseFEM::VectorPartitioned FrictionQPotFEM::Generic2d::System::m_vector_elas
protected

m_vector for elastic elements only.

Definition at line 2091 of file Generic2d.h.

◆ m_vector_plas

GooseFEM::VectorPartitioned FrictionQPotFEM::Generic2d::System::m_vector_plas
protected

m_vector for plastic elements only.

Definition at line 2092 of file Generic2d.h.


The documentation for this class was generated from the following file: