v0.14.0
SmallStrainPlasticity.hpp

Problem formulation

Theory

Introduction

In the following document, we define functions, \(\psi(\varepsilon^e,\alpha)\), \(f(\sigma,\mathbf{q})\) and \(\Psi(\sigma,\mathbf{q})\) which are free Helmholtz energy, yield function and flow potential, respectively. Note that free Helmholtz energy is defined by stains and kinematic internal state variables, whereas yield function and flow potential is defined by fluxes (stress tensor and fluxes work conjugated to internal kinematic state variables). In this implantation, we assume that all functions are convex. Also, is assumed that energy potential and flow potential is positive and takes zero value at zero strain and zero internal variable states. It can be shown that for such defined problem, the first law of thermodynamics and non-negative dissipation of energy is automatically satisfied, see [23] [66]. Moreover, in the following implementation is assumed that all functions are smooth enough that appropriate derivatives could be calculated.

In the following implementation, the backward-Euler difference scheme is applied to solve set nonlinear system of equations with unilateral constrains. Admissible stress state is enforced by Lagrange multiplier method and (Khun-Tucker) loading/unloading conditions, see [66]. That lead to well known Closet point projection algorithm.

Following [66] can be shown that plastic loading unloading conditions can be exclusively characterised in terms of trial state

\[ \varepsilon^{e,\textrm{trial}}_{n+1} := \varepsilon_{n+1} - \varepsilon^p_n \]

\[ \alpha^\textrm{trial}_{n+1} := \alpha_n \]

\[ \sigma^\textrm{trial}_{n+1} := \frac{\partial \psi}{\partial \varepsilon^e}(\varepsilon^{e,\textrm{trial}}_{n+1},\alpha_n) \]

\[ \mathbf{q}^\textrm{trial}_{n+1} := \frac{\partial \psi}{\partial \alpha}(\varepsilon^{e,\textrm{trial}}_{n+1},\alpha_n) \]

\[ f^\textrm{trial}_{n+1} := f(\sigma^\textrm{trial}_{n+1},\mathbf{q}_n) \]

then loading/unloading conditions are given

\[ \begin{split} f^\textrm{trial}_{n+1} < 0,\quad\textrm{elastic step}\,\to\,\Delta\gamma_{n+1} = 0 \\ f^\textrm{trial}_{n+1} \geq 0,\quad\textrm{plastic step}\,\to\,\Delta\gamma_{n+1} > 0 \end{split} \]

If first condition is true \(f^\textrm{trial}_{n+1} < 0\), trial stress is solution at integration point. If second condition is true, plastic strains and values of internal state variables are calculated using Closet point projection algorithm.

Closet point projection algorithm

Additive strain decomposition

\[ \varepsilon^e_{n+1} := \varepsilon_{n+1} - \varepsilon^p_{n+1} \]

Constitutive relations

\[ \sigma_{n+1} := \frac{\partial \psi}{\partial \varepsilon^e}(\varepsilon^e_{n+1},\alpha_{n+1}) \]

\[ \mathbf{q}_{n+1} := \frac{\partial \psi}{\partial \alpha}(\varepsilon^e_{n+1},\alpha_{n+1}) \]

Constrains

\[ f_{n+1}:=f(\sigma_{n+1},\mathbf{q}_{n+1}) \]

Evolution functions

Flow vector

\[ \mathbf{n}_{n+1}:=\frac{\partial \Psi}{\partial \sigma}(\sigma_{n+1},\mathbf{q}_{n+1}) \]

Hardening law:

\[ \mathbf{h}_{n+1}:=-\frac{\partial \Psi}{\partial \mathbf{q}}(\sigma_{n+1},\mathbf{q}_{n+1}) \]

Residual

\[ \mathbf{R}_{n+1}:= \left\{ \begin{array}{c} -\varepsilon^p_{n+1}+\varepsilon^p_n \\ -\alpha_{n+1}+\alpha_n \\ f_{n+1} \end{array} \right\} + \Delta\gamma_{n+1} \left\{ \begin{array}{c} \mathbf{n}_{n+1} \\ \mathbf{h}_{n+1} \\ 0 \end{array} \right\} \]

\[ \left[ \begin{array}{ccc} -1+\gamma \left( \frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p} + \frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial\alpha\varepsilon^p} \right) & \gamma \left( \frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2} + \frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \alpha} \right) & \mathbf{n} \\ \gamma \left( \frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p} + \frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial\alpha\partial\varepsilon^p} \right) & -1+\gamma \left( \frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2} + \frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \alpha} \right) & \mathbf{h}\\ \frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p} + \frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha\partial\varepsilon^p} & \frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2} + \frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial\alpha} & 0 \end{array} \right]_n \left\{ \begin{array}{c} \delta \varepsilon^p \\ \delta \alpha \\ \delta \gamma \end{array} \right\}_{n+1} = \mathbf{R}_n \]

where \(\varepsilon^p_{n+1}=\varepsilon^p_n + \delta\varepsilon^p_{n+1}\), \(\alpha^p_{n+1}=\alpha_n + \delta\alpha_{n+1}\) and \(\Delta\gamma_{n+1} = \Delta\gamma_n+\delta\gamma_{n+1}\). Usually it is assumed free energy split, i.e.

\[ \psi(\varepsilon^e,\alpha) = \psi^{\varepsilon}(\varepsilon^e) + \psi^\alpha(\alpha) \]

in that case linearized residual takes the form

\[ \left[ \begin{array}{ccc} -1+\Delta\gamma \frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p} & \Delta\gamma \frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2} & \mathbf{n} \\ \Delta\gamma \frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p} & -1+\Delta\gamma \frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2} & \mathbf{h}\\ \frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p} & \frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2} & 0 \end{array} \right]_n \left\{ \begin{array}{c} \delta \varepsilon^p \\ \delta \alpha \\ \delta \gamma \end{array} \right\}_{n+1} = -\mathbf{R}_n \]

\[ \chi_{n+1} := \left\{ \begin{array}{c} \delta \varepsilon^p \\ \delta \alpha \\ \delta \gamma \end{array} \right\}_{n+1} = - \mathbf{A}^{-1} \mathbf{R}_n \]

Algorithmic tangent matrix

\[ \delta \sigma := \mathbf{C}^e(\delta \varepsilon - \delta \varepsilon^p) \]

where \(\mathbf{C}^e = \frac{\partial^2 \psi^e}{\partial {\varepsilon^e}^2}\) and \(\delta \varepsilon = \textrm{sym}[ \nabla \delta \mathbf{u}]\)

\[ \left\{ \begin{array}{c} \mathbf{E}^p \\ \mathbf{L}^p \\ \mathbf{Y}^p \end{array} \right\} \delta \varepsilon = - \mathbf{A}^{-1} \left\{ \begin{array}{c} \Delta\gamma \frac{\partial \mathbf{n}}{\partial \sigma}\mathbf{C}^e \\ \Delta\gamma \frac{\partial \mathbf{h}}{\partial \sigma}\mathbf{C}^e \\ \frac{\partial f}{\partial \sigma}\mathbf{C}^e \end{array} \right\} \delta \varepsilon \]

\[ \delta \varepsilon^p = \mathbf{E}^p \delta \varepsilon \]

\[ \delta \alpha = \mathbf{L}^p \delta \varepsilon \]

\[ \delta \gamma = \mathbf{Y}^p \delta \varepsilon \]

\[ \delta \sigma = \mathbf{C}^e(\mathbf{1} - \mathbf{E}^p)\delta \varepsilon \]

\[ \delta \sigma = \mathbf{C}^{ep}\delta \varepsilon \]

where \(\mathbf{C}^{ep} = \mathbf{C}^e-\mathbf{C}^e\mathbf{E}^p\)

For more details and easy access to all source files related to this see group Small Strain Plasticity.

