FrictionQPotFEM 0.23.3
|
System with elastic elements and plastic elements (GMatElastoPlasticQPot). More...
#include <FrictionQPotFEM/Generic2d.h>
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::VectorPartitioned & | vector () const |
GooseFEM vector definition. More... | |
const GooseFEM::Element::Quad4::Quadrature & | quad () 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... | |
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.
|
inlinevirtual |
Definition at line 189 of file Generic2d.h.
|
inline |
Define the geometry, including boundary conditions and element sets.
C | Type of nodal coordinates, e.g. array_type::tensor<double, 2> |
E | Type of connectivity and DOFs, e.g. array_type::tensor<size_t, 2> |
L | Type of node/element lists, e.g. array_type::tensor<size_t, 1> |
coor | Nodal coordinates. |
conn | Connectivity. |
dofs | DOFs per node. |
iip | DOFs whose displacement is fixed. |
elastic_elem | Elastic elements. |
elastic_K | Bulk modulus per quad. point of each elastic element, see setElastic(). |
elastic_G | Shear modulus per quad. point of each elastic element, see setElastic(). |
plastic_elem | Plastic elements. |
plastic_K | Bulk modulus per quad. point of each plastic element, see Plastic(). |
plastic_G | Shear modulus per quad. point of each plastic element, see Plastic(). |
plastic_epsy | Yield strain per quad. point of each plastic element, see Plastic(). |
dt | Time step, set setDt(). |
rho | Mass density, see setMassMatrix(). |
alpha | Background damping density, see setDampingMatrix(). |
eta | Damping at the interface (homogeneous), see setEta(). |
Definition at line 215 of file Generic2d.h.
|
inline |
|
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.
delta_gamma | Strain to add (the shear component of deformation gradient is twice that). |
Definition at line 2041 of file Generic2d.h.
|
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).
delta_gamma | Strain to add (the shear component of deformation gradient is twice that). |
Definition at line 2061 of file Generic2d.h.
|
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.
Definition at line 507 of file Generic2d.h.
|
inline |
Modify the external force such that a shear stress is applied to the top boundary.
sigxy | The shear stress to be applied. |
Definition at line 789 of file Generic2d.h.
|
inlineprotected |
Compute m_Sig and m_Eps using m_u.
Definition at line 1134 of file Generic2d.h.
|
inlineprotected |
Update m_fmaterial based on the current displacement field m_u.
Definition at line 1179 of file Generic2d.h.
|
inlineprotectedvirtual |
Compute:
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.
|
inline |
Connectivity.
Shape: [System::m_nelem, System::m_nne].
Definition at line 698 of file Generic2d.h.
|
inline |
Nodal coordinates.
Shape: [System::m_nnode, System::m_ndim].
Definition at line 708 of file Generic2d.h.
|
inline |
Damping matrix, see setAlpha() and setDampingMatrix().
Note that there can be second source of damping, see setEta().
Definition at line 769 of file Generic2d.h.
|
inline |
|
inline |
Get time step.
Definition at line 566 of file Generic2d.h.
|
inline |
Integration point volume (of each integration point)
Definition at line 928 of file Generic2d.h.
|
inline |
Elastic material mode, used for all elastic elements.
Definition at line 959 of file Generic2d.h.
|
inline |
List of elastic elements.
Shape: [System::m_nelem_elas].
Definition at line 678 of file Generic2d.h.
|
inline |
Strain tensor of each integration point.
[nelem, nip, 2, 2]
. Definition at line 1054 of file Generic2d.h.
|
inline |
Symmetric gradient of the nodal accelerations of each integration point.
[nelem, nip, 2, 2]
. Definition at line 1076 of file Generic2d.h.
|
inline |
Strain-rate tensor (the symmetric gradient of the nodal velocities) of each integration point.
[nelem, nip, 2, 2]
. Definition at line 1066 of file Generic2d.h.
|
inline |
|
inline |
Get displacement perturbation used for event driven code, see eventDriven_setDeltaU().
Definition at line 1261 of file Generic2d.h.
|
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().
delta_u | Nodal displacement field. |
autoscale | Scale such that the perturbation is small enough. |
Definition at line 1234 of file Generic2d.h.
|
inlinevirtual |
Add event driven step for the boundary conditions that correspond to the displacement perturbation set in eventDriven_setDeltaU().
FRICTIONQPOTFEM_WIP
assertion is made.deps | Size of the local stain kick to apply. |
kick | If false , increment displacements to deps / 2 of yielding again. If true , increment displacements by a affine simple shear increment deps . |
direction | If +1 : apply shear to the right. If -1 applying shear to the left. |
yield_element | If 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. |
iterative | If true the step is iteratively searched. This is more costly by recommended, if the perturbation is non-analytical (and contains numerical errors). |
Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System, and FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System.
Definition at line 1300 of file Generic2d.h.
|
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).
[nnode, ndim]
. Definition at line 854 of file Generic2d.h.
|
inline |
External force.
[nnode, ndim]
. Definition at line 779 of file Generic2d.h.
|
inline |
Internal force.
[nnode, ndim]
. Definition at line 831 of file Generic2d.h.
|
inline |
Make a number of steps with the following protocol.
n | Number of steps to make. |
v | Nodal velocity [nnode, ndim] .s |
nmargin | Number of potentials to leave as margin.
|
Definition at line 1803 of file Generic2d.h.
|
inline |
Force deriving from elasticity.
[nnode, ndim]
. Definition at line 842 of file Generic2d.h.
|
inline |
Shear modulus per integration point.
Convenience function: assembles output from elastic() and plastic().
[nelem, nip]
. Definition at line 1012 of file Generic2d.h.
|
inline |
|
inline |
Check if elasticity is homogeneous.
true
is elasticity is homogeneous (false
otherwise). Definition at line 554 of file Generic2d.h.
|
inline |
Bulk modulus per integration point.
Convenience function: assembles output from elastic() and plastic().
[nelem, nip]
. Definition at line 980 of file Generic2d.h.
|
inline |
Mass matrix, see System::m_M.
Definition at line 758 of file Generic2d.h.
|
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)
nmargin | Number of potentials to leave as margin.
|
tol | Relative force tolerance for equilibrium. See System::residual for definition. |
niter_tol | Enforce the residual check for niter_tol consecutive increments. |
max_iter | Maximum number of time-steps. Throws std::runtime_error otherwise. |
time_activity | If true plastic activity is timed. After this function you can find:
|
max_iter_is_error | If true an error is thrown when the maximum number of time-steps is reached. If false the function simply returns max_iter . |
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
.Definition at line 1868 of file Generic2d.h.
|
inline |
Minimise energy while doing a high-frequency measurement of:
nodes | List of nodes to sample. |
top | Nodes of the top boundary. f_x = sum(fext[top, 0]) / (top.size - 1) (accounting for periodicity). |
t_step | Sample every t_step time steps. |
nmargin | Number of potentials to have as margin initially. |
tol | Relative force tolerance for equilibrium. See System::residual for definition. |
niter_tol | Enforce the residual check for niter_tol consecutive increments. |
max_iter | Maximum number of iterations. Throws std::runtime_error otherwise. |
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.
|
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).
A_truncate
and S_truncate
are defined on the first integration point.A_truncate | Truncate if A_truncate blocks have yielded at least once. |
S_truncate | Truncate 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. |
tol | Relative force tolerance for equilibrium. See System::residual for definition. |
niter_tol | Enforce the residual check for niter_tol consecutive increments. |
max_iter | Maximum 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.
idx_n | Reference potential index of the first integration point. |
Definition at line 1962 of file Generic2d.h.
|
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).
A_truncate
and S_truncate
are defined on the first integration point.A_truncate | Truncate if A_truncate blocks have yielded at least once. |
S_truncate | Truncate 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. |
tol | Relative force tolerance for equilibrium. See System::residual for definition. |
niter_tol | Enforce the residual check for niter_tol consecutive increments. |
max_iter | Maximum number of time-steps. Throws std::runtime_error otherwise. |
Definition at line 2023 of file Generic2d.h.
|
inlinevirtual |
Return the linear system size (in number of blocks).
Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System.
Definition at line 385 of file Generic2d.h.
|
inline |
Elastic material mode, used for all elastic elements.
Definition at line 969 of file Generic2d.h.
|
inline |
List of plastic elements.
Shape: [System::m_nelem_plas].
Definition at line 688 of file Generic2d.h.
|
inline |
Strain-rate tensor of integration points of plastic elements only, see System::plastic.
Definition at line 1119 of file Generic2d.h.
|
inline |
GooseFEM quadrature definition.
Takes case of interpolation, integration, and differentiation.
Definition at line 949 of file Generic2d.h.
|
inline |
Increment with the first plastic event.
This value is only relevant if time_activity = true
was used in the last call of minimise().
Definition at line 1936 of file Generic2d.h.
|
inline |
Increment with the last plastic event.
This value is only relevant if time_activity = true
was used in the last call of minimise().
Definition at line 1947 of file Generic2d.h.
|
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.
|
inline |
Recompute all forces.
Under normal circumstances, this function is called automatically when needed.
Definition at line 1215 of file Generic2d.h.
|
inline |
Norm of the relative residual force (the external force at the top/bottom boundaries is used for normalisation).
Definition at line 865 of file Generic2d.h.
|
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.
Definition at line 406 of file Generic2d.h.
|
inline |
Overwrite nodal accelerations.
a | [nnode, ndim] . |
Definition at line 615 of file Generic2d.h.
|
inline |
Overwrite background damping density (proportional to the velocity), To use a non-homogeneous density use setDampingMatrix().
alpha | Damping parameter. |
Definition at line 489 of file Generic2d.h.
|
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().
val_elem | Damping per element. |
Definition at line 521 of file Generic2d.h.
|
inline |
Overwrite time step.
Using for example in System::timeStep and System::minimise.
Definition at line 574 of file Generic2d.h.
|
inline |
Overwrite the value of the damping at the interface.
Note that you can specify either setEta() or setDampingMatrix() or both.
eta | Damping parameter |
Definition at line 461 of file Generic2d.h.
|
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.
fext | [nnode, ndim] . |
Definition at line 667 of file Generic2d.h.
|
inlinevirtual |
Set increment.
arg | size_t. |
Reimplemented in FrictionQPotFEM::UniformSingleLayerThermal2d::System.
Definition at line 918 of file Generic2d.h.
|
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().
T | e.g. array_type::tensor<double, 1> . |
val_elem | Density per element. |
Definition at line 431 of file Generic2d.h.
|
inline |
Overwrite the mass density to a homogeneous quantity.
To use a non-homogeneous density use setMassMatrix().
rho | Mass density. |
Definition at line 417 of file Generic2d.h.
|
inline |
Overwrite time.
Definition at line 899 of file Generic2d.h.
|
inline |
Overwrite nodal displacements.
Internally, this updates the relevant forces using updated_u().
u | [nnode, ndim] . |
Definition at line 586 of file Generic2d.h.
|
inline |
Overwrite nodal velocities.
Internally, this updates the relevant forces using updated_v().
v | [nnode, ndim] . |
Definition at line 601 of file Generic2d.h.
|
inline |
Stress tensor of each integration point.
[nelem, nip, 2, 2]
. Definition at line 1043 of file Generic2d.h.
|
inline |
Stiffness tensor of each integration point.
Convenience function: assembles output from elastic() and plastic().
[nelem, nip, 2, 2, 2, 2]
. Definition at line 1087 of file Generic2d.h.
|
inline |
Current time.
Definition at line 891 of file Generic2d.h.
|
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.
|
inline |
Make a number of time steps, see timeStep().
n | Number of steps to make. |
nmargin | Number of potentials to leave as margin.
|
== n
.Definition at line 1601 of file Generic2d.h.
|
inline |
Perform a series of time-steps until the next plastic event, or equilibrium.
nmargin | Number of potentials to have as margin initially. |
tol | Relative force tolerance for equilibrium. See System::residual for definition. |
niter_tol | Enforce the residual check for niter_tol consecutive increments. |
max_iter | Maximum number of iterations. Throws std::runtime_error otherwise. |
0
if there was no plastic activity and the residual was reached. Definition at line 1638 of file Generic2d.h.
|
inlinevirtual |
Return the type of system.
Reimplemented in FrictionQPotFEM::UniformMultiLayerIndividualDrive2d::System, FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System, FrictionQPotFEM::UniformSingleLayer2d::System, and FrictionQPotFEM::UniformSingleLayerThermal2d::System.
Definition at line 394 of file Generic2d.h.
|
inline |
|
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.
|
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.
|
inline |
|
inline |
|
protected |
Nodal accelerations.
Definition at line 2114 of file Generic2d.h.
|
protected |
Nodal accelerations last time-step.
Definition at line 2116 of file Generic2d.h.
|
protected |
Background damping density (non-zero only if homogeneous).
Definition at line 2142 of file Generic2d.h.
|
protected |
Nodal coordinates, see coor().
Definition at line 2076 of file Generic2d.h.
|
protected |
Damping matrix (diagonal).
Definition at line 2094 of file Generic2d.h.
|
protected |
Time-step.
Definition at line 2139 of file Generic2d.h.
|
protected |
Material for elastic el.
Definition at line 2095 of file Generic2d.h.
|
protected |
Elastic elements.
Definition at line 2085 of file Generic2d.h.
|
protected |
Plastic elements.
Definition at line 2086 of file Generic2d.h.
|
protected |
Quad-point tensor: strain.
Definition at line 2132 of file Generic2d.h.
|
protected |
Quad-point tensor: strain-rate for plastic el.
Definition at line 2134 of file Generic2d.h.
|
protected |
Damping at the interface.
Definition at line 2140 of file Generic2d.h.
|
protected |
Nodal force from damping.
Definition at line 2126 of file Generic2d.h.
|
protected |
Element vector (used for forces).
Definition at line 2118 of file Generic2d.h.
|
protected |
m_fe for elastic element only
Definition at line 2120 of file Generic2d.h.
|
protected |
m_fe for plastic element only
Definition at line 2122 of file Generic2d.h.
|
protected |
Nodal force from elasticity of elastic elements.
Definition at line 2124 of file Generic2d.h.
|
protected |
Nodal force: total external force (reaction force)
Definition at line 2130 of file Generic2d.h.
|
protected |
Nodal force: total internal force.
Definition at line 2129 of file Generic2d.h.
|
protected |
Nodal force, deriving from elasticity.
Definition at line 2123 of file Generic2d.h.
|
protected |
Nodal force from elasticity of plastic elements.
Definition at line 2125 of file Generic2d.h.
|
protected |
Nodal force: residual force.
Definition at line 2131 of file Generic2d.h.
|
protected |
Temporary for internal use.
Definition at line 2128 of file Generic2d.h.
|
protected |
Keep track of the need to recompute fields on full geometry.
Definition at line 2145 of file Generic2d.h.
|
protected |
Nodal force from damping at the interface.
Definition at line 2127 of file Generic2d.h.
|
protected |
Current increment (current time = m_dt * m_inc).
Definition at line 2138 of file Generic2d.h.
|
protected |
Stiffness matrix for elastic elements only.
Definition at line 2135 of file Generic2d.h.
|
protected |
Mass matrix (diagonal).
Definition at line 2093 of file Generic2d.h.
|
protected |
Number of plastic elements, alias of m_nelem_plas.
Definition at line 2077 of file Generic2d.h.
|
protected |
Number of spatial dimensions.
Definition at line 2083 of file Generic2d.h.
|
protected |
Number of elements.
Definition at line 2078 of file Generic2d.h.
|
protected |
Number of elastic elements.
Definition at line 2079 of file Generic2d.h.
|
protected |
Number of plastic elements.
Definition at line 2080 of file Generic2d.h.
|
protected |
Number of integration points.
Definition at line 2084 of file Generic2d.h.
|
protected |
Number of nodes per element.
Definition at line 2081 of file Generic2d.h.
|
protected |
Number of nodes.
Definition at line 2082 of file Generic2d.h.
|
protected |
Strain deviator for m_pert_u.
Definition at line 2147 of file Generic2d.h.
|
protected |
Definition at line 2146 of file Generic2d.h.
|
protected |
Material for plastic el.
Definition at line 2096 of file Generic2d.h.
|
protected |
First increment with plastic activity during minimisation.
Definition at line 2136 of file Generic2d.h.
|
protected |
Last increment with plastic activity during minimisation.
Definition at line 2137 of file Generic2d.h.
|
protected |
Quadrature for all elements.
Definition at line 2087 of file Generic2d.h.
|
protected |
m_quad for elastic elements only.
Definition at line 2088 of file Generic2d.h.
|
protected |
m_quad for plastic elements only.
Definition at line 2089 of file Generic2d.h.
|
protected |
Mass density (non-zero only if homogeneous).
Definition at line 2141 of file Generic2d.h.
|
protected |
Use m_D only if it is non-zero.
Definition at line 2143 of file Generic2d.h.
|
protected |
Use m_eta only if it is non-zero.
Definition at line 2144 of file Generic2d.h.
|
protected |
Quad-point tensor: stress.
Definition at line 2133 of file Generic2d.h.
|
protected |
Nodal displacements.
Definition at line 2104 of file Generic2d.h.
|
protected |
Element vector (used for displacements).
Definition at line 2117 of file Generic2d.h.
|
protected |
m_ue for elastic element only
Definition at line 2119 of file Generic2d.h.
|
protected |
m_ue for plastic element only
Definition at line 2121 of file Generic2d.h.
|
protected |
Nodal velocities.
Definition at line 2112 of file Generic2d.h.
|
protected |
Nodal velocities last time-step.
Definition at line 2115 of file Generic2d.h.
|
protected |
Convert vectors between 'nodevec', 'elemvec', ....
Definition at line 2090 of file Generic2d.h.
|
protected |
m_vector for elastic elements only.
Definition at line 2091 of file Generic2d.h.
|
protected |
m_vector for plastic elements only.
Definition at line 2092 of file Generic2d.h.