FrictionQPotFEM 0.23.3
|
System that comprises several layers (elastic or plastic). More...
#include <FrictionQPotFEM/UniformMultiLayerIndividualDrive2d.h>
Public Member Functions | |
System (const array_type::tensor< double, 2 > &coor, const array_type::tensor< size_t, 2 > &conn, const array_type::tensor< size_t, 2 > &dofs, const array_type::tensor< size_t, 1 > &iip, const std::vector< array_type::tensor< size_t, 1 > > &elem, const std::vector< array_type::tensor< size_t, 1 > > &node, const array_type::tensor< bool, 1 > &layer_is_plastic, const array_type::tensor< double, 2 > &elastic_K, const array_type::tensor< double, 2 > &elastic_G, const array_type::tensor< double, 2 > &plastic_K, const array_type::tensor< double, 2 > &plastic_G, const array_type::tensor< double, 3 > &plastic_epsy, double dt, double rho, double alpha, double eta, const array_type::tensor< bool, 2 > &drive_is_active, double k_drive) | |
Define basic geometry. More... | |
size_t | N () const override |
Return the linear system size (in number of blocks). More... | |
std::string | type () const override |
Return the type of system. More... | |
size_t | nlayer () const |
Return number of layers. More... | |
array_type::tensor< size_t, 1 > | layerElements (size_t i) const |
Return the elements belonging to the i-th layer. More... | |
array_type::tensor< size_t, 1 > | layerNodes (size_t i) const |
Return the nodes belonging to the i-th layer. More... | |
const array_type::tensor< bool, 1 > & | layerIsPlastic () const |
Return if a layer is elastic (false ) or plastic (true ). More... | |
void | layerSetDriveStiffness (double k, bool symmetric=true) |
Set the stiffness of the springs connecting the average displacement of a layer ("ubar") to its set target value. More... | |
template<class S , class T > | |
double | initEventDriven (const S &ubar, const T &active) |
Initialise the event driven protocol by applying a perturbation to the loading springs and computing and storing the linear, purely elastic, response. More... | |
template<class S , class T , class U > | |
double | initEventDriven (const S &ubar, const T &active, const U &delta_u) |
Restore perturbation used from event driven protocol. More... | |
const array_type::tensor< double, 2 > & | eventDriven_deltaUbar () const |
Get target average position perturbation used for event driven code. More... | |
const array_type::tensor< bool, 2 > & | eventDriven_targetActive () const |
Get if the target average position is prescribed in the event driven code. More... | |
double | eventDrivenStep (double deps, bool kick, int direction=+1, bool yield_element=false, bool fallback=false) override |
Add event driven step for the boundary conditions that correspond to the displacement perturbation set in eventDriven_setDeltaU(). More... | |
template<class T > | |
void | layerSetTargetActive (const T &active) |
Turn on (or off) springs connecting the average displacement of a layer ("ubar") to its set target value. More... | |
array_type::tensor< double, 2 > | layerUbar () |
List the average displacement of each layer. More... | |
const array_type::tensor< double, 2 > & | layerTargetUbar () const |
List the target average displacement per layer. More... | |
const array_type::tensor< bool, 2 > & | layerTargetActive () const |
List if the driving spring is activate. More... | |
template<class T > | |
void | layerSetTargetUbar (const T &ubar) |
Overwrite the target average displacement per layer. More... | |
template<class S , class T > | |
void | layerSetUbar (const S &ubar, const T &prescribe) |
Move the layer such that the average displacement is exactly equal to its input value. More... | |
template<class T > | |
array_type::tensor< double, 2 > | layerTargetUbar_affineSimpleShear (double delta_gamma, const T &height) const |
Get target average displacements that correspond to affine simple shear (with the bottom fixed). More... | |
const array_type::tensor< double, 2 > & | fdrive () const |
Nodal force induced by the driving springs. More... | |
array_type::tensor< double, 2 > | layerFdrive () const |
Force of each of the driving springs. More... | |
![]() | |
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 | computeLayerUbarActive () |
Compute the average displacement of all layers with an active driving spring. More... | |
void | computeForceFromTargetUbar () |
Compute force deriving from the activate springs between the average displacement of the layer and its target value. More... | |
void | updated_u () override |
Evaluate relevant forces when m_u is updated. More... | |
void | computeInternalExternalResidualForce () override |
Compute: More... | |
![]() | |
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 | |
size_t | m_N |
Linear system size. More... | |
size_t | m_n_layer |
Number of layers. More... | |
std::vector< array_type::tensor< size_t, 1 > > | m_layer_node |
Nodes per layer. More... | |
std::vector< array_type::tensor< size_t, 1 > > | m_layer_elem |
Elements per layer. More... | |
array_type::tensor< bool, 1 > | m_layer_is_plastic |
Per layer true if the layer is plastic. More... | |
array_type::tensor< size_t, 1 > | m_slice_index |
Per layer the index in m_slice_plas or m_slice_elas. More... | |
array_type::tensor< size_t, 1 > | m_slice_plas |
How to slice elastic_elem(): start & end index. More... | |
array_type::tensor< size_t, 1 > | m_slice_elas |
How to slice plastic_elem(): start & end index. More... | |
bool | m_drive_spring_symmetric |
If false the drive spring buckles in compression. More... | |
double | m_drive_k |
Stiffness of the drive control frame. More... | |
array_type::tensor< bool, 2 > | m_layerdrive_active |
See prescribe in layerSetTargetUbar() More... | |
array_type::tensor< double, 2 > | m_layerdrive_targetubar |
Per layer, the prescribed average position. More... | |
array_type::tensor< double, 2 > | m_layer_ubar |
See computeLayerUbarActive(). More... | |
array_type::tensor< double, 2 > | m_fdrive |
Force related to driving frame. More... | |
array_type::tensor< double, 2 > | m_layer_dV1 |
volume per layer (as vector, same for all dimensions) More... | |
array_type::tensor< double, 2 > | m_dV |
copy of m_quad.dV() More... | |
array_type::tensor< double, 3 > | m_uq |
qvector More... | |
array_type::tensor< double, 1 > | m_ud |
dofval More... | |
array_type::tensor< bool, 2 > | m_pert_layerdrive_active |
Event driven: applied lever setting. More... | |
array_type::tensor< double, 2 > | m_pert_layerdrive_targetubar |
Event driven: applied lever setting. More... | |
![]() | |
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 that comprises several layers (elastic or plastic).
The average displacement of each layer can be coupled to a prescribed target value using a linear spring (one spring per spatial dimension):
Terminology:
ubar
: the average displacement per layer [nlayer, 2], see layerUbar() and layerSetUbar().target_ubar
: the target average displacement per layer [nlayer, 2], see layerTargetUbar() and layerSetTargetUbar().target_active
: the average displacement per layer/DOF is only enforced if the spring is active, see layerTargetActive() and layerSetTargetActive(). Definition at line 57 of file UniformMultiLayerIndividualDrive2d.h.
|
inlinevirtual |
Reimplemented from FrictionQPotFEM::Generic2d::System.
Definition at line 62 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Define basic geometry.
coor | Nodal coordinates. |
conn | Connectivity. |
dofs | DOFs per node. |
iip | DOFs whose displacement is fixed. |
elem | Elements per layer. |
node | Nodes per layer. |
layer_is_plastic | Per layer set if elastic (= 0) or plastic (= 1). |
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_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(). |
drive_is_active | Per layer (de)activate drive for each DOF, see layerSetTargetActive(). |
k_drive | Stiffness of driving spring, see layerSetDriveStiffness(). |
Definition at line 86 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineprotected |
Compute force deriving from the activate springs between the average displacement of the layer and its target value.
The force is applied as a force density.
Internal rule: computeLayerUbarActive() is called before this function, if the displacement changed since the last time the average was computed.
Definition at line 721 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineoverrideprotectedvirtual |
Compute:
Internal rule: all relevant forces are expected to be updated before this function is called.
Reimplemented from FrictionQPotFEM::Generic2d::System.
Definition at line 765 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineprotected |
Compute the average displacement of all layers with an active driving spring.
Definition at line 690 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Get target average position perturbation used for event driven code.
Definition at line 475 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Get if the target average position is prescribed in the event driven code.
Definition at line 484 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineoverridevirtual |
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 from FrictionQPotFEM::Generic2d::System.
Reimplemented in FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System.
Definition at line 489 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Nodal force induced by the driving springs.
The only non-zero contribution comes from:
Definition at line 658 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Initialise the event driven protocol by applying a perturbation to the loading springs and computing and storing the linear, purely elastic, response.
The system can thereafter be moved forward to the next event. Note that this function itself does not change the system in any way, it just stores the relevant perturbations.
S | array_type::tensor<double, 2> |
ubar | The perturbation in the target average position of each layer [nlayer, 2]. |
T | array_type::tensor<bool, 2> |
active | For each layer and each degree-of-freedom specify if the spring is active (true ) or not (false ) [nlayer, 2]. |
Definition at line 388 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Restore perturbation used from event driven protocol.
ubar | See eventDriven_deltaUbar(). |
active | See eventDriven_targetActive(). |
delta_u | See eventDriven_deltaU(). |
Definition at line 460 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Return the elements belonging to the i-th layer.
i | Index of the layer. |
Definition at line 325 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Force of each of the driving springs.
The only non-zero contribution comes from:
Definition at line 671 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Return if a layer is elastic (false
) or plastic (true
).
Definition at line 346 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Return the nodes belonging to the i-th layer.
i | Index of the layer. |
Definition at line 336 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Set the stiffness of the springs connecting the average displacement of a layer ("ubar") to its set target value.
Note that the stiffness of all springs is taken the same.
k | The stiffness. |
symmetric | If true the spring is a normal spring. If false the spring has no stiffness under compression. |
Definition at line 361 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Turn on (or off) springs connecting the average displacement of a layer ("ubar") to its set target value.
T | e.g. array_type::tensor<bool, 2> |
active | For each layer and each degree-of-freedom specify if the spring is active (true ) or not (false ) [nlayer, 2]. |
Definition at line 512 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Overwrite the target average displacement per layer.
Only layers that have an active driving spring will feel a force (if its average displacement is different from the target displacement), see layerSetTargetActive().
T | e.g. array_type::tensor<double, 2> |
ubar | The target average position of each layer [nlayer, 2]. |
Definition at line 582 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Move the layer such that the average displacement is exactly equal to its input value.
S | e.g. array_type::tensor<double, 2> |
T | e.g. array_type::tensor<bool, 2> |
ubar | The target average position of each layer [nlayer, 2]. |
prescribe | Per layers/degree-of-freedom specify if its average is modified [nlayer, 2]. Note that this not modify which of the driving springs is active or not, that can only be changed using layerSetTargetActive(). |
Definition at line 604 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
List if the driving spring is activate.
Definition at line 567 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
List the target average displacement per layer.
Definition at line 558 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Get target average displacements that correspond to affine simple shear (with the bottom fixed).
In particular \( \bar{u}_x^i = 2 \Delta \gamma h_i \) with \( \bar{u}_x^i \) the \( x \)-component of the target average displacement of layer \( i \).
T | e.g. array_type::tensor<double, 1> |
delta_gamma | Affine strain to add. |
height | The height \( h_i \) of the loading frame of each layer [nlayer]. |
Definition at line 637 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
List the average displacement of each layer.
Requires to recompute the average displacements (as they are normally only computed on the driven DOFs).
Definition at line 527 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineoverridevirtual |
Return the linear system size (in number of blocks).
Reimplemented from FrictionQPotFEM::Generic2d::System.
Definition at line 301 of file UniformMultiLayerIndividualDrive2d.h.
|
inline |
Return number of layers.
Definition at line 315 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineoverridevirtual |
Return the type of system.
Reimplemented from FrictionQPotFEM::Generic2d::System.
Reimplemented in FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System.
Definition at line 306 of file UniformMultiLayerIndividualDrive2d.h.
|
inlineoverrideprotectedvirtual |
Evaluate relevant forces when m_u is updated.
Reimplemented from FrictionQPotFEM::Generic2d::System.
Reimplemented in FrictionQPotFEM::UniformMultiLayerLeverDrive2d::System.
Definition at line 749 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Stiffness of the drive control frame.
Definition at line 784 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
If false the drive spring buckles in compression.
Definition at line 783 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
copy of m_quad.dV()
Definition at line 792 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Force related to driving frame.
Definition at line 789 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
volume per layer (as vector, same for all dimensions)
Definition at line 791 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Elements per layer.
Definition at line 776 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Per layer true
if the layer is plastic.
Definition at line 777 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Nodes per layer.
Definition at line 775 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Definition at line 788 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
See prescribe
in layerSetTargetUbar()
Definition at line 785 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Per layer, the prescribed average position.
Definition at line 787 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Linear system size.
Definition at line 773 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Number of layers.
Definition at line 774 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Event driven: applied lever setting.
Definition at line 796 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Event driven: applied lever setting.
Definition at line 798 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
How to slice plastic_elem(): start & end index.
Definition at line 781 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
Per layer the index in m_slice_plas or m_slice_elas.
Definition at line 779 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
How to slice elastic_elem(): start & end index.
Definition at line 780 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
dofval
Definition at line 794 of file UniformMultiLayerIndividualDrive2d.h.
|
protected |
qvector
Definition at line 793 of file UniformMultiLayerIndividualDrive2d.h.