v0.14.0
HMHNeohookean.cpp
Go to the documentation of this file.
1 /**
2  * @file HMHNeohookean.cpp
3  * @brief Implementation of NeoHookean material
4  * @date 2024-08-31
5  *
6  * @copyright Copyright (c) 2024
7  *
8  */
9 
10 namespace EshelbianPlasticity {
11 
13 
14  static constexpr int numberOfActiveVariables = 9;
15  static constexpr int numberOfDependentVariables = 9;
16 
17  HMHNeohookean(const double c10, const double K)
19  c10(c10), K(K) {}
20 
23  MoFEMErrorCode evaluateRhs(EntData &data) { return 0; }
24  MoFEMErrorCode evaluateLhs(EntData &data) { return 0; }
25  };
26 
28  returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
29  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
30  boost::shared_ptr<PhysicalEquations> physics_ptr) {
31  return (new OpJacobian(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
32  }
33 
36  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "neo_hookean_", "", "none");
37 
38  c10 = 1;
39  CHKERR PetscOptionsScalar("-c10", "C10", "", c10, &c10, PETSC_NULL);
40  K = 1;
41  CHKERR PetscOptionsScalar("-K", "Bulk modulus K", "", K, &K, PETSC_NULL);
42 
43  ierr = PetscOptionsEnd();
44  CHKERRG(ierr);
46  }
47 
48  double c10;
49  double K;
50 
51  MoFEMErrorCode recordTape(const int tape, DTensor2Ptr *t_h_ptr) {
55  }
56 
58 
59  OpSpatialPhysical(const std::string &field_name,
60  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
61  const double alpha_u);
62 
64 
65  private:
66  const double alphaU;
67  };
68 
69  virtual VolUserDataOperator *
70  returnOpSpatialPhysical(const std::string &field_name,
71  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
72  const double alpha_u) {
73  return new OpSpatialPhysical(field_name, data_ptr, alpha_u);
74  }
75 
77  OpSpatialPhysical_du_du(std::string row_field, std::string col_field,
78  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
79  const double alpha);
80  MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
81 
82  private:
83  const double alphaU;
84  };
85 
87  std::string row_field, std::string col_field,
88  boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha) {
89  return new OpSpatialPhysical_du_du(row_field, col_field, data_ptr,
90  alpha);
91  }
92 };
93 
95  const std::string &field_name,
96  boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha_u)
97  : OpAssembleVolume(field_name, data_ptr, OPROW), alphaU(alpha_u) {}
98 
101 
102  auto neohookean_ptr =
103  boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
104  if (!neohookean_ptr) {
105  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
106  "Pointer to HMHNeohookean is null");
107  }
108 
109 #ifdef NDEBUG
111  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
112  "Stretch selector is not equal to LOG");
113  } else {
114  if (EshelbianCore::exponentBase != exp(1)) {
115  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
116  "Exponent base is not equal to exp(1)");
117  }
118  }
119 #endif // NDEBUG
120 
121  const auto c10 = neohookean_ptr->c10;
122  const auto K = neohookean_ptr->K;
123 
125  auto t_L = symm_L_tensor();
126 
127  constexpr auto t_kd_sym = FTensor::Kronecker_Delta_symmetric<int>();
128  constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
129 
130  int nb_dofs = data.getIndices().size();
131  int nb_integration_pts = data.getN().size1();
132  auto v = getVolume();
133  auto t_w = getFTensor0IntegrationWeight();
134  auto t_approx_P_adjont_log_du =
135  getFTensor1FromMat<size_symm>(dataAtPts->adjointPdUAtPts);
136  auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
137  auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
138  auto t_total_log_u =
139  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTotalTensorAtPts);
140  auto t_dot_log_u =
141  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
142  auto t_diff_u =
143  getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
144  auto t_log_u =
145  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTensorAtPts);
146  auto t_log_u2_h1 =
147  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretch2H1AtPts);
148 
149  auto t_diff = diff_tensor();
150 
157 
158  auto get_ftensor2 = [](auto &v) {
160  &v[0], &v[1], &v[2], &v[3], &v[4], &v[5]);
161  };
162 
163  int nb_base_functions = data.getN().size2();
164  auto t_row_base_fun = data.getFTensor0N();
165  for (int gg = 0; gg != nb_integration_pts; ++gg) {
166  double a = v * t_w;
167  ++t_w;
168 
171 
173  case NO_H1_CONFIGURATION:
174  t_h1(i, j) = t_kd(i, j);
175  // t_h1_du(i, j, m, n) = t_kd(i, m) * t_kd(j, n);
176  t_Ldiff_u(i, j, L) = t_diff_u(i, j, k, l) * t_L(k, l, L);
177  break;
178  case LARGE_ROT:
179  case MODERATE_ROT:
180  t_h1(i, j) = t_grad_h1(i, j) + t_kd(i, j);
181  // t_h1_du(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_h1(k, j));
182  t_Ldiff_u(i, j, L) = (t_diff_u(i, m, k, l) * t_h1(m, j)) * t_L(k, l, L);
183  break;
184  case SMALL_ROT:
185  t_h1(i, j) = t_kd(i, j);
186  // t_h1_du(i, j, m, n) = t_kd(i, m) * t_kd(j, n);
187  t_Ldiff_u(i, j, L) = t_diff_u(i, j, k, l) * t_L(k, l, L);
188  break;
189  };
190  ++t_diff_u;
191  ++t_grad_h1;
192 
194  t_Sigma_u(i, j) = 2.0 * c10 * (t_u(i, m) * t_h1(m, j));
195  const double tr = t_total_log_u(i, j) * t_kd_sym(i, j);
196  const double J = EshelbianCore::f(tr);
197  const double Simga_J = -2 * c10 + K * (J - 1) * J;
198 
200  t_viscous_P(i, j) =
201  alphaU *
202  t_dot_log_u(i,
203  j); // That is cheat, should be split on axial and deviator
204 
206  t_residual(L) = t_approx_P_adjont_log_du(L);
207  t_residual(L) -= t_Ldiff_u(i, j, L) * t_Sigma_u(i, j);
208  t_residual(L) -= (t_L(i, j, L) * t_kd(i, j)) * Simga_J;
209  t_residual(L) -= t_L(i, j, L) * t_viscous_P(i, j);
210  t_residual(L) *= a;
211 
212  ++t_approx_P_adjont_log_du;
213  ++t_u;
214  ++t_total_log_u;
215  ++t_dot_log_u;
216 
217  auto t_nf = getFTensor1FromPtr<size_symm>(&*nF.data().begin());
218  int bb = 0;
219  for (; bb != nb_dofs / size_symm; ++bb) {
220  t_nf(L) += t_row_base_fun * t_residual(L);
221  ++t_nf;
222  ++t_row_base_fun;
223  }
224  for (; bb != nb_base_functions; ++bb)
225  ++t_row_base_fun;
226  }
228 }
229 
231  std::string row_field, std::string col_field,
232  boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha)
233  : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
234  alphaU(alpha) {
235  sYmm = false;
236 }
237 
240  EntData &col_data) {
242 
243  auto neohookean_ptr =
244  boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
245  if (!neohookean_ptr) {
246  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
247  "Pointer to HMHNeohookean is null");
248  }
249 
250 #ifdef NDEBUG
252  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
253  "Stretch selector is not equal to LOG");
254  } else {
255  if (EshelbianCore::exponentBase != exp(1)) {
256  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
257  "Exponent base is not equal to exp(1)");
258  }
259  }
260 #endif // NDEBUG
261 
262  const auto c10 = neohookean_ptr->c10;
263  const auto lambda = neohookean_ptr->K;
264 
267 
268  constexpr auto t_kd_sym = FTensor::Kronecker_Delta_symmetric<int>();
269  constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
270 
271  auto t_L = symm_L_tensor();
272  auto t_diff = diff_tensor();
273 
274  int nb_integration_pts = row_data.getN().size1();
275  int row_nb_dofs = row_data.getIndices().size();
276  int col_nb_dofs = col_data.getIndices().size();
277 
278  auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
280  size_symm>(
281 
282  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 0, c + 3),
283  &m(r + 0, c + 4), &m(r + 0, c + 5),
284 
285  &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 1, c + 3),
286  &m(r + 1, c + 4), &m(r + 1, c + 5),
287 
288  &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2), &m(r + 2, c + 3),
289  &m(r + 2, c + 4), &m(r + 2, c + 5),
290 
291  &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2), &m(r + 3, c + 3),
292  &m(r + 3, c + 4), &m(r + 3, c + 5),
293 
294  &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2), &m(r + 4, c + 3),
295  &m(r + 4, c + 4), &m(r + 4, c + 5),
296 
297  &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2), &m(r + 5, c + 3),
298  &m(r + 5, c + 4), &m(r + 5, c + 5)
299 
300  );
301  };
302 
309 
310  auto v = getVolume();
311  auto ts_a = getTSa();
312  auto t_w = getFTensor0IntegrationWeight();
313 
314  int row_nb_base_functions = row_data.getN().size2();
315  auto t_row_base_fun = row_data.getFTensor0N();
316 
317  auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
318  auto t_diff_u =
319  getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
320  auto t_total_log_u =
321  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTotalTensorAtPts);
322  auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
323  auto t_approx_P_adjont_dstretch =
324  getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
325  auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
326  auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
327  auto &nbUniq = dataAtPts->nbUniq;
328 
329  auto t_P = getFTensor2FromMat<3, 3>(dataAtPts->PAtPts);
330  auto r_P_du = getFTensor4FromMat<3, 3, 3, 3>(dataAtPts->P_du);
331 
332  for (int gg = 0; gg != nb_integration_pts; ++gg) {
333  double a = v * t_w;
334  ++t_w;
335 
338  t_h1(i, j) = t_grad_h1(i, j) + t_kd(i, j);
339 
341  case NO_H1_CONFIGURATION:
342  t_Ldiff_u(i, j, L) = t_diff_u(i, j, k, l) * t_L(k, l, L);
343  break;
344  case LARGE_ROT:
345  case MODERATE_ROT:
346  t_Ldiff_u(i, j, L) = (t_diff_u(i, m, k, l) * t_h1(m, j)) * t_L(k, l, L);
347  break;
348  case SMALL_ROT:
349  t_Ldiff_u(i, j, L) = t_diff_u(i, j, k, l) * t_L(k, l, L);
350  break;
351  };
352  ++t_diff_u;
353  ++t_grad_h1;
354 
355  const double tr = t_total_log_u(i, j) * t_kd_sym(i, j);
356  const double det_u = EshelbianCore::f(tr);
358  Simga_J_dtr(L) =
359  (+lambda * (2 * det_u * det_u - det_u)) * (t_kd(i, j) * t_L(i, j, L));
360  ++t_total_log_u;
361 
363  t_dP(L, J) = (-2.0 * c10) * (t_Ldiff_u(i, j, L) * t_Ldiff_u(i, j, J));
364  t_dP(L, J) -= (t_L(i, j, L) * t_kd(i, j)) * Simga_J_dtr(J);
365  t_dP(L, J) -=
366  (alphaU * ts_a) * (t_L(i, j, L) * (t_diff(i, j, k, l) * t_L(k, l, J)));
367 
369  t_Sigma_u(i, j) = 2.0 * c10 * (t_u(i, m) * t_h1(m, j));
370  ++t_u;
371 
373  t_deltaP(i, j) =
374  t_approx_P_adjont_dstretch(i, j) - t_h1(j, n) * t_Sigma_u(i, n);
375  ++t_approx_P_adjont_dstretch;
376 
379  t_deltaP_sym(i, j) = (t_deltaP(i, j) || t_deltaP(j, i));
380  t_deltaP_sym(i, j) /= 2.0;
381  auto t_diff2_uP2 = EigenMatrix::getDiffDiffMat(
382  t_eigen_vals, t_eigen_vecs, EshelbianCore::f, EshelbianCore::d_f,
383  EshelbianCore::dd_f, t_deltaP_sym, nbUniq[gg]);
384  t_dP(L, J) += t_L(i, j, L) * (t_diff2_uP2(i, j, k, l) * t_L(k, l, J));
385  }
386  ++t_eigen_vals;
387  ++t_eigen_vecs;
388 
389  int rr = 0;
390  for (; rr != row_nb_dofs / size_symm; ++rr) {
391  auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
392  auto t_m = get_ftensor2(K, 6 * rr, 0);
393  for (int cc = 0; cc != col_nb_dofs / size_symm; ++cc) {
394  double b = a * t_row_base_fun * t_col_base_fun;
395  t_m(L, J) += b * t_dP(L, J);
396  ++t_m;
397  ++t_col_base_fun;
398  }
399  ++t_row_base_fun;
400  }
401 
402  for (; rr != row_nb_base_functions; ++rr) {
403  ++t_row_base_fun;
404  }
405 
406  ++t_P;
407  ++r_P_du;
408  }
409 
411 }
412 
413 } // namespace EshelbianPlasticity
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical_du_du::OpSpatialPhysical_du_du
OpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
Definition: HMHNeohookean.cpp:230
EshelbianPlasticity::LINEAR
@ LINEAR
Definition: EshelbianPlasticity.hpp:44
EshelbianPlasticity::size_symm
constexpr static auto size_symm
Definition: EshelbianAux.hpp:39
EshelbianPlasticity::EshelbianCore::dd_f
static boost::function< double(const double)> dd_f
Definition: EshelbianPlasticity.hpp:905
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical_du_du::integrate
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
Definition: HMHNeohookean.cpp:239
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EshelbianPlasticity::EshelbianCore::exponentBase
static double exponentBase
Definition: EshelbianPlasticity.hpp:902
EshelbianPlasticity::HMHNeohookean::numberOfActiveVariables
static constexpr int numberOfActiveVariables
Definition: HMHNeohookean.cpp:14
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical::OpSpatialPhysical
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
Definition: HMHNeohookean.cpp:94
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
FTensor::Kronecker_Delta
Kronecker Delta class.
Definition: Kronecker_Delta.hpp:15
EshelbianPlasticity
Definition: CGGTonsorialBubbleBase.hpp:11
EshelbianPlasticity::diff_tensor
auto diff_tensor()
Definition: EshelbianAux.hpp:41
EshelbianPlasticity::PhysicalEquations
Definition: EshelbianPlasticity.hpp:219
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical
Definition: HMHNeohookean.cpp:57
J
FTensor::Index< 'J', DIM1 > J
Definition: level_set.cpp:30
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
EshelbianPlasticity::EshelbianCore::d_f
static boost::function< double(const double)> d_f
Definition: EshelbianPlasticity.hpp:904
FTENSOR_INDEX
#define FTENSOR_INDEX(DIM, I)
Definition: Templates.hpp:2011
EshelbianPlasticity::NO_H1_CONFIGURATION
@ NO_H1_CONFIGURATION
Definition: EshelbianPlasticity.hpp:43
sdf.r
int r
Definition: sdf.py:8
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
EshelbianPlasticity::HMHNeohookean::OpJacobian::evaluateRhs
MoFEMErrorCode evaluateRhs(EntData &data)
Definition: HMHNeohookean.cpp:23
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
FTensor::Tensor3
Definition: Tensor3_value.hpp:12
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical_du_du
Definition: HMHNeohookean.cpp:76
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
EshelbianPlasticity::EshelbianCore::gradApproximator
static enum RotSelector gradApproximator
Definition: EshelbianPlasticity.hpp:897
a
constexpr double a
Definition: approx_sphere.cpp:30
EshelbianPlasticity::HMHNeohookean::OpJacobian::evaluateLhs
MoFEMErrorCode evaluateLhs(EntData &data)
Definition: HMHNeohookean.cpp:24
EshelbianPlasticity::EshelbianCore::stretchSelector
static enum StretchSelector stretchSelector
Definition: EshelbianPlasticity.hpp:898
EshelbianPlasticity::OpAssembleVolume
Definition: EshelbianPlasticity.hpp:517
EshelbianPlasticity::HMHNeohookean::K
double K
Definition: HMHNeohookean.cpp:49
EshelbianPlasticity::EshelbianCore::f
static boost::function< double(const double)> f
Definition: EshelbianPlasticity.hpp:903
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
EshelbianPlasticity::HMHNeohookean
Definition: HMHNeohookean.cpp:12
EshelbianPlasticity::HMHNeohookean::numberOfDependentVariables
static constexpr int numberOfDependentVariables
Definition: HMHNeohookean.cpp:15
EshelbianPlasticity::HMHNeohookean::recordTape
MoFEMErrorCode recordTape(const int tape, DTensor2Ptr *t_h_ptr)
Definition: HMHNeohookean.cpp:51
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical::integrate
MoFEMErrorCode integrate(EntData &data)
Definition: HMHNeohookean.cpp:99
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical::alphaU
const double alphaU
Definition: HMHNeohookean.cpp:66
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
EshelbianPlasticity::symm_L_tensor
auto symm_L_tensor()
Definition: EshelbianAux.hpp:52
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index
Definition: Index.hpp:23
convert.n
n
Definition: convert.py:82
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
EshelbianPlasticity::HMHNeohookean::returnOpSpatialPhysical
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
Definition: HMHNeohookean.cpp:70
EshelbianPlasticity::HMHNeohookean::getOptions
MoFEMErrorCode getOptions()
Definition: HMHNeohookean.cpp:34
EshelbianPlasticity::SMALL_ROT
@ SMALL_ROT
Definition: EshelbianPlasticity.hpp:43
EshelbianPlasticity::OpJacobian::OpJacobian
OpJacobian()
Definition: EshelbianPlasticity.hpp:394
EshelbianPlasticity::LARGE_ROT
@ LARGE_ROT
Definition: EshelbianPlasticity.hpp:43
EshelbianPlasticity::LOG
@ LOG
Definition: EshelbianPlasticity.hpp:44
EshelbianPlasticity::HMHNeohookean::returnOpJacobian
virtual EshelbianPlasticity::OpJacobian * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
Definition: HMHNeohookean.cpp:28
EshelbianPlasticity::HMHNeohookean::HMHNeohookean
HMHNeohookean(const double c10, const double K)
Definition: HMHNeohookean.cpp:17
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
EshelbianPlasticity::OpJacobian
Definition: EshelbianPlasticity.hpp:380
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
EshelbianPlasticity::HMHNeohookean::OpSpatialPhysical_du_du::alphaU
const double alphaU
Definition: HMHNeohookean.cpp:83
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
EigenMatrix::getDiffDiffMat
FTensor::Ddg< double, 3, 3 > getDiffDiffMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, FTensor::Tensor2< double, 3, 3 > &t_S, const int nb)
Definition: MatrixFunction.cpp:78
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
EshelbianPlasticity::HMHNeohookean::OpJacobian
Definition: HMHNeohookean.cpp:21
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
EshelbianPlasticity::HMHNeohookean::c10
double c10
Definition: HMHNeohookean.cpp:48
EshelbianPlasticity::MODERATE_ROT
@ MODERATE_ROT
Definition: EshelbianPlasticity.hpp:43
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
EshelbianPlasticity::HMHNeohookean::returnOpSpatialPhysical_du_du
VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
Definition: HMHNeohookean.cpp:86