v0.14.0
ADOLCPlasticity.hpp
Go to the documentation of this file.
1 /** \file ADOLCPlasticity.hpp
2  * \ingroup adoc_plasticity
3  * \example ADOLCPlasticity.hpp
4  *
5  * \brief Operators and data structures for small strain plasticity
6  *
7  * \defgroup adoc_plasticity ADOL-C plasticity
8  * \ingroup user_modules
9  * \defgroup user_modules User modules
10  *
11 **/
12 
13 #ifndef __ADOLCPLASTICITY_HPP_
14 #define __ADOLCPLASTICITY_HPP_
15 
16 #ifndef WITH_ADOL_C
17 #error "MoFEM need to be compiled with ADOL-C"
18 #endif
19 
20 /**
21  * \ingroup adoc_plasticity
22  *
23  */
24 namespace ADOLCPlasticity {
25 
26 /**
27  * \brief Op convert Vight strain vector to strain tensor
28 */
30  return FTensor::Dg<double, 3, 6>{1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
31  0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0,
32  0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
33  0.0, 0.0, 0.5, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0};
34 };
35 
36 /**
37  * \brief Op convert strain tensor to Vight strain vector
38 */
40  return FTensor::Dg<double, 3, 6>{1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
41  1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0,
42  0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
43  0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0};
44 };
45 
46 /**
47  * \brief Op convert Vight stress vector to stress tensor
48 */
50  return FTensor::Dg<double, 3, 6>{1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
51  1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0,
52  0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
53  0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0};
54 };
55 
56 /** \brief common data used by volume elements
57  * \ingroup nonlinear_elastic_elem
58  */
59 struct CommonData : boost::enable_shared_from_this<CommonData> {
60 
63 
64  inline auto getFTensor1StressAtGaussPts(int gg = 0) {
66  &gradientW(gg, 7),
67  &gradientW(gg, 8),
68  &gradientW(gg, 9),
69  &gradientW(gg, 10),
70  &gradientW(gg, 11),
71  static_cast<int>(gradientW.size2())};
72  };
73 
74  inline auto getFTensor1PlasticStrainAtGaussPts(int gg = 0) {
76  &activeVariablesW(gg, 0),
77  &activeVariablesW(gg, 1),
78  &activeVariablesW(gg, 2),
79  &activeVariablesW(gg, 3),
80  &activeVariablesW(gg, 4),
81  &activeVariablesW(gg, 5),
82  static_cast<int>(activeVariablesW.size2())};
83  };
84 
87 
88  inline auto getPlasticStrain(int gg = 0) {
89  return getVectorAdaptor(&(activeVariablesW(gg, 0)), 6);
90  }
91 
92  inline auto getInternalVariables(int gg = 0) {
93  return getVectorAdaptor(&(activeVariablesW(gg, 12)),
94  activeVariablesW.size2() - 12);
95  }
96 
97  vector<double> deltaGamma; //< Lagrange plastic multiplier
102 
103  inline auto getGardAtGaussPtsPtr() {
104  return boost::shared_ptr<MatrixDouble>(shared_from_this(), &gradAtGaussPts);
105  }
106 
107  bool bBar = true;
108 
111  PetscBool b_bar = bBar ? PETSC_TRUE : PETSC_FALSE;
112  CHKERR PetscOptionsGetBool(PETSC_NULL, PETSC_NULL, "-b_bar", &b_bar,
113  PETSC_NULL);
114  bBar = b_bar;
116  }
117 
118  boost::shared_ptr<MatrixDouble> getStressMatrixPtr() {
119  return boost::shared_ptr<MatrixDouble>(shared_from_this(), &stressMatrix);
120  }
121  boost::shared_ptr<MatrixDouble> getPlasticStrainMatrixPtr() {
122  return boost::shared_ptr<MatrixDouble>(shared_from_this(),
124  }
125 
128 
130 
131  CHK_THROW_MESSAGE(getDefaultMaterialParameters(), "get parameters failed");
132  }
133 };
134 
135 /** \brief Closest point projection algorithm
136 * \ingroup small_strain_plasticity
137 
138 */
140 
142  boost::function<int(int, int, int)> integrationRule = [](int, int, int p) {
143  return 2 * (p - 1);
144  };
145 
148 
150  return getVectorAdaptor(&(activeVariablesW[0]), 6);
151  }
153  return getVectorAdaptor(&(activeVariablesW[6]), 6);
154  }
157  }
160  }
162  return getVectorAdaptor(&(gradientW[6]), 6);
163  }
166  }
167 
168  enum TypesTags { W = 0, Y, H, LAST_TAPE }; //< W - energy, Y - yield, H - flow
169  std::array<int, LAST_TAPE> tapesTags; //< Tapes nmbers
170 
172 
175 
176  double w;
177  double y;
178  double h;
179 
180  double deltaGamma; // Increment of plastic multiplier
182  ublas::symmetric_matrix<double, ublas::lower> C, D;
183 
184  /**
185  * \brief Record strain energy
186  */
188 
189  /**
190  * \brief Record yield function
191  */
193 
194  /**
195  * \brief Record flow potential
196  */
198 
199  MatrixDouble hessianW; //< Hessian of energy
202 
203  VectorDouble gradientY; //< Gradient of yield function
206 
207  VectorDouble gradientH; //< Gradient of flow potential
208  MatrixDouble hessianH; //< Hessian of flow potential
211 
212  MoFEMErrorCode createMatAVecR(); //< For integration point SNES solver
213 
214  MoFEMErrorCode evaluatePotentials(); //< Evaluate potentials
215 
217  playPotentials(); //< Play potentials from recorded ADOl-C tapes
218 
220  playPotentials_NoHessian(); //< Play potentials from recorded ADOl-C tapes
221 
222  MoFEMErrorCode calculateR(Vec R); //< Calculate residual
223 
224  MoFEMErrorCode calculateA(); //< Calculate tangent matrix
225 
226  /**
227  * \brief Function executed by nested SENES at evaluationg residual
228  */
229  friend MoFEMErrorCode ADOLCPlasticityRes(SNES, Vec, Vec, void *ctx);
230 
231  /**
232  * \brief Function executed by nested SENES at evaluationg tangent matrix
233  */
234  friend MoFEMErrorCode ADOLCPlasticityJac(SNES, Vec, Mat, Mat, void *ctx);
235 
236  /**
237  * \brief Create nested snes
238  */
240 
241  /**
242  * \brief Solve nonlinear system of equations to find stress, internal
243  * fluxes, and Lagrange plastic multiplier
244  */
246 
247  /**
248  * \brief Calculate consistent tangent matrix
249  */
251 
252  /**
253  * \brief Record tapes
254  */
256 
257  /**
258  * \brief Get model parameters from blocks
259  */
260  virtual MoFEMErrorCode
262  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
263  std::string block_name, Sev sev) {
264  return 0;
265  }
266 
267  /**
268  * \brief Set parameters for ADOL-C of constitutive relations
269  *
270  * \param tag [in] - tag of the tape
271  * \param recalculate_elastic_tangent [out] - if setParam set it to true,
272  * tangent matrix for elastic domain should be recalculated
273  */
274  virtual MoFEMErrorCode setParams(short tag,
275  bool &recalculate_elastic_tangent) {
276  return 0;
277  }
278 
279  /**
280  * \brief Set Hemholtz energy
281  */
282  virtual MoFEMErrorCode freeHemholtzEnergy() = 0;
283 
284  /**
285  * \brief Set yield function
286  */
287  virtual MoFEMErrorCode yieldFunction() = 0;
288 
289  /**
290  * \brief Set flow potential
291  */
292  virtual MoFEMErrorCode flowPotential() = 0;
293 
294  // protected:
295 
296  MatrixDouble partialWStrainPlasticStrain; //< Partial derivative of free energy
297  // with respect to plastic strain
298  VectorAdaptor partialYSigma; //< Partial derivative of yield function with
299  // respect to stress
300  VectorAdaptor partialYFlux; //< Partial derivative of yield function with
301  // respect to internal fluxes
302  VectorAdaptor partialHSigma; //< Partial derivative of flow potential with
303  // respect to stress
304  VectorAdaptor partialHFlux; //< Partial derivative of flow potential with
305  // respect to internal fluxes
306 
307  /// Second partial derivative of flow potential with respect to stresses
308  /// and internal
309  ublas::symmetric_matrix<double, ublas::lower> partial2HSigma;
310  /// Second partial derivative of flow potential with respect to internal
311  /// fluxes
312  ublas::symmetric_matrix<double, ublas::lower> partial2HFlux;
313  /// Mixed decond partial derivative of flow potential with respect to stresses
314  /// and internal fluxes
316 
317  SmartPetscObj<Mat> A; //< Nested SNES tangent matrix
318  SmartPetscObj<Vec> R; //< Nested SNES residual
319  SmartPetscObj<Vec> Chi; //< Nested SNES unknown vector
320  SmartPetscObj<Vec> dChi; //< Nested SNES increment vector
321  ublas::matrix<double, ublas::column_major> dataA;
322  SmartPetscObj<SNES> sNes; //< Nested SNES solver
323 
324  ublas::vector<adouble> a_sTrain;
325  ublas::vector<adouble> a_plasticStrain;
326  ublas::vector<adouble> a_internalVariables;
327  ublas::vector<adouble> a_sTress, a_internalFluxes;
331 };
332 
333 /**
334  * \brief Internal SNES function used at integration points to calulate stress
335 */
336 MoFEMErrorCode ADOLCPlasticityRes(SNES snes, Vec chi, Vec r, void *ctx);
337 
338 /**
339  * \brief Internal SNES function used at integration points to calulate
340  * tangent matrix
341 */
342 MoFEMErrorCode ADOLCPlasticityJac(SNES, Vec, Mat, Mat, void *ctx);
343 
344 /**
345  * \brief Get opreator to calulate stress
346 */
347 template <int DIM>
349  MoFEM::Interface &m_field, boost::shared_ptr<CommonData> common_data_ptr,
350  boost::shared_ptr<ClosestPointProjection> cp_ptr, bool calc_lhs);
351 
352 template <>
354  MoFEM::Interface &m_field, boost::shared_ptr<CommonData> common_data_ptr,
355  boost::shared_ptr<ClosestPointProjection> cp_ptr, bool calc_lhs);
356 
357 template <>
359  MoFEM::Interface &m_field, boost::shared_ptr<CommonData> common_data_ptr,
360  boost::shared_ptr<ClosestPointProjection> cp_ptr, bool calc_lhs);
361 
362 /**
363  * \brief Assemble right hand side
364  *
365  * \param field_name [in] - name of field
366  * \param common_data_ptr [in] - shared pointer to common data
367  *
368  * \ingroup small_strain_plasticity
369  */
370 template <int DIM, IntegrationType I, typename AssemblyDomainEleOp>
371 struct OpRhsImpl;
372 
373 /**
374  * \brief Assemble left hand side
375  *
376  * \param field_name [in] - name of field
377  * \param common_data_ptr [in] - shared pointer to common data
378  *
379  * \ingroup small_strain_plasticity
380  */
381 template <int DIM, IntegrationType I, typename AssemblyDomainEleOp>
382 struct OpLhsImpl;
383 
384 template <typename DomainEleOp> struct ADOLCPlasticityIntegrators {
385 
386  template <AssemblyType A> struct Assembly {
387 
388  using AssemblyDomainEleOp =
389  typename FormsIntegrators<DomainEleOp>::template Assembly<A>::OpBase;
390 
391  template <int DIM, IntegrationType I>
393 
394  template <int DIM, IntegrationType I>
396  };
397 };
398 
399 using Pip = boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator>;
400 
401 /**
402  * @brief Assemble the left hand side, i.e. tangent matrix
403  * @ingroup adoc_plasticity
404  *
405  * @tparam DIM dimension of the problem
406  * @tparam A assembly type
407  * @tparam I integration type
408  * @tparam DomainEleOp operator type
409  * @param m_field
410  * @param field_name
411  * @param pip
412  * @param block_name esh block name caring material parameters
413  * @param common_data_ptr
414  * @param cp_ptr
415  * @return MoFEMErrorCode
416  */
417 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
419  MoFEM::Interface &m_field, string field_name, Pip &pip,
420  std::string block_name, boost::shared_ptr<CommonData> common_data_ptr,
421  boost::shared_ptr<ClosestPointProjection> cp_ptr, Sev sev = Sev::inform) {
424  CHKERR cp_ptr->addMatBlockOps(m_field, pip, block_name, sev);
425  pip.push_back(
426  getRawPtrOpCalculateStress<DIM>(m_field, common_data_ptr, cp_ptr, false));
427  pip.push_back(new typename P::template Assembly<A>::template OpRhs<DIM, I>(
428  field_name, common_data_ptr));
430 }
431 /**
432  * @brief Assemble the left hand side, i.e. tangent matrix
433  * @ingroup adoc_plasticity
434  *
435  * @tparam DIM dimension of the problem
436  * @tparam A assembly type
437  * @tparam I integration type
438  * @tparam DomainEleOp operator type
439  * @param m_field
440  * @param field_name
441  * @param pip
442  * @param block_name esh block name caring material parameters
443  * @param common_data_ptr
444  * @param cp_ptr
445  * @return MoFEMErrorCode
446  */
447 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
450  std::string block_name,
451  boost::shared_ptr<CommonData> common_data_ptr,
452  boost::shared_ptr<ClosestPointProjection> cp_ptr) {
455  CHKERR cp_ptr->addMatBlockOps(m_field, pip, block_name, Sev::noisy);
456  pip.push_back(
457  getRawPtrOpCalculateStress<DIM>(m_field, common_data_ptr, cp_ptr, true));
458  pip.push_back(new typename P::template Assembly<A>::template OpLhs<DIM, I>(
459  field_name, common_data_ptr));
461 }
462 
463 /**
464  * @brief Push operators to update history variables
465  * @ingroup adoc_plasticity
466  *
467  * @tparam DIM dimension of the problem
468  * @param m_field core interface
469  * @param pip
470  * @param block_name mesh block name caring material parameters
471  * @param common_data_ptr
472  * @param cp_ptr
473  * @return MoFEMErrorCode
474  */
475 template <int DIM>
478  std::string block_name,
479  boost::shared_ptr<CommonData> common_data_ptr,
480  boost::shared_ptr<ClosestPointProjection> cp_ptr);
481 
482 /**
483  * @copydoc ADOLCPlasticity::opFactoryDomainUpdate
484  */
485 template <>
488  std::string block_name,
489  boost::shared_ptr<CommonData> common_data_ptr,
490  boost::shared_ptr<ClosestPointProjection> cp_ptr);
491 
492 /**
493  * @copydoc ADOLCPlasticity::opFactoryDomainUpdate
494  */
495 template <>
498  std::string block_name,
499  boost::shared_ptr<CommonData> common_data_ptr,
500  boost::shared_ptr<ClosestPointProjection> cp_ptr);
501 
502 
503 
504 /**
505  * \brief Update internal fluxes (update history variables)
506  * \ingroup adoc_plasticity
507  */
508 struct TSUpdate {
509  virtual MoFEMErrorCode postProcess(TS ts) = 0;
510 };
511 boost::shared_ptr<TSUpdate> createTSUpdate(std::string fe_name,
512  boost::shared_ptr<FEMethod> fe_ptr);
513 
514 /**
515  * \brief Calculate tensorial base functions. Apply bBar method when needed
516  */
517 struct MakeB {
518 
521  MatrixDouble &storage, const bool b_bar,
522  const int nb_integration_pts, double *w_ptr,
524 
527  MatrixDouble &storage, const bool b_bar,
528  const int nb_integration_pts, double *w_ptr,
530 };
531 
532 
533 template <int DIM, typename AssemblyDomainEleOp>
535  OpRhsImpl(std::string field_name,
536  boost::shared_ptr<CommonData> common_data_ptr);
538 
539 private:
540  boost::shared_ptr<CommonData> commonDataPtr;
541  MatrixDouble baseStorage; ///< Store tensorial base functions
542 };
543 
544 template <int DIM, typename AssemblyDomainEleOp>
546  OpLhsImpl(std::string field_name,
547  boost::shared_ptr<CommonData> common_data_ptr);
550 
551 protected:
552  boost::shared_ptr<CommonData> commonDataPtr;
553  MatrixDouble baseRowStorage; ///< Store tensorial base functions
554  MatrixDouble baseColStorage; ///< Store tensorial base functions
555 };
556 
557 template <int DIM, typename AssemblyDomainEleOp>
559  string field_name, boost::shared_ptr<CommonData> common_data_ptr)
561  commonDataPtr(common_data_ptr) {}
562 
563 //! [OpRhsImpl integrate]
564 template <int DIM, typename AssemblyDomainEleOp>
568  using OP = AssemblyDomainEleOp;
571  double *w_ptr = &(OP::getGaussPts()(DIM, 0));
572  // Calulate tensorial base functions. Apply bBar method when needed
574  data, baseStorage, commonDataPtr->bBar, OP::getGaussPts().size2(), w_ptr,
576  auto t_stress = getFTensor2SymmetricFromMat<3>(commonDataPtr->stressMatrix);
577  const auto vol = OP::getMeasure();
578  auto t_w = OP::getFTensor0IntegrationWeight();
579  for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
580  double alpha = vol * t_w;
581  for (int bb = 0; bb != OP::nbRows; ++bb) {
582  OP::locF[bb] += alpha * (t_stress(i, j) * t_diff(i, j));
583  ++t_diff;
584  }
585  ++t_w;
586  ++t_stress;
587  }
589 }
590 //! [OpRhsImpl integrate]
591 
592 template <int DIM, typename AssemblyDomainEleOp>
594  string field_name, boost::shared_ptr<CommonData> common_data_ptr)
596  AssemblyDomainEleOp::OPROWCOL),
597  commonDataPtr(common_data_ptr) {
598  this->sYmm = false; // It has to be false for not-associtive plasticity
599 }
600 
601 template <int DIM, typename AssemblyDomainEleOp>
606  using OP = AssemblyDomainEleOp;
607 
608  double *w_ptr = &(OP::getGaussPts()(DIM, 0));
609  auto t_diff_row = MakeB::getFTensor2SymmetricDiffBase(
610  row_data, baseRowStorage, commonDataPtr->bBar, OP::getGaussPts().size2(),
611  w_ptr, FTensor::Number<DIM>());
612  auto t_diff_col = MakeB::getFTensor2SymmetricDiffBase(
613  col_data, baseColStorage, commonDataPtr->bBar, OP::getGaussPts().size2(),
614  w_ptr, FTensor::Number<DIM>());
615 
616  auto get_ftensor2_symmetric = [&](auto &storage, const auto gg,
617  const auto rr) {
619  &storage(gg, 6 * rr + 0), &storage(gg, 6 * rr + 1),
620  &storage(gg, 6 * rr + 2), &storage(gg, 6 * rr + 3),
621  &storage(gg, 6 * rr + 4), &storage(gg, 6 * rr + 5)};
622  };
623 
628 
629  auto t_Cep =
630  getFTensor4DdgFromMat<3, 3, 1>(commonDataPtr->materialTangentOperator);
631  const auto vol = OP::getMeasure();
632  auto t_w = OP::getFTensor0IntegrationWeight();
633  for (int gg = 0; gg != OP::nbIntegrationPts; ++gg) {
634  const double alpha = vol * t_w;
635  ++t_w;
636 
637  for (auto rr = 0; rr != OP::nbRows; ++rr) {
639  t_rowCep(k, l) = t_Cep(i, j, k, l) * t_diff_row(i, j);
640  auto t_diff_col = get_ftensor2_symmetric(baseColStorage, gg, 0);
641  for (auto cc = 0; cc != OP::nbCols; ++cc) {
642  OP::locMat(rr, cc) += alpha * (t_rowCep(k, l) * t_diff_col(k, l));
643  ++t_diff_col;
644  }
645  ++t_diff_row;
646  }
647 
648  ++t_Cep;
649  }
651 }
652 
653 } // namespace ADOLCPlasticity
654 
655 #endif //__ADOLCPLASTICITY_HPP_
656 
ADOLCPlasticity::ClosestPointProjection::getActiveVariablesYH
VectorAdaptor getActiveVariablesYH()
Definition: ADOLCPlasticity.hpp:158
ADOLCPlasticity::CommonData::deltaGamma
vector< double > deltaGamma
Definition: ADOLCPlasticity.hpp:97
ADOLCPlasticity::ClosestPointProjection::playW
MoFEMErrorCode playW()
Definition: ClosetPointProjection.cpp:122
ADOLCPlasticity::ClosestPointProjection::playW_NoHessian
MoFEMErrorCode playW_NoHessian()
Definition: ClosetPointProjection.cpp:172
ADOLCPlasticity::ClosestPointProjection::getInternalVariables
VectorAdaptor getInternalVariables()
Definition: ADOLCPlasticity.hpp:155
ADOLCPlasticity::ClosestPointProjection::flowPotential
virtual MoFEMErrorCode flowPotential()=0
Set flow potential.
ADOLCPlasticity::ClosestPointProjection::partial2HSigmaFlux
MatrixDouble partial2HSigmaFlux
Definition: ADOLCPlasticity.hpp:315
ADOLCPlasticity::ClosestPointProjection::plasticStrain0
VectorDouble plasticStrain0
Definition: ADOLCPlasticity.hpp:147
ADOLCPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: ADOLCPlasticity.hpp:330
ADOLCPlasticity::getRawPtrOpCalculateStress< 3 >
ForcesAndSourcesCore::UserDataOperator * getRawPtrOpCalculateStress< 3 >(MoFEM::Interface &m_field, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr, bool calc_lhs)
Definition: ADOLCPlasticity.cpp:547
ADOLCPlasticity::ClosestPointProjection::playH_NoHessian
MoFEMErrorCode playH_NoHessian()
Definition: ClosetPointProjection.cpp:287
ADOLCPlasticity::Pip
boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > Pip
Definition: ADOLCPlasticity.hpp:399
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
ADOLCPlasticity::ClosestPointProjection::Ep
MatrixDouble Ep
Definition: ADOLCPlasticity.hpp:181
ADOLCPlasticity::ClosestPointProjection::deltaGamma
double deltaGamma
Definition: ADOLCPlasticity.hpp:180
ADOLCPlasticity::createTSUpdate
boost::shared_ptr< TSUpdate > createTSUpdate(std::string fe_name, boost::shared_ptr< FEMethod > fe_ptr)
ADOLCPlasticity::TSUpdate
Update internal fluxes (update history variables)
Definition: ADOLCPlasticity.hpp:508
ADOLCPlasticity::ClosestPointProjection::ClosestPointProjection
ClosestPointProjection()
Definition: ClosetPointProjection.cpp:15
ADOLCPlasticity::CommonData::getDefaultMaterialParameters
MoFEMErrorCode getDefaultMaterialParameters()
Definition: ADOLCPlasticity.hpp:109
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:596
ADOLCPlasticity::strain_to_voight_op
FTensor::Dg< double, 3, 6 > strain_to_voight_op()
Op convert strain tensor to Vight strain vector.
Definition: ADOLCPlasticity.hpp:39
ADOLCPlasticity::ClosestPointProjection::TypesTags
TypesTags
Definition: ADOLCPlasticity.hpp:168
ADOLCPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: ADOLCPlasticity.hpp:325
ADOLCPlasticity::ClosestPointProjection::evaluatePotentials
MoFEMErrorCode evaluatePotentials()
Definition: ClosetPointProjection.cpp:335
ADOLCPlasticity::ClosestPointProjection::A
SmartPetscObj< Mat > A
Definition: ADOLCPlasticity.hpp:317
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
ADOLCPlasticity::CommonData::getFTensor1StressAtGaussPts
auto getFTensor1StressAtGaussPts(int gg=0)
Definition: ADOLCPlasticity.hpp:64
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
ADOLCPlasticity::MakeB::getFTensor2SymmetricDiffBase
static FTensor::Tensor2_symmetric< FTensor::PackPtr< double *, 6 >, 3 > getFTensor2SymmetricDiffBase(DataForcesAndSourcesCore::EntData &data, MatrixDouble &storage, const bool b_bar, const int nb_integration_pts, double *w_ptr, FTensor::Number< 2 >)
Definition: ADOLCPlasticity.cpp:186
ADOLCPlasticity::ADOLCPlasticityIntegrators
Definition: ADOLCPlasticity.hpp:384
ADOLCPlasticity::CommonData::activeVariablesW
MatrixDouble activeVariablesW
Definition: ADOLCPlasticity.hpp:61
ADOLCPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: ADOLCPlasticity.hpp:324
ADOLCPlasticity::ClosestPointProjection::partialYFlux
VectorAdaptor partialYFlux
Definition: ADOLCPlasticity.hpp:300
ADOLCPlasticity::ClosestPointProjection::setParams
virtual MoFEMErrorCode setParams(short tag, bool &recalculate_elastic_tangent)
Set parameters for ADOL-C of constitutive relations.
Definition: ADOLCPlasticity.hpp:274
ADOLCPlasticity::ClosestPointProjection::partialWStrainPlasticStrain
MatrixDouble partialWStrainPlasticStrain
Definition: ADOLCPlasticity.hpp:296
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
ADOLCPlasticity::CommonData::gradientW
MatrixDouble gradientW
Definition: ADOLCPlasticity.hpp:62
ADOLCPlasticity::CommonData::plasticStrainSize
int plasticStrainSize
Definition: ADOLCPlasticity.hpp:100
ADOLCPlasticity::ClosestPointProjection::Cep
MatrixDouble Cep
Definition: ADOLCPlasticity.hpp:181
sdf.r
int r
Definition: sdf.py:8
ADOLCPlasticity::OpRhsImpl< DIM, GAUSS, AssemblyDomainEleOp >::baseStorage
MatrixDouble baseStorage
Store tensorial base functions.
Definition: ADOLCPlasticity.hpp:541
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
ADOLCPlasticity::ClosestPointProjection::dChi
SmartPetscObj< Vec > dChi
Definition: ADOLCPlasticity.hpp:320
ADOLCPlasticity::ClosestPointProjection::hessianW
MatrixDouble hessianW
Definition: ADOLCPlasticity.hpp:199
ADOLCPlasticity::CommonData::getGardAtGaussPtsPtr
auto getGardAtGaussPtsPtr()
Definition: ADOLCPlasticity.hpp:103
ADOLCPlasticity::CommonData::internalVariablesPtr
double * internalVariablesPtr
Definition: ADOLCPlasticity.hpp:99
ADOLCPlasticity::ClosestPointProjection::calculateR
MoFEMErrorCode calculateR(Vec R)
Definition: ClosetPointProjection.cpp:360
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:170
ADOLCPlasticity::ClosestPointProjection::playY
MoFEMErrorCode playY()
Definition: ClosetPointProjection.cpp:191
ADOLCPlasticity::ClosestPointProjection::recordY
MoFEMErrorCode recordY()
Record yield function.
Definition: ClosetPointProjection.cpp:70
ADOLCPlasticity::ClosestPointProjection::gradientW
VectorAdaptor gradientW
Definition: ADOLCPlasticity.hpp:174
FTensor::Number< 2 >
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
ADOLCPlasticity::ClosestPointProjection
Closest point projection algorithm.
Definition: ADOLCPlasticity.hpp:139
ADOLCPlasticity::ClosestPointProjection::H
@ H
Definition: ADOLCPlasticity.hpp:168
ADOLCPlasticity::ClosestPointProjection::integrationRule
boost::function< int(int, int, int)> integrationRule
Definition: ADOLCPlasticity.hpp:142
ADOLCPlasticity::ClosestPointProjection::getPlasticStrain
VectorAdaptor getPlasticStrain()
Definition: ADOLCPlasticity.hpp:149
OpLhs
Definition: approx_sphere.cpp:134
AssemblyDomainEleOp
ADOLCPlasticity::ClosestPointProjection::playPotentials
MoFEMErrorCode playPotentials()
Definition: ClosetPointProjection.cpp:343
ADOLCPlasticity::OpRhsImpl
Assemble right hand side.
Definition: ADOLCPlasticity.hpp:371
MoFEM::getVectorAdaptor
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition: Templates.hpp:31
ADOLCPlasticity::CommonData::CommonData
CommonData()
Definition: ADOLCPlasticity.hpp:129
ADOLCPlasticity::CommonData::stressMatrix
MatrixDouble stressMatrix
Definition: ADOLCPlasticity.hpp:126
ADOLCPlasticity::CommonData::plasticStrainPtr
double * plasticStrainPtr
Definition: ADOLCPlasticity.hpp:98
ADOLCPlasticity::ClosestPointProjection::getStress
VectorAdaptor getStress()
Definition: ADOLCPlasticity.hpp:161
ADOLCPlasticity::opFactoryDomainLhs
MoFEMErrorCode opFactoryDomainLhs(MoFEM::Interface &m_field, string field_name, Pip &pip, std::string block_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr)
Assemble the left hand side, i.e. tangent matrix.
Definition: ADOLCPlasticity.hpp:449
ADOLCPlasticity::ClosestPointProjection::solveClosetProjection
MoFEMErrorCode solveClosetProjection()
Solve nonlinear system of equations to find stress, internal fluxes, and Lagrange plastic multiplier.
Definition: ClosetPointProjection.cpp:459
ADOLCPlasticity::getRawPtrOpCalculateStress
ForcesAndSourcesCore::UserDataOperator * getRawPtrOpCalculateStress(MoFEM::Interface &m_field, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr, bool calc_lhs)
Get opreator to calulate stress.
ADOLCPlasticity::ClosestPointProjection::getInternalFluxes
VectorAdaptor getInternalFluxes()
Definition: ADOLCPlasticity.hpp:164
ADOLCPlasticity::OpLhsImpl< DIM, GAUSS, AssemblyDomainEleOp >::baseRowStorage
MatrixDouble baseRowStorage
Store tensorial base functions.
Definition: ADOLCPlasticity.hpp:553
ADOLCPlasticity::CommonData::getPlasticStrain
auto getPlasticStrain(int gg=0)
Definition: ADOLCPlasticity.hpp:88
ADOLCPlasticity::ClosestPointProjection::h
double h
Definition: ADOLCPlasticity.hpp:178
ADOLCPlasticity::CommonData::internalVariablesSize
int internalVariablesSize
Definition: ADOLCPlasticity.hpp:101
ADOLCPlasticity::ClosestPointProjection::partial2HFlux
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
Definition: ADOLCPlasticity.hpp:312
ADOLCPlasticity::CommonData::getPlasticStrainMatrixPtr
boost::shared_ptr< MatrixDouble > getPlasticStrainMatrixPtr()
Definition: ADOLCPlasticity.hpp:121
EshelbianPlasticity::P
@ P
Definition: EshelbianContact.cpp:193
ADOLCPlasticity::ClosestPointProjection::playPotentials_NoHessian
MoFEMErrorCode playPotentials_NoHessian()
Definition: ClosetPointProjection.cpp:351
ADOLCPlasticity::opFactoryDomainRhs
MoFEMErrorCode opFactoryDomainRhs(MoFEM::Interface &m_field, string field_name, Pip &pip, std::string block_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr, Sev sev=Sev::inform)
Assemble the left hand side, i.e. tangent matrix.
Definition: ADOLCPlasticity.hpp:418
ADOLCPlasticity::ClosestPointProjection::partialYSigma
VectorAdaptor partialYSigma
Definition: ADOLCPlasticity.hpp:298
ADOLCPlasticity::voight_to_strain_op
FTensor::Dg< double, 3, 6 > voight_to_strain_op()
Op convert Vight strain vector to strain tensor.
Definition: ADOLCPlasticity.hpp:29
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:128
ADOLCPlasticity::CommonData::gradAtGaussPts
MatrixDouble gradAtGaussPts
Definition: ADOLCPlasticity.hpp:83
ADOLCPlasticity::CommonData::plasticStrainMatrix
MatrixDouble plasticStrainMatrix
Definition: ADOLCPlasticity.hpp:127
ADOLCPlasticity::ClosestPointProjection::Chi
SmartPetscObj< Vec > Chi
Definition: ADOLCPlasticity.hpp:319
ADOLCPlasticity::OpLhsImpl< DIM, GAUSS, AssemblyDomainEleOp >::baseColStorage
MatrixDouble baseColStorage
Store tensorial base functions.
Definition: ADOLCPlasticity.hpp:554
ADOLCPlasticity::ClosestPointProjection::R
SmartPetscObj< Vec > R
Definition: ADOLCPlasticity.hpp:318
ADOLCPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: ADOLCPlasticity.hpp:328
ADOLCPlasticity::OpLhsImpl
Assemble left hand side.
Definition: ADOLCPlasticity.hpp:382
ADOLCPlasticity::ClosestPointProjection::dataA
ublas::matrix< double, ublas::column_major > dataA
Definition: ADOLCPlasticity.hpp:321
ADOLCPlasticity::TSUpdate::postProcess
virtual MoFEMErrorCode postProcess(TS ts)=0
ADOLCPlasticity::ClosestPointProjection::partialHSigma
VectorAdaptor partialHSigma
Definition: ADOLCPlasticity.hpp:302
ADOLCPlasticity
Definition: ADOLCPlasticity.hpp:24
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
ADOLCPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: ADOLCPlasticity.hpp:327
ADOLCPlasticity::ClosestPointProjection::nbInternalVariables
int nbInternalVariables
Definition: ADOLCPlasticity.hpp:141
ADOLCPlasticity::ClosestPointProjection::activeVariablesW
VectorAdaptor activeVariablesW
Definition: ADOLCPlasticity.hpp:173
ADOLCPlasticity::ClosestPointProjection::tapesTags
std::array< int, LAST_TAPE > tapesTags
Definition: ADOLCPlasticity.hpp:169
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
ADOLCPlasticity::opFactoryDomainUpdate
MoFEMErrorCode opFactoryDomainUpdate(MoFEM::Interface &m_field, Pip &pip, std::string block_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr)
Push operators to update history variables.
ADOLCPlasticity::ADOLCPlasticityRes
MoFEMErrorCode ADOLCPlasticityRes(SNES snes, Vec chi, Vec r, void *ctx)
Internal SNES function used at integration points to calulate stress.
Definition: ClosetPointProjection.cpp:558
ADOLCPlasticity::ClosestPointProjection::partialHFlux
VectorAdaptor partialHFlux
Definition: ADOLCPlasticity.hpp:304
OpRhs
Definition: approx_sphere.cpp:68
adouble
ADOLCPlasticity::ClosestPointProjection::ADOLCPlasticityRes
friend MoFEMErrorCode ADOLCPlasticityRes(SNES, Vec, Vec, void *ctx)
Function executed by nested SENES at evaluationg residual.
Definition: ClosetPointProjection.cpp:558
ADOLCPlasticity::ClosestPointProjection::w
double w
Definition: ADOLCPlasticity.hpp:176
ADOLCPlasticity::opFactoryDomainUpdate< 2 >
MoFEMErrorCode opFactoryDomainUpdate< 2 >(MoFEM::Interface &m_field, Pip &pip, std::string block_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr)
Push operators to update history variables.
Definition: ADOLCPlasticity.cpp:586
ADOLCPlasticity::ClosestPointProjection::sNes
SmartPetscObj< SNES > sNes
Definition: ADOLCPlasticity.hpp:322
FTensor::Dg
Definition: Dg_value.hpp:9
ADOLCPlasticity::ClosestPointProjection::getTotalStrain
VectorAdaptor getTotalStrain()
Definition: ADOLCPlasticity.hpp:152
ADOLCPlasticity::OpRhsImpl< DIM, GAUSS, AssemblyDomainEleOp >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: ADOLCPlasticity.hpp:540
ADOLCPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: ADOLCPlasticity.hpp:326
ADOLCPlasticity::ClosestPointProjection::recordW
MoFEMErrorCode recordW()
Record strain energy.
Definition: ClosetPointProjection.cpp:40
ADOLCPlasticity::ClosestPointProjection::playH
MoFEMErrorCode playH()
Definition: ClosetPointProjection.cpp:230
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
ADOLCPlasticity::ClosestPointProjection::consistentTangent
MoFEMErrorCode consistentTangent()
Calculate consistent tangent matrix.
Definition: ClosetPointProjection.cpp:509
ADOLCPlasticity::MakeB
Calculate tensorial base functions. Apply bBar method when needed.
Definition: ADOLCPlasticity.hpp:517
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
ADOLCPlasticity::ClosestPointProjection::addMatBlockOps
virtual MoFEMErrorCode addMatBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, Sev sev)
Get model parameters from blocks.
Definition: ADOLCPlasticity.hpp:261
ADOLCPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: ADOLCPlasticity.hpp:329
ADOLCPlasticity::ClosestPointProjection::createMatAVecR
MoFEMErrorCode createMatAVecR()
Definition: ClosetPointProjection.cpp:314
ADOLCPlasticity::OpLhsImpl< DIM, GAUSS, AssemblyDomainEleOp >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: ADOLCPlasticity.hpp:552
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
ADOLCPlasticity::CommonData::getInternalVariables
auto getInternalVariables(int gg=0)
Definition: ADOLCPlasticity.hpp:92
ADOLCPlasticity::ClosestPointProjection::recordTapes
MoFEMErrorCode recordTapes()
Record tapes.
Definition: ClosetPointProjection.cpp:610
ADOLCPlasticity::voight_to_stress_op
FTensor::Dg< double, 3, 6 > voight_to_stress_op()
Op convert Vight stress vector to stress tensor.
Definition: ADOLCPlasticity.hpp:49
ADOLCPlasticity::ClosestPointProjection::calculateA
MoFEMErrorCode calculateA()
Definition: ClosetPointProjection.cpp:382
ADOLCPlasticity::ClosestPointProjection::LAST_TAPE
@ LAST_TAPE
Definition: ADOLCPlasticity.hpp:168
ADOLCPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: ADOLCPlasticity.hpp:327
ADOLCPlasticity::CommonData
common data used by volume elements
Definition: ADOLCPlasticity.hpp:59
ADOLCPlasticity::ClosestPointProjection::gradientY
VectorDouble gradientY
Definition: ADOLCPlasticity.hpp:203
ADOLCPlasticity::ClosestPointProjection::D
ublas::symmetric_matrix< double, ublas::lower > D
Definition: ADOLCPlasticity.hpp:182
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
ADOLCPlasticity::ClosestPointProjection::gradientH
VectorDouble gradientH
Definition: ADOLCPlasticity.hpp:207
ADOLCPlasticity::ClosestPointProjection::C
ublas::symmetric_matrix< double, ublas::lower > C
Definition: ADOLCPlasticity.hpp:182
AssemblyDomainEleOp
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp
Definition: tensor_divergence_operator.cpp:59
ADOLCPlasticity::ClosestPointProjection::snesCreate
MoFEMErrorCode snesCreate()
Create nested snes.
Definition: ClosetPointProjection.cpp:433
ADOLCPlasticity::ClosestPointProjection::ADOLCPlasticityJac
friend MoFEMErrorCode ADOLCPlasticityJac(SNES, Vec, Mat, Mat, void *ctx)
Function executed by nested SENES at evaluationg tangent matrix.
Definition: ClosetPointProjection.cpp:584
ADOLCPlasticity::CommonData::getStressMatrixPtr
boost::shared_ptr< MatrixDouble > getStressMatrixPtr()
Definition: ADOLCPlasticity.hpp:118
ADOLCPlasticity::getRawPtrOpCalculateStress< 2 >
ForcesAndSourcesCore::UserDataOperator * getRawPtrOpCalculateStress< 2 >(MoFEM::Interface &m_field, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr, bool calc_lhs)
Definition: ADOLCPlasticity.cpp:554
ADOLCPlasticity::ClosestPointProjection::partial2HSigma
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
Definition: ADOLCPlasticity.hpp:309
ADOLCPlasticity::ClosestPointProjection::freeHemholtzEnergy
virtual MoFEMErrorCode freeHemholtzEnergy()=0
Set Hemholtz energy.
ADOLCPlasticity::opFactoryDomainUpdate< 3 >
MoFEMErrorCode opFactoryDomainUpdate< 3 >(MoFEM::Interface &m_field, Pip &pip, std::string block_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< ClosestPointProjection > cp_ptr)
Push operators to update history variables.
Definition: ADOLCPlasticity.cpp:576
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
ADOLCPlasticity::ClosestPointProjection::playY_NoGradient
MoFEMErrorCode playY_NoGradient()
Definition: ClosetPointProjection.cpp:217
ADOLCPlasticity::CommonData::getFTensor1PlasticStrainAtGaussPts
auto getFTensor1PlasticStrainAtGaussPts(int gg=0)
Definition: ADOLCPlasticity.hpp:74
OP
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
ADOLCPlasticity::ClosestPointProjection::internalVariables0
VectorDouble internalVariables0
Definition: ADOLCPlasticity.hpp:146
ADOLCPlasticity::CommonData::bBar
bool bBar
Definition: ADOLCPlasticity.hpp:107
ADOLCPlasticity::ADOLCPlasticityIntegrators::Assembly
Definition: ADOLCPlasticity.hpp:386
ADOLCPlasticity::ClosestPointProjection::Cp
MatrixDouble Cp
Definition: ADOLCPlasticity.hpp:181
ADOLCPlasticity::CommonData::materialTangentOperator
MatrixDouble materialTangentOperator
Definition: ADOLCPlasticity.hpp:86
ADOLCPlasticity::ClosestPointProjection::Y
@ Y
Definition: ADOLCPlasticity.hpp:168
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
ADOLCPlasticity::ClosestPointProjection::yieldFunction
virtual MoFEMErrorCode yieldFunction()=0
Set yield function.
ADOLCPlasticity::ClosestPointProjection::hessianH
MatrixDouble hessianH
Definition: ADOLCPlasticity.hpp:208
ADOLCPlasticity::ClosestPointProjection::y
double y
Definition: ADOLCPlasticity.hpp:177
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
ADOLCPlasticity::ClosestPointProjection::recordH
MoFEMErrorCode recordH()
Record flow potential.
Definition: ClosetPointProjection.cpp:96
ADOLCPlasticity::ClosestPointProjection::W
@ W
Definition: ADOLCPlasticity.hpp:168
ADOLCPlasticity::ADOLCPlasticityJac
MoFEMErrorCode ADOLCPlasticityJac(SNES, Vec, Mat, Mat, void *ctx)
Internal SNES function used at integration points to calulate tangent matrix.
Definition: ClosetPointProjection.cpp:584