/** \file SmallStrainPlasticity.hpp
* \ingroup small_strain_plasticity
* \brief Operators and data structures for small strain plasticity
* \example SmallStrainPlasticity.hpp
\section problem_formulation Problem formulation
\subsection theory_small_strain Theory
\subsubsection introduction Introduction
In the following document, we define functions,
\f$\psi(\varepsilon^e,\alpha)\f$, \f$f(\sigma,\mathbf{q})\f$ and
\f$\Psi(\sigma,\mathbf{q})\f$ which are free Helmholtz energy, yield function
and flow potential, respectively. Note that free Helmholtz energy is defined
by stains and kinematic internal state variables, whereas yield function and
flow potential is defined by fluxes (stress tensor and fluxes work conjugated
to internal kinematic state variables). In this implantation, we assume that all
functions are convex. Also, is assumed that energy potential and flow
potential is positive and takes zero value at zero strain and zero internal
variable states. It can be shown that for such defined problem, the first law
of thermodynamics and non-negative dissipation of energy is automatically
satisfied, see \cite de2011computational \cite simo2006computational.
Moreover, in the following implementation is assumed that all functions are
smooth enough that appropriate derivatives could be calculated.
In the following implementation, the backward-Euler difference scheme is
applied to solve set nonlinear system of equations with unilateral constrains.
Admissible stress state is enforced by Lagrange multiplier method and
(Khun-Tucker) loading/unloading conditions, see \cite simo2006computational.
That lead to well known \em Closet \em point \em projection algorithm.
Following \cite simo2006computational can be shown that plastic loading
unloading conditions can be exclusively characterised in terms of trial
state
\f[
\varepsilon^{e,\textrm{trial}}_{n+1} := \varepsilon_{n+1} - \varepsilon^p_n
\f]
\f[
\alpha^\textrm{trial}_{n+1} := \alpha_n
\f]
\f[
\sigma^\textrm{trial}_{n+1} := \frac{\partial \psi}{\partial \varepsilon^e}(\varepsilon^{e,\textrm{trial}}_{n+1},\alpha_n)
\f]
\f[
\mathbf{q}^\textrm{trial}_{n+1} := \frac{\partial \psi}{\partial \alpha}(\varepsilon^{e,\textrm{trial}}_{n+1},\alpha_n)
\f]
\f[
f^\textrm{trial}_{n+1} := f(\sigma^\textrm{trial}_{n+1},\mathbf{q}_n)
\f]
then loading/unloading conditions are given
\f[
\begin{split}
f^\textrm{trial}_{n+1} < 0,\quad\textrm{elastic step}\,\to\,\Delta\gamma_{n+1} = 0 \\
f^\textrm{trial}_{n+1} \geq 0,\quad\textrm{plastic step}\,\to\,\Delta\gamma_{n+1} > 0
\end{split}
\f]
If first condition is true \f$f^\textrm{trial}_{n+1} < 0\f$, trial stress is solution
at integration point. If second condition is true, plastic strains and values
of internal state variables are calculated using \em Closet \em point \em
projection algorithm.
\subsection cloaset_point_projection Closet point projection algorithm
\subsubsection elastic_strain Additive strain decomposition
\f[
\varepsilon^e_{n+1} := \varepsilon_{n+1} - \varepsilon^p_{n+1}
\f]
\subsubsection constitutive_equation Constitutive relations
\f[
\sigma_{n+1} := \frac{\partial \psi}{\partial \varepsilon^e}(\varepsilon^e_{n+1},\alpha_{n+1})
\f]
\f[
\mathbf{q}_{n+1} := \frac{\partial \psi}{\partial \alpha}(\varepsilon^e_{n+1},\alpha_{n+1})
\f]
\subsubsection constrains Constrains
\f[
f_{n+1}:=f(\sigma_{n+1},\mathbf{q}_{n+1})
\f]
\subsubsection evolution Evolution functions
Flow vector
\f[
\mathbf{n}_{n+1}:=\frac{\partial \Psi}{\partial \sigma}(\sigma_{n+1},\mathbf{q}_{n+1})
\f]
Hardening law:
\f[
\mathbf{h}_{n+1}:=-\frac{\partial \Psi}{\partial \mathbf{q}}(\sigma_{n+1},\mathbf{q}_{n+1})
\f]
\subsubsection residual Residual
\f[
\mathbf{R}_{n+1}:=
\left\{
\begin{array}{c}
-\varepsilon^p_{n+1}+\varepsilon^p_n \\
-\alpha_{n+1}+\alpha_n \\
f_{n+1}
\end{array}
\right\}
+
\Delta\gamma_{n+1}
\left\{
\begin{array}{c}
\mathbf{n}_{n+1} \\
\mathbf{h}_{n+1} \\
0
\end{array}
\right\}
\f]
\f[
\left[
\begin{array}{ccc}
-1+\gamma
\left(
\frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
+
\frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial\alpha\varepsilon^p}
\right)
&
\gamma
\left(
\frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
+
\frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \alpha}
\right)
&
\mathbf{n}
\\
\gamma
\left(
\frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
+
\frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial\alpha\partial\varepsilon^p}
\right)
&
-1+\gamma
\left(
\frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
+
\frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \alpha}
\right)
&
\mathbf{h}\\
\frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
+
\frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha\partial\varepsilon^p}
&
\frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
+
\frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial\alpha}
& 0
\end{array}
\right]_n
\left\{
\begin{array}{c}
\delta \varepsilon^p \\
\delta \alpha \\
\delta \gamma
\end{array}
\right\}_{n+1}
=
\mathbf{R}_n
\f]
where
\f$\varepsilon^p_{n+1}=\varepsilon^p_n + \delta\varepsilon^p_{n+1}\f$,
\f$\alpha^p_{n+1}=\alpha_n + \delta\alpha_{n+1}\f$
and
\f$\Delta\gamma_{n+1} = \Delta\gamma_n+\delta\gamma_{n+1}\f$.
Usually it is assumed free energy split, i.e.
\f[
\psi(\varepsilon^e,\alpha) = \psi^{\varepsilon}(\varepsilon^e) + \psi^\alpha(\alpha)
\f]
in that case linearized residual takes the form
\f[
\left[
\begin{array}{ccc}
-1+\Delta\gamma
\frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
&
\Delta\gamma
\frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
&
\mathbf{n}
\\
\Delta\gamma
\frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
&
-1+\Delta\gamma
\frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
&
\mathbf{h}\\
\frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
&
\frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
& 0
\end{array}
\right]_n
\left\{
\begin{array}{c}
\delta \varepsilon^p \\
\delta \alpha \\
\delta \gamma
\end{array}
\right\}_{n+1}
=
-\mathbf{R}_n
\f]
\f[
\chi_{n+1} :=
\left\{
\begin{array}{c}
\delta \varepsilon^p \\
\delta \alpha \\
\delta \gamma
\end{array}
\right\}_{n+1}
=
-
\mathbf{A}^{-1}
\mathbf{R}_n
\f]
\subsubsection small_strain_algebraic_tangent Algorithmic tangent matrix
\f[
\delta \sigma := \mathbf{C}^e(\delta \varepsilon - \delta \varepsilon^p)
\f]
where \f$\mathbf{C}^e = \frac{\partial^2 \psi^e}{\partial {\varepsilon^e}^2}\f$ and
\f$\delta \varepsilon = \textrm{sym}[ \nabla \delta \mathbf{u}]\f$
\f[
\left\{
\begin{array}{c}
\mathbf{E}^p \\
\mathbf{L}^p \\
\mathbf{Y}^p
\end{array}
\right\}
\delta \varepsilon
=
-
\mathbf{A}^{-1}
\left\{
\begin{array}{c}
\Delta\gamma \frac{\partial \mathbf{n}}{\partial \sigma}\mathbf{C}^e \\
\Delta\gamma \frac{\partial \mathbf{h}}{\partial \sigma}\mathbf{C}^e \\
\frac{\partial f}{\partial \sigma}\mathbf{C}^e
\end{array}
\right\}
\delta \varepsilon
\f]
\f[
\delta \varepsilon^p = \mathbf{E}^p \delta \varepsilon
\f]
\f[
\delta \alpha = \mathbf{L}^p \delta \varepsilon
\f]
\f[
\delta \gamma = \mathbf{Y}^p \delta \varepsilon
\f]
\f[
\delta \sigma = \mathbf{C}^e(\mathbf{1} - \mathbf{E}^p)\delta \varepsilon
\f]
\f[
\delta \sigma = \mathbf{C}^{ep}\delta \varepsilon
\f]
where \f$\mathbf{C}^{ep} = \mathbf{C}^e-\mathbf{C}^e\mathbf{E}^p\f$
For more details and easy access to all source files related to this see group
\ref small_strain_plasticity.
*/
/*
* This file is part of MoFEM.
* MoFEM is free software: you can redistribute it and/or modify it under
* the terms of the GNU Lesser General Public License as published by the
* Free Software Foundation, either version 3 of the License, or (at your
* option) any later version.
*
* MoFEM is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
#ifndef __SMALL_STRAIN_PLASTICITY_HPP
#define __SMALL_STRAIN_PLASTICITY_HPP
#ifndef WITH_ADOL_C
#error "MoFEM need to be compiled with ADOL-C"
#endif
PetscErrorCode SmallStrainPlasticityfRes(SNES snes,Vec chi,Vec r,void *ctx);
PetscErrorCode SmallStrainPlasticityfJac(SNES,Vec,Mat,Mat,void *ctx);
/** \brief Small strain finite element implementation
* \ingroup small_strain_plasticity
*/
/// \brief definition of volume element
struct MyVolumeFE: public MoFEM::VolumeElementForcesAndSourcesCore {
int addToRule;
addToRule(0) { }
/** \brief it is used to calculate nb. of Gauss integration points
*
* for more details pleas look
* Reference:
*
* Albert Nijenhuis, Herbert Wilf,
* Combinatorial Algorithms for Computers and Calculators,
* Second Edition,
* Academic Press, 1978,
* ISBN: 0-12-519260-6,
* LC: QA164.N54.
*
* More details about algorithm
* http://people.sc.fsu.edu/~jburkardt/cpp_src/gm_rule/gm_rule.html
**/
int getRule(int order) {
// cerr << "order " << order << endl;
return 2*(order-1)+addToRule;
}
};
MyVolumeFE feRhs; ///< calculate right hand side for tetrahedral elements
MyVolumeFE feLhs; //< calculate left hand side for tetrahedral elements
MyVolumeFE feUpdate; //< update history variables
mField(m_field),
feRhs(m_field),
feLhs(m_field),
feUpdate(m_field) {
}
PetscErrorCode createInternalVariablesTag();
/** \brief common data used by volume elements
* \ingroup nonlinear_elastic_elem
*/
struct CommonData {
map<string,vector<VectorDouble> > dataAtGaussPts;
map<string,vector<MatrixDouble> > gradAtGaussPts;
vector<VectorDouble> sTress;
vector<MatrixDouble> materialTangentOperator;
vector<VectorDouble> plasticStrain;
vector<VectorDouble> internalVariables;
vector<VectorDouble> internalFluxes;
vector<double> deltaGamma;
bool bBar;
bBar(true) {
}
};
vector<VectorDouble > &valuesAtGaussPts;
vector<MatrixDouble > &gradientAtGaussPts;
const EntityType zeroAtType;
vector<VectorDouble > &values_at_gauss_pts,
vector<MatrixDouble > &gardient_at_gauss_pts
);
/** \brief operator calculating field data and gradients
*/
PetscErrorCode doWork(
int side,EntityType type,DataForcesAndSourcesCore::EntData &data
);
};
struct OpGetCommonDataAtGaussPts: public OpGetDataAtGaussPts {
OpGetCommonDataAtGaussPts(const string field_name,CommonData &common_data);
};
string field_name,
CommonData &common_data
);
PetscErrorCode doWork(
int side,EntityType type,DataForcesAndSourcesCore::EntData &data
);
};
struct MakeB {
PetscErrorCode makeB(const MatrixAdaptor &diffN,MatrixDouble &B);
PetscErrorCode addVolumetricB(const MatrixAdaptor &diffN,MatrixDouble &volB,double alpha);
};
string field_name,
CommonData &common_data
);
PetscErrorCode doWork(
int side,EntityType type,DataForcesAndSourcesCore::EntData &data
);
};
string field_name,
CommonData &common_data
);
PetscErrorCode doWork(
int row_side,int col_side,
EntityType row_type,EntityType col_type,
);
};
/** \brief Closest point projection algorithm
* \ingroup small_strain_plasticity
*/
struct ClosestPointProjection {
double deltaGamma;
double tOl;
int gG;
ublas::symmetric_matrix<double,ublas::lower> C,D;
ublas::symmetric_matrix<double,ublas::lower> partial2HSigma;
ublas::symmetric_matrix<double,ublas::lower> partial2HFlux;
vector<int> tAgs;
virtual PetscErrorCode setActiveVariablesW();
virtual PetscErrorCode setActiveVariablesYH();
double w,y,h;
ublas::vector<adouble> a_sTress,a_internalFluxes;
virtual PetscErrorCode rEcordW();
virtual PetscErrorCode rEcordY();
virtual PetscErrorCode rEcordH();
virtual PetscErrorCode pLayW();
virtual PetscErrorCode pLayW_NoHessian();
virtual PetscErrorCode pLayY();
virtual PetscErrorCode pLayY_NoGradient();
virtual PetscErrorCode pLayH();
Mat A;
ublas::matrix<double,ublas::column_major> dataA;
virtual PetscErrorCode createMatAVecR();
virtual PetscErrorCode destroyMatAVecR();
virtual PetscErrorCode evaluatePotentials();
virtual PetscErrorCode cAlculateR(Vec R);
virtual PetscErrorCode cAlculateA();
friend PetscErrorCode SmallStrainPlasticityfRes(SNES,Vec,Vec,void *ctx);
friend PetscErrorCode SmallStrainPlasticityfJac(SNES,Vec,Mat,Mat,void *ctx);
SNES sNes;
KSP kSp;
PC pC;
SNESLineSearch lineSearch;
PetscErrorCode snesCreate();
PetscErrorCode snesDestroy();
PetscErrorCode solveColasetProjection();
PetscErrorCode consistentTangent();
virtual PetscErrorCode freeHemholtzEnergy() {
PetscFunctionBegin;
SETERRQ(
PETSC_COMM_SELF,
"Not implemented (overload class and implement this function)"
);
PetscFunctionReturn(0);
}
virtual PetscErrorCode yieldFunction() {
PetscFunctionBegin;
SETERRQ(
PETSC_COMM_SELF,
"Not implemented (overload class and implement this function)"
);
PetscFunctionReturn(0);
}
virtual PetscErrorCode flowPotential() {
PetscFunctionBegin;
SETERRQ(
PETSC_COMM_SELF,
"Not implemented (overload class and implement this function)"
);
PetscFunctionReturn(0);
}
};
ClosestPointProjection &cP;
bool initCp;
MoFEM::Interface &m_field,
string field_name,
CommonData &common_data,
ClosestPointProjection &cp,
bool is_linear = true
):
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator(field_name,UserDataOperator::OPROW),
mField(m_field),
commonData(common_data),
cP(cp),
initCp(false),
isLinear(is_linear),
}
PetscErrorCode getTags();
PetscErrorCode setTagsData(
const EntityHandle tet,const int nb_gauss_pts,const int nb_internal_variables
);
PetscErrorCode doWork(
int side,EntityType type,DataForcesAndSourcesCore::EntData &data
);
};
};
#endif //__SMALL_STRAIN_PLASTICITY_HPP
/***************************************************************************//**
* \defgroup small_strain_plasticity Small Strain Plasticity
* \ingroup user_modules
******************************************************************************/
/***************************************************************************//**
* \defgroup user_modules User modules
******************************************************************************/
SmallStrainPlasticity::ClosestPointProjection::SmallStrainPlasticityfRes
friend PetscErrorCode SmallStrainPlasticityfRes(SNES, Vec, Vec, void *ctx)
Definition: SmallStrainPlasticity.cpp:1154
SmallStrainPlasticity::MyVolumeFE::getRule
int getRule(int order)
it is used to calculate nb. of Gauss integration points
Definition: SmallStrainPlasticity.hpp:353
SmallStrainPlasticity::OpAssembleRhs::commonData
CommonData & commonData
Definition: SmallStrainPlasticity.hpp:440
SmallStrainPlasticity::ClosestPointProjection::Cp
MatrixDouble Cp
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::ClosestPointProjection::evaluatePotentials
virtual PetscErrorCode evaluatePotentials()
Definition: SmallStrainPlasticity.cpp:945
SmallStrainPlasticity::OpCalculateStress::cP
ClosestPointProjection & cP
Definition: SmallStrainPlasticity.hpp:607
SmallStrainPlasticity::OpGetDataAtGaussPts::gradientAtGaussPts
vector< MatrixDouble > & gradientAtGaussPts
Definition: SmallStrainPlasticity.hpp:398
SmallStrainPlasticity::ClosestPointProjection::plasticStrain0
VectorDouble plasticStrain0
Definition: SmallStrainPlasticity.hpp:488
SmallStrainPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: SmallStrainPlasticity.hpp:517
SmallStrainPlasticity::OpAssembleLhs::colVolB
MatrixDouble colVolB
Definition: SmallStrainPlasticity.hpp:461
SmallStrainPlasticity::OpGetDataAtGaussPts::zeroAtType
const EntityType zeroAtType
Definition: SmallStrainPlasticity.hpp:399
SmallStrainPlasticity::ClosestPointProjection::lineSearch
SNESLineSearch lineSearch
Definition: SmallStrainPlasticity.hpp:559
SmallStrainPlasticity::feUpdate
MyVolumeFE feUpdate
Definition: SmallStrainPlasticity.hpp:362
EntityHandle
SmallStrainPlasticity::ClosestPointProjection::internalVariables0
VectorDouble internalVariables0
Definition: SmallStrainPlasticity.hpp:487
SmallStrainPlasticity::ClosestPointProjection::yieldFunction
virtual PetscErrorCode yieldFunction()
Definition: SmallStrainPlasticity.hpp:581
SmallStrainPlasticity::CommonData::sTress
vector< VectorDouble > sTress
Definition: SmallStrainPlasticity.hpp:379
SmallStrainPlasticity::OpAssembleLhs::transK
MatrixDouble transK
Definition: SmallStrainPlasticity.hpp:459
SmallStrainPlasticity::OpAssembleLhs::colB
MatrixDouble colB
Definition: SmallStrainPlasticity.hpp:460
SmallStrainPlasticity::ClosestPointProjection::setActiveVariablesW
virtual PetscErrorCode setActiveVariablesW()
Definition: SmallStrainPlasticity.cpp:665
SmallStrainPlasticity::OpAssembleRhs::doWork
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Definition: SmallStrainPlasticity.cpp:227
SmallStrainPlasticity::ClosestPointProjection::partialHFlux
VectorAdaptor partialHFlux
Definition: SmallStrainPlasticity.hpp:501
SmallStrainPlasticity::ClosestPointProjection::pLayH
virtual PetscErrorCode pLayH()
Definition: SmallStrainPlasticity.cpp:870
SmallStrainPlasticity::ClosestPointProjection::destroyMatAVecR
virtual PetscErrorCode destroyMatAVecR()
Definition: SmallStrainPlasticity.cpp:936
SmallStrainPlasticity::ClosestPointProjection::plasticStrain
VectorDouble plasticStrain
Definition: SmallStrainPlasticity.hpp:489
SmallStrainPlasticity::OpGetCommonDataAtGaussPts::OpGetCommonDataAtGaussPts
OpGetCommonDataAtGaussPts(const string field_name, CommonData &common_data)
Definition: SmallStrainPlasticity.cpp:110
SmallStrainPlasticity::OpAssembleRhs::nF
VectorDouble nF
Definition: SmallStrainPlasticity.hpp:441
SmallStrainPlasticity::ClosestPointProjection::rEcordH
virtual PetscErrorCode rEcordH()
Definition: SmallStrainPlasticity.cpp:740
SmallStrainPlasticity::MyVolumeFE::addToRule
int addToRule
Definition: SmallStrainPlasticity.hpp:333
SmallStrainPlasticity::CommonData::plasticStrainSize
int plasticStrainSize
Definition: SmallStrainPlasticity.hpp:387
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
SmallStrainPlasticity::ClosestPointProjection::R
Vec R
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::sTrain
VectorDouble sTrain
Definition: SmallStrainPlasticity.hpp:485
SmallStrainPlasticity::ClosestPointProjection::A
Mat A
Definition: SmallStrainPlasticity.hpp:540
SmallStrainPlasticity::CommonData::deltaGamma
vector< double > deltaGamma
Definition: SmallStrainPlasticity.hpp:384
SmallStrainPlasticity::SmallStrainPlasticity
SmallStrainPlasticity(MoFEM::Interface &m_field)
Definition: SmallStrainPlasticity.hpp:364
SmallStrainPlasticity::ClosestPointProjection::solveColasetProjection
PetscErrorCode solveColasetProjection()
Definition: SmallStrainPlasticity.cpp:1050
SmallStrainPlasticity::ClosestPointProjection::partialHSigma
VectorAdaptor partialHSigma
Definition: SmallStrainPlasticity.hpp:500
SmallStrainPlasticity::ClosestPointProjection::rEcordY
virtual PetscErrorCode rEcordY()
Definition: SmallStrainPlasticity.cpp:718
SmallStrainPlasticity::ClosestPointProjection::Ep
MatrixDouble Ep
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::ClosestPointProjection::hessianH
MatrixDouble hessianH
Definition: SmallStrainPlasticity.hpp:537
SmallStrainPlasticity::ClosestPointProjection::h
double h
Definition: SmallStrainPlasticity.hpp:516
SmallStrainPlasticity::ClosestPointProjection::tAgs
vector< int > tAgs
Definition: SmallStrainPlasticity.hpp:506
SmallStrainPlasticity::ClosestPointProjection::hessianW
MatrixDouble hessianW
Definition: SmallStrainPlasticity.hpp:528
SmallStrainPlasticity::ClosestPointProjection::D
ublas::symmetric_matrix< double, ublas::lower > D
Definition: SmallStrainPlasticity.hpp:496
SmallStrainPlasticity::OpUpdate::OpUpdate
OpUpdate(string field_name, CommonData &common_data)
Definition: SmallStrainPlasticity.cpp:118
SmallStrainPlasticity::ClosestPointProjection::dataA
ublas::matrix< double, ublas::column_major > dataA
Definition: SmallStrainPlasticity.hpp:541
SmallStrainPlasticity::ClosestPointProjection::partialYSigma
VectorAdaptor partialYSigma
Definition: SmallStrainPlasticity.hpp:498
SmallStrainPlasticity::ClosestPointProjection::pLayY_NoGradient
virtual PetscErrorCode pLayY_NoGradient()
Definition: SmallStrainPlasticity.cpp:860
SmallStrainPlasticity::ClosestPointProjection::partial2HFlux
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
Definition: SmallStrainPlasticity.hpp:503
SmallStrainPlasticity::ClosestPointProjection::snesCreate
PetscErrorCode snesCreate()
Definition: SmallStrainPlasticity.cpp:1027
SmallStrainPlasticity::ClosestPointProjection::gG
int gG
Definition: SmallStrainPlasticity.hpp:492
sdf.r
int r
Definition: sdf.py:8
SmallStrainPlasticity::ClosestPointProjection::createMatAVecR
virtual PetscErrorCode createMatAVecR()
Definition: SmallStrainPlasticity.cpp:924
order
constexpr int order
Definition: dg_projection.cpp:18
SmallStrainPlasticity::OpAssembleLhs::K
MatrixDouble K
Definition: SmallStrainPlasticity.hpp:459
SmallStrainPlasticity::ClosestPointProjection::SmallStrainPlasticityfJac
friend PetscErrorCode SmallStrainPlasticityfJac(SNES, Vec, Mat, Mat, void *ctx)
Definition: SmallStrainPlasticity.cpp:1207
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
SmallStrainPlasticity::OpAssembleLhs::doWork
PetscErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
Definition: SmallStrainPlasticity.cpp:311
SmallStrainPlasticity::CommonData::internalVariablesPtr
double * internalVariablesPtr
Definition: SmallStrainPlasticity.hpp:386
SmallStrainPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: SmallStrainPlasticity.hpp:518
MoFEM::VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore
VolumeElementForcesAndSourcesCore(Interface &m_field, const EntityType type=MBTET)
Definition: VolumeElementForcesAndSourcesCore.cpp:9
SmallStrainPlasticity::CommonData::internalVariablesSize
int internalVariablesSize
Definition: SmallStrainPlasticity.hpp:388
SmallStrainPlasticity::ClosestPointProjection::internalFluxes
VectorAdaptor internalFluxes
Definition: SmallStrainPlasticity.hpp:495
SmallStrainPlasticity::CommonData::plasticStrain
vector< VectorDouble > plasticStrain
Definition: SmallStrainPlasticity.hpp:381
SmallStrainPlasticity::ClosestPointProjection::sNes
SNES sNes
Definition: SmallStrainPlasticity.hpp:556
SmallStrainPlasticity::ClosestPointProjection::gradientY
VectorDouble gradientY
Definition: SmallStrainPlasticity.hpp:532
SmallStrainPlasticity::OpCalculateStress::thInternalVariables
Tag thInternalVariables
Definition: SmallStrainPlasticity.hpp:631
SmallStrainPlasticity::ClosestPointProjection::tOl
double tOl
Definition: SmallStrainPlasticity.hpp:491
SmallStrainPlasticity::OpCalculateStress::mField
MoFEM::Interface & mField
Definition: SmallStrainPlasticity.hpp:605
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
SmallStrainPlasticity::ClosestPointProjection::partial2HSigmaFlux
MatrixDouble partial2HSigmaFlux
Definition: SmallStrainPlasticity.hpp:504
SmallStrainPlasticity::ClosestPointProjection::partialWStrainPlasticStrain
MatrixDouble partialWStrainPlasticStrain
Definition: SmallStrainPlasticity.hpp:497
SmallStrainPlasticity::ClosestPointProjection::snesDestroy
PetscErrorCode snesDestroy()
Definition: SmallStrainPlasticity.cpp:1044
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
SmallStrainPlasticity::CommonData::gradAtGaussPts
map< string, vector< MatrixDouble > > gradAtGaussPts
Definition: SmallStrainPlasticity.hpp:378
SmallStrainPlasticity::OpAssembleLhs::rowVolB
MatrixDouble rowVolB
Definition: SmallStrainPlasticity.hpp:461
convert.type
type
Definition: convert.py:64
SmallStrainPlasticity::CommonData::internalFluxes
vector< VectorDouble > internalFluxes
Definition: SmallStrainPlasticity.hpp:383
SmallStrainPlasticity::feLhs
MyVolumeFE feLhs
Definition: SmallStrainPlasticity.hpp:361
SmallStrainPlasticityfJac
PetscErrorCode SmallStrainPlasticityfJac(SNES, Vec, Mat, Mat, void *ctx)
Definition: SmallStrainPlasticity.cpp:1207
SmallStrainPlasticity::ClosestPointProjection::Lp
MatrixDouble Lp
Definition: SmallStrainPlasticity.hpp:566
SmallStrainPlasticity::OpCalculateStress::hessianWCalculated
bool hessianWCalculated
Definition: SmallStrainPlasticity.hpp:610
SmallStrainPlasticity::ClosestPointProjection::Chi
Vec Chi
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::y
double y
Definition: SmallStrainPlasticity.hpp:516
SmallStrainPlasticity::ClosestPointProjection::partialYFlux
VectorAdaptor partialYFlux
Definition: SmallStrainPlasticity.hpp:499
SmallStrainPlasticity::OpGetDataAtGaussPts::OpGetDataAtGaussPts
OpGetDataAtGaussPts(const string field_name, vector< VectorDouble > &values_at_gauss_pts, vector< MatrixDouble > &gardient_at_gauss_pts)
Definition: SmallStrainPlasticity.cpp:48
SmallStrainPlasticity::OpAssembleLhs::rowB
MatrixDouble rowB
Definition: SmallStrainPlasticity.hpp:460
SmallStrainPlasticity::OpCalculateStress::thPlasticStrain
Tag thPlasticStrain
Definition: SmallStrainPlasticity.hpp:630
SmallStrainPlasticity::commonData
CommonData commonData
Definition: SmallStrainPlasticity.hpp:394
SmallStrainPlasticity::MyVolumeFE::MyVolumeFE
MyVolumeFE(MoFEM::Interface &m_field)
Definition: SmallStrainPlasticity.hpp:334
SmallStrainPlasticity::ClosestPointProjection::activeVariablesYH
VectorDouble activeVariablesYH
Definition: SmallStrainPlasticity.hpp:513
SmallStrainPlasticity::OpCalculateStress::commonData
CommonData & commonData
Definition: SmallStrainPlasticity.hpp:606
SmallStrainPlasticity::OpAssembleLhs::OpAssembleLhs
OpAssembleLhs(string field_name, CommonData &common_data)
Definition: SmallStrainPlasticity.cpp:302
SmallStrainPlasticity::ClosestPointProjection::activeVariablesW
VectorDouble activeVariablesW
Definition: SmallStrainPlasticity.hpp:510
SmallStrainPlasticity::CommonData::dataAtGaussPts
map< string, vector< VectorDouble > > dataAtGaussPts
Definition: SmallStrainPlasticity.hpp:377
SmallStrainPlasticity::CommonData::materialTangentOperator
vector< MatrixDouble > materialTangentOperator
Definition: SmallStrainPlasticity.hpp:380
SmallStrainPlasticity::ClosestPointProjection::rEcordW
virtual PetscErrorCode rEcordW()
Definition: SmallStrainPlasticity.cpp:692
SmallStrainPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: SmallStrainPlasticity.hpp:517
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
SmallStrainPlasticity
Small strain finite element implementation.
Definition: SmallStrainPlasticity.hpp:320
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
SmallStrainPlasticity::ClosestPointProjection::partial2HSigma
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
Definition: SmallStrainPlasticity.hpp:502
SmallStrainPlasticity::ClosestPointProjection::freeHemholtzEnergy
virtual PetscErrorCode freeHemholtzEnergy()
Definition: SmallStrainPlasticity.hpp:571
SmallStrainPlasticity::OpGetDataAtGaussPts::valuesAtGaussPts
vector< VectorDouble > & valuesAtGaussPts
Definition: SmallStrainPlasticity.hpp:397
SmallStrainPlasticity::MakeB::addVolumetricB
PetscErrorCode addVolumetricB(const MatrixAdaptor &diffN, MatrixDouble &volB, double alpha)
Definition: SmallStrainPlasticity.cpp:200
SmallStrainPlasticity::thPlasticStrain
Tag thPlasticStrain
Definition: SmallStrainPlasticity.hpp:327
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
SmallStrainPlasticity::OpAssembleRhs::OpAssembleRhs
OpAssembleRhs(string field_name, CommonData &common_data)
Definition: SmallStrainPlasticity.cpp:219
SmallStrainPlasticity::ClosestPointProjection::sTress
VectorAdaptor sTress
Definition: SmallStrainPlasticity.hpp:494
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
SmallStrainPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: SmallStrainPlasticity.hpp:518
SmallStrainPlasticity::ClosestPointProjection::w
double w
Definition: SmallStrainPlasticity.hpp:516
SmallStrainPlasticity::ClosestPointProjection::dChi
Vec dChi
Definition: SmallStrainPlasticity.hpp:542
SmallStrainPlasticity::ClosestPointProjection::internalVariables
VectorDouble internalVariables
Definition: SmallStrainPlasticity.hpp:486
SmallStrainPlasticity::CommonData::plasticStrainPtr
double * plasticStrainPtr
Definition: SmallStrainPlasticity.hpp:385
SmallStrainPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: SmallStrainPlasticity.hpp:519
adouble
SmallStrainPlasticity::OpCalculateStress::getTags
PetscErrorCode getTags()
Definition: SmallStrainPlasticity.cpp:432
SmallStrainPlasticity::ClosestPointProjection::pLayW
virtual PetscErrorCode pLayW()
Definition: SmallStrainPlasticity.cpp:762
SmallStrainPlasticity::OpCalculateStress::initCp
bool initCp
Definition: SmallStrainPlasticity.hpp:609
SmallStrainPlasticity::ClosestPointProjection::flowPotential
virtual PetscErrorCode flowPotential()
Definition: SmallStrainPlasticity.hpp:591
SmallStrainPlasticity::ClosestPointProjection::ClosestPointProjection
ClosestPointProjection()
Definition: SmallStrainPlasticity.hpp:508
SmallStrainPlasticity::ClosestPointProjection::pLayW_NoHessian
virtual PetscErrorCode pLayW_NoHessian()
Definition: SmallStrainPlasticity.cpp:813
SmallStrainPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: SmallStrainPlasticity.hpp:519
CommonData
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:22
SmallStrainPlasticity::createInternalVariablesTag
PetscErrorCode createInternalVariablesTag()
Definition: SmallStrainPlasticity.cpp:34
SmallStrainPlasticity::OpAssembleRhs::volB
MatrixDouble volB
Definition: SmallStrainPlasticity.hpp:442
SmallStrainPlasticity::ClosestPointProjection::kSp
KSP kSp
Definition: SmallStrainPlasticity.hpp:557
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
SmallStrainPlasticity::OpAssembleRhs::B
MatrixDouble B
Definition: SmallStrainPlasticity.hpp:442
SmallStrainPlasticity::ClosestPointProjection::cAlculateR
virtual PetscErrorCode cAlculateR(Vec R)
Definition: SmallStrainPlasticity.cpp:962
SmallStrainPlasticity::OpUpdate::commonData
CommonData & commonData
Definition: SmallStrainPlasticity.hpp:420
SmallStrainPlasticity::mField
MoFEM::Interface & mField
Definition: SmallStrainPlasticity.hpp:325
SmallStrainPlasticity::thInternalVariables
Tag thInternalVariables
Definition: SmallStrainPlasticity.hpp:328
SmallStrainPlasticity::OpGetDataAtGaussPts::doWork
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
operator calculating field data and gradients
Definition: SmallStrainPlasticity.cpp:58
SmallStrainPlasticity::ClosestPointProjection::setActiveVariablesYH
virtual PetscErrorCode setActiveVariablesYH()
Definition: SmallStrainPlasticity.cpp:680
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
SmallStrainPlasticity::OpCalculateStress::setTagsData
PetscErrorCode setTagsData(const EntityHandle tet, const int nb_gauss_pts, const int nb_internal_variables)
Definition: SmallStrainPlasticity.cpp:444
OpCalculateStress
Definition: HookeElement.hpp:153
SmallStrainPlasticity::OpCalculateStress::isLinear
bool isLinear
Definition: SmallStrainPlasticity.hpp:610
SmallStrainPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: SmallStrainPlasticity.hpp:517
MoFEM::Types::MatrixAdaptor
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
SmallStrainPlasticity::CommonData::bBar
bool bBar
Definition: SmallStrainPlasticity.hpp:389
SmallStrainPlasticity::ClosestPointProjection::consistentTangent
PetscErrorCode consistentTangent()
Definition: SmallStrainPlasticity.cpp:1095
SmallStrainPlasticity::ClosestPointProjection::deltaGamma
double deltaGamma
Definition: SmallStrainPlasticity.hpp:490
SmallStrainPlasticity::feRhs
MyVolumeFE feRhs
calculate right hand side for tetrahedral elements
Definition: SmallStrainPlasticity.hpp:360
SmallStrainPlasticity::ClosestPointProjection::gradientW
VectorDouble gradientW
Definition: SmallStrainPlasticity.hpp:527
SmallStrainPlasticity::ClosestPointProjection::pLayY
virtual PetscErrorCode pLayY()
Definition: SmallStrainPlasticity.cpp:837
SmallStrainPlasticity::ClosestPointProjection::pC
PC pC
Definition: SmallStrainPlasticity.hpp:558
SmallStrainPlasticity::MakeB::makeB
PetscErrorCode makeB(const MatrixAdaptor &diffN, MatrixDouble &B)
Definition: SmallStrainPlasticity.cpp:172
SmallStrainPlasticity::ClosestPointProjection::gradientH
VectorDouble gradientH
Definition: SmallStrainPlasticity.hpp:536
SmallStrainPlasticity::ClosestPointProjection::cAlculateA
virtual PetscErrorCode cAlculateA()
Definition: SmallStrainPlasticity.cpp:978
SmallStrainPlasticity::OpCalculateStress::OpCalculateStress
OpCalculateStress(MoFEM::Interface &m_field, string field_name, CommonData &common_data, ClosestPointProjection &cp, bool is_linear=true)
Definition: SmallStrainPlasticity.hpp:614
SmallStrainPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: SmallStrainPlasticity.hpp:519
SmallStrainPlasticity::OpAssembleLhs::CB
MatrixDouble CB
Definition: SmallStrainPlasticity.hpp:460
SmallStrainPlasticity::CommonData::internalVariables
vector< VectorDouble > internalVariables
Definition: SmallStrainPlasticity.hpp:382
SmallStrainPlasticity::ClosestPointProjection::C
ublas::symmetric_matrix< double, ublas::lower > C
Definition: SmallStrainPlasticity.hpp:496
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
SmallStrainPlasticity::OpAssembleLhs::commonData
CommonData & commonData
Definition: SmallStrainPlasticity.hpp:458
SmallStrainPlasticity::OpUpdate::doWork
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Definition: SmallStrainPlasticity.cpp:126
SmallStrainPlasticity::CommonData::CommonData
CommonData()
Definition: SmallStrainPlasticity.hpp:390
SmallStrainPlasticityfRes
PetscErrorCode SmallStrainPlasticityfRes(SNES snes, Vec chi, Vec r, void *ctx)
Definition: SmallStrainPlasticity.cpp:1154
SmallStrainPlasticity::OpCalculateStress::doWork
PetscErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Definition: SmallStrainPlasticity.cpp:489
SmallStrainPlasticity::ClosestPointProjection::Cep
MatrixDouble Cep
Definition: SmallStrainPlasticity.hpp:566
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567
SmallStrainPlasticity::ClosestPointProjection::Yp
VectorDouble Yp
Definition: SmallStrainPlasticity.hpp:567