v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
PlasticOps::OpCalculatePlasticity Struct Reference

#include <users_modules/tutorials/adv-0/src/PlasticOps.hpp>

Inheritance diagram for PlasticOps::OpCalculatePlasticity:
[legend]
Collaboration diagram for PlasticOps::OpCalculatePlasticity:
[legend]

Public Member Functions

 OpCalculatePlasticity (const std::string field_name, MoFEM::Interface &m_field, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< MatrixDouble > m_D_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Protected Attributes

MoFEM::InterfacemField
 
boost::shared_ptr< CommonDatacommonDataPtr
 
boost::shared_ptr< MatrixDouble > mDPtr
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Detailed Description

Examples
plastic.cpp.

Definition at line 119 of file PlasticOps.hpp.

Constructor & Destructor Documentation

◆ OpCalculatePlasticity()

PlasticOps::OpCalculatePlasticity::OpCalculatePlasticity ( const std::string  field_name,
MoFEM::Interface m_field,
boost::shared_ptr< CommonData common_data_ptr,
boost::shared_ptr< MatrixDouble >  m_D_ptr 
)
Examples
PlasticOpsGeneric.hpp.

Definition at line 45 of file PlasticOpsGeneric.hpp.

50 commonDataPtr(common_data_ptr), mDPtr(m_D_ptr) {
51 // Opetor is only executed for vertices
52 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
53}
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
constexpr auto field_name
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
@ OPROW
operator doWork function is executed on FE rows
boost::shared_ptr< CommonData > commonDataPtr
Definition: PlasticOps.hpp:127
boost::shared_ptr< MatrixDouble > mDPtr
Definition: PlasticOps.hpp:128

Member Function Documentation

◆ doWork()

MoFEMErrorCode PlasticOps::OpCalculatePlasticity::doWork ( int  side,
EntityType  type,
EntData data 
)
virtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Examples
PlasticOps.hpp, and PlasticOpsGeneric.hpp.

Definition at line 55 of file PlasticOpsGeneric.hpp.

56 {
58
59 const size_t nb_gauss_pts = getGaussPts().size2();
61 auto t_tau = getFTensor0FromVec(commonDataPtr->plasticTau);
62 auto t_tau_dot = getFTensor0FromVec(commonDataPtr->plasticTauDot);
63 auto t_f = getFTensor0FromVec(commonDataPtr->plasticSurface);
64 auto t_flow =
65 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->plasticFlow);
66 auto t_plastic_strain_dot =
67 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->plasticStrainDot);
68 auto t_stress =
69 getFTensor2SymmetricFromMat<SPACE_DIM>(*(commonDataPtr->mStressPtr));
70
71 auto t_D =
72 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*commonDataPtr->mDPtr);
73 auto t_D_Op = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mDPtr);
74
75 auto t_diff_plastic_strain = diff_tensor();
76 auto t_diff_deviator = diff_deviator(diff_tensor());
77
80 t_flow_dir_dstress(i, j, k, l) =
81 1.5 * (t_diff_deviator(M, N, i, j) * t_diff_deviator(M, N, k, l));
82 t_flow_dir_dstrain(i, j, k, l) =
83 t_flow_dir_dstress(i, j, m, n) * t_D_Op(m, n, k, l);
84
85 commonDataPtr->resC.resize(nb_gauss_pts, false);
86 commonDataPtr->resCdTau.resize(nb_gauss_pts, false);
87 commonDataPtr->resCdStrain.resize(size_symm, nb_gauss_pts, false);
88 commonDataPtr->resCdStrainDot.resize(size_symm, nb_gauss_pts, false);
89 commonDataPtr->resFlow.resize(size_symm, nb_gauss_pts, false);
90 commonDataPtr->resFlowDtau.resize(size_symm, nb_gauss_pts, false);
91 commonDataPtr->resFlowDstrain.resize(size_symm * size_symm, nb_gauss_pts,
92 false);
93 commonDataPtr->resFlowDstrainDot.resize(size_symm * size_symm, nb_gauss_pts,
94 false);
95
96 commonDataPtr->resC.clear();
97 commonDataPtr->resCdTau.clear();
98 commonDataPtr->resCdStrain.clear();
99 commonDataPtr->resCdStrainDot.clear();
100 commonDataPtr->resFlow.clear();
101 commonDataPtr->resFlowDtau.clear();
102 commonDataPtr->resFlowDstrain.clear();
103 commonDataPtr->resFlowDstrainDot.clear();
104
105 auto t_res_c = getFTensor0FromVec(commonDataPtr->resC);
106 auto t_res_c_dtau = getFTensor0FromVec(commonDataPtr->resCdTau);
107 auto t_res_c_dstrain =
108 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->resCdStrain);
109 auto t_res_c_dstrain_dot =
110 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->resCdStrainDot);
111 auto t_res_flow =
112 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->resFlow);
113 auto t_res_flow_dtau =
114 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->resFlowDtau);
115 auto t_res_flow_dstrain = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM>(
116 commonDataPtr->resFlowDstrain);
117 auto t_res_flow_dstrain_dot = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM>(
118 commonDataPtr->resFlowDstrainDot);
119
120 auto next = [&]() {
121 ++t_tau;
122 ++t_tau_dot;
123 ++t_f;
124 ++t_flow;
125 ++t_plastic_strain_dot;
126 ++t_stress;
127 ++t_res_c;
128 ++t_res_c_dtau;
129 ++t_res_c_dstrain;
130 ++t_res_c_dstrain_dot;
131 ++t_res_flow;
132 ++t_res_flow_dtau;
133 ++t_res_flow_dstrain;
134 ++t_res_flow_dstrain_dot;
135 ++t_w;
136 };
137
138 auto get_avtive_pts = [&]() {
139 int nb_points_avtive_on_elem = 0;
140 int nb_points_on_elem = 0;
141
142 auto t_tau = getFTensor0FromVec(commonDataPtr->plasticTau);
143 auto t_tau_dot = getFTensor0FromVec(commonDataPtr->plasticTauDot);
144 auto t_f = getFTensor0FromVec(commonDataPtr->plasticSurface);
145 auto t_plastic_strain_dot =
146 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->plasticStrainDot);
147
148 for (auto &f : commonDataPtr->plasticSurface) {
149 auto eqiv = equivalent_strain_dot(t_plastic_strain_dot);
150 const auto ww = w(eqiv, t_tau_dot, t_f, hardening(t_tau));
151 const auto sign_ww = constrian_sign(ww);
152
153 ++nb_points_on_elem;
154 if (sign_ww > 0) {
155 ++nb_points_avtive_on_elem;
156 }
157
158 ++t_tau;
159 ++t_tau_dot;
160 ++t_f;
161 ++t_plastic_strain_dot;
162 }
163
164 int &active_points = commonDataPtr->activityData[0];
165 int &avtive_full_elems = commonDataPtr->activityData[1];
166 int &avtive_elems = commonDataPtr->activityData[2];
167 int &nb_points = commonDataPtr->activityData[3];
168 int &nb_elements = commonDataPtr->activityData[4];
169
170 ++nb_elements;
171 nb_points += nb_points_on_elem;
172 if (nb_points_avtive_on_elem > 0) {
173 ++avtive_elems;
174 active_points += nb_points_avtive_on_elem;
175 if (nb_points_avtive_on_elem == nb_points_on_elem) {
176 ++avtive_full_elems;
177 }
178 }
179
180 if (nb_points_avtive_on_elem != nb_points_on_elem)
181 return 1;
182 else
183 return 0;
184 };
185
186 if (getTSCtx() == TSMethod::TSContext::CTX_TSSETIJACOBIAN) {
187 get_avtive_pts();
188 }
189
190 for (auto &f : commonDataPtr->plasticSurface) {
191
192 auto eqiv = equivalent_strain_dot(t_plastic_strain_dot);
193 auto t_diff_eqiv = diff_equivalent_strain_dot(eqiv, t_plastic_strain_dot,
194 t_diff_plastic_strain);
195
196 const auto sigma_y = hardening(t_tau);
197 const auto d_sigma_y = hardening_dtau(t_tau);
198
199 auto ww = w(eqiv, t_tau_dot, t_f, sigma_y);
200 auto abs_ww = constrain_abs(ww);
201 auto sign_ww = constrian_sign(ww);
202
203 auto c = constraint(eqiv, t_tau_dot, t_f, sigma_y, abs_ww);
204 auto c_dot_tau = diff_constrain_ddot_tau(sign_ww, eqiv);
205 auto c_equiv = diff_constrain_eqiv(sign_ww, eqiv, t_tau_dot);
206 auto c_sigma_y = diff_constrain_dsigma_y(sign_ww);
207 auto c_f = diff_constrain_df(sign_ww);
208
209 auto t_dev_stress = deviator(t_stress, trace(t_stress));
211 t_flow_dir(k, l) = 1.5 * (t_dev_stress(I, J) * t_diff_deviator(I, J, k, l));
213 t_flow_dstrain(i, j) = t_flow(k, l) * t_D_Op(k, l, i, j);
214
215 auto get_res_c = [&]() { return c; };
216
217 auto get_res_c_dstrain = [&](auto &t_diff_res) {
218 t_diff_res(i, j) = c_f * t_flow_dstrain(i, j);
219 };
220
221 auto get_res_c_dstrain_dot = [&](auto &t_diff_res) {
222 t_diff_res(i, j) = (getTSa() * c_equiv) * t_diff_eqiv(i, j);
223 };
224
225 auto get_res_c_dtau = [&]() {
226 return getTSa() * c_dot_tau + c_sigma_y * d_sigma_y;
227 };
228
229 auto get_res_flow = [&](auto &t_res_flow) {
230 const auto a = sigma_y;
231 const auto b = t_tau_dot;
232 t_res_flow(k, l) = a * t_plastic_strain_dot(k, l) - b * t_flow_dir(k, l);
233 };
234
235 auto get_res_flow_dtau = [&](auto &t_res_flow_dtau) {
236 const auto da = d_sigma_y;
237 const auto db = getTSa();
238 t_res_flow_dtau(k, l) =
239 da * t_plastic_strain_dot(k, l) - db * t_flow_dir(k, l);
240 };
241
242 auto get_res_flow_dstrain = [&](auto &t_res_flow_dstrain) {
243 const auto b = t_tau_dot;
244 t_res_flow_dstrain(m, n, k, l) = -t_flow_dir_dstrain(m, n, k, l) * b;
245 };
246
247 auto get_res_flow_dstrain_dot = [&](auto &t_res_flow_dstrain_dot) {
248 const auto a = sigma_y;
249 t_res_flow_dstrain_dot(m, n, k, l) =
250 (a * getTSa()) * t_diff_plastic_strain(m, n, k, l);
251 };
252
253 t_res_c = get_res_c();
254 get_res_flow(t_res_flow);
255
256 if (getTSCtx() == TSMethod::TSContext::CTX_TSSETIJACOBIAN) {
257 t_res_c_dtau = get_res_c_dtau();
258 get_res_c_dstrain(t_res_c_dstrain);
259 get_res_c_dstrain_dot(t_res_c_dstrain_dot);
260 get_res_flow_dtau(t_res_flow_dtau);
261 get_res_flow_dstrain(t_res_flow_dstrain);
262 get_res_flow_dstrain_dot(t_res_flow_dstrain_dot);
263 }
264
265 next();
266 }
267
269}
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
const double c
speed of light (cm/ns)
auto deviator(FTensor::Tensor2_symmetric< T, SPACE_DIM > &t_stress, double trace)
Definition: PlasticOps.hpp:373
double constrian_sign(double x)
Definition: PlasticOps.hpp:488
FTensor::Index< 'j', SPACE_DIM > j
Definition: PlasticOps.hpp:93
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:103
double diff_constrain_eqiv(double sign, double eqiv, double dot_tau)
Definition: PlasticOps.hpp:540
FTensor::Index< 'l', SPACE_DIM > l
Definition: PlasticOps.hpp:95
FTensor::Index< 'k', SPACE_DIM > k
Definition: PlasticOps.hpp:94
auto diff_tensor()
[Operators definitions]
Definition: PlasticOps.hpp:308
auto diff_constrain_df(double sign)
Definition: PlasticOps.hpp:545
double w(double eqiv, double dot_tau, double f, double sigma_y)
Definition: PlasticOps.hpp:513
FTensor::Index< 'N', 3 > N
Definition: PlasticOps.hpp:104
auto diff_constrain_dsigma_y(double sign)
Definition: PlasticOps.hpp:547
FTensor::Index< 'I', 3 > I
Definition: PlasticOps.hpp:101
auto diff_deviator(FTensor::Ddg< double, SPACE_DIM, SPACE_DIM > &&t_diff_stress)
Definition: PlasticOps.hpp:387
FTensor::Index< 'i', SPACE_DIM > i
[Common data]
Definition: PlasticOps.hpp:92
FTensor::Index< 'm', SPACE_DIM > m
Definition: PlasticOps.hpp:96
double constrain_abs(double x)
Definition: PlasticOps.hpp:502
FTensor::Index< 'J', 3 > J
Definition: PlasticOps.hpp:102
auto diff_equivalent_strain_dot(const T1 eqiv, T2 &t_plastic_strain_dot, T3 &t_diff_plastic_strain)
Definition: PlasticOps.hpp:576
double trace(FTensor::Tensor2_symmetric< T, SPACE_DIM > &t_stress)
Definition: PlasticOps.hpp:364
double constraint(double eqiv, double dot_tau, double f, double sigma_y, double abs_w)
Definition: PlasticOps.hpp:528
double diff_constrain_ddot_tau(double sign, double eqiv)
Definition: PlasticOps.hpp:536
auto equivalent_strain_dot(FTensor::Tensor2_symmetric< T, SPACE_DIM > &t_plastic_strain_dot)
Definition: PlasticOps.hpp:567
FTensor::Index< 'n', SPACE_DIM > n
Definition: PlasticOps.hpp:97
constexpr auto size_symm
Definition: plastic.cpp:33
double hardening(double tau)
Definition: plastic.cpp:122
double hardening_dtau(double tau)
Definition: plastic.cpp:126
auto getFTensor0IntegrationWeight()
Get integration weights.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<CommonData> PlasticOps::OpCalculatePlasticity::commonDataPtr
protected
Examples
PlasticOps.hpp, and PlasticOpsGeneric.hpp.

Definition at line 127 of file PlasticOps.hpp.

◆ mDPtr

boost::shared_ptr<MatrixDouble> PlasticOps::OpCalculatePlasticity::mDPtr
protected
Examples
PlasticOps.hpp, and PlasticOpsGeneric.hpp.

Definition at line 128 of file PlasticOps.hpp.

◆ mField

MoFEM::Interface& PlasticOps::OpCalculatePlasticity::mField
protected
Examples
PlasticOps.hpp.

Definition at line 126 of file PlasticOps.hpp.


The documentation for this struct was generated from the following files: