v0.14.0
SmallStrainPlasticity.hpp
Go to the documentation of this file.
1 /** \file SmallStrainPlasticity.hpp
2  * \ingroup small_strain_plasticity
3  * \brief Operators and data structures for small strain plasticity
4  * \example SmallStrainPlasticity.hpp
5 
6  \section problem_formulation Problem formulation
7 
8  \subsection theory_small_strain Theory
9 
10  \subsubsection introduction Introduction
11 
12  In the following document, we define functions,
13  \f$\psi(\varepsilon^e,\alpha)\f$, \f$f(\sigma,\mathbf{q})\f$ and
14  \f$\Psi(\sigma,\mathbf{q})\f$ which are free Helmholtz energy, yield function
15  and flow potential, respectively. Note that free Helmholtz energy is defined
16  by stains and kinematic internal state variables, whereas yield function and
17  flow potential is defined by fluxes (stress tensor and fluxes work conjugated
18  to internal kinematic state variables). In this implantation, we assume that all
19  functions are convex. Also, is assumed that energy potential and flow
20  potential is positive and takes zero value at zero strain and zero internal
21  variable states. It can be shown that for such defined problem, the first law
22  of thermodynamics and non-negative dissipation of energy is automatically
23  satisfied, see \cite de2011computational \cite simo2006computational.
24  Moreover, in the following implementation is assumed that all functions are
25  smooth enough that appropriate derivatives could be calculated.
26 
27  In the following implementation, the backward-Euler difference scheme is
28  applied to solve set nonlinear system of equations with unilateral constrains.
29  Admissible stress state is enforced by Lagrange multiplier method and
30  (Khun-Tucker) loading/unloading conditions, see \cite simo2006computational.
31  That lead to well known \em Closet \em point \em projection algorithm.
32 
33  Following \cite simo2006computational can be shown that plastic loading
34  unloading conditions can be exclusively characterised in terms of trial
35  state
36  \f[
37  \varepsilon^{e,\textrm{trial}}_{n+1} := \varepsilon_{n+1} - \varepsilon^p_n
38  \f]
39  \f[
40  \alpha^\textrm{trial}_{n+1} := \alpha_n
41  \f]
42  \f[
43  \sigma^\textrm{trial}_{n+1} := \frac{\partial \psi}{\partial \varepsilon^e}(\varepsilon^{e,\textrm{trial}}_{n+1},\alpha_n)
44  \f]
45  \f[
46  \mathbf{q}^\textrm{trial}_{n+1} := \frac{\partial \psi}{\partial \alpha}(\varepsilon^{e,\textrm{trial}}_{n+1},\alpha_n)
47  \f]
48  \f[
49  f^\textrm{trial}_{n+1} := f(\sigma^\textrm{trial}_{n+1},\mathbf{q}_n)
50  \f]
51  then loading/unloading conditions are given
52  \f[
53  \begin{split}
54  f^\textrm{trial}_{n+1} < 0,\quad\textrm{elastic step}\,\to\,\Delta\gamma_{n+1} = 0 \\
55  f^\textrm{trial}_{n+1} \geq 0,\quad\textrm{plastic step}\,\to\,\Delta\gamma_{n+1} > 0
56  \end{split}
57  \f]
58  If first condition is true \f$f^\textrm{trial}_{n+1} < 0\f$, trial stress is solution
59  at integration point. If second condition is true, plastic strains and values
60  of internal state variables are calculated using \em Closet \em point \em
61  projection algorithm.
62 
63  \subsection cloaset_point_projection Closet point projection algorithm
64 
65  \subsubsection elastic_strain Additive strain decomposition
66 
67  \f[
68  \varepsilon^e_{n+1} := \varepsilon_{n+1} - \varepsilon^p_{n+1}
69  \f]
70 
71  \subsubsection constitutive_equation Constitutive relations
72 
73  \f[
74  \sigma_{n+1} := \frac{\partial \psi}{\partial \varepsilon^e}(\varepsilon^e_{n+1},\alpha_{n+1})
75  \f]
76 
77  \f[
78  \mathbf{q}_{n+1} := \frac{\partial \psi}{\partial \alpha}(\varepsilon^e_{n+1},\alpha_{n+1})
79  \f]
80 
81  \subsubsection constrains Constrains
82 
83  \f[
84  f_{n+1}:=f(\sigma_{n+1},\mathbf{q}_{n+1})
85  \f]
86 
87  \subsubsection evolution Evolution functions
88 
89  Flow vector
90  \f[
91  \mathbf{n}_{n+1}:=\frac{\partial \Psi}{\partial \sigma}(\sigma_{n+1},\mathbf{q}_{n+1})
92  \f]
93 
94  Hardening law:
95  \f[
96  \mathbf{h}_{n+1}:=-\frac{\partial \Psi}{\partial \mathbf{q}}(\sigma_{n+1},\mathbf{q}_{n+1})
97  \f]
98 
99  \subsubsection residual Residual
100 
101  \f[
102  \mathbf{R}_{n+1}:=
103  \left\{
104  \begin{array}{c}
105  -\varepsilon^p_{n+1}+\varepsilon^p_n \\
106  -\alpha_{n+1}+\alpha_n \\
107  f_{n+1}
108  \end{array}
109  \right\}
110  +
111  \Delta\gamma_{n+1}
112  \left\{
113  \begin{array}{c}
114  \mathbf{n}_{n+1} \\
115  \mathbf{h}_{n+1} \\
116  0
117  \end{array}
118  \right\}
119  \f]
120 
121  \f[
122  \left[
123  \begin{array}{ccc}
124  -1+\gamma
125  \left(
126  \frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
127  +
128  \frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial\alpha\varepsilon^p}
129  \right)
130  &
131  \gamma
132  \left(
133  \frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
134  +
135  \frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \alpha}
136  \right)
137  &
138  \mathbf{n}
139  \\
140  \gamma
141  \left(
142  \frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
143  +
144  \frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial\alpha\partial\varepsilon^p}
145  \right)
146  &
147  -1+\gamma
148  \left(
149  \frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
150  +
151  \frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \alpha}
152  \right)
153  &
154  \mathbf{h}\\
155  \frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
156  +
157  \frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha\partial\varepsilon^p}
158  &
159  \frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
160  +
161  \frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial\alpha}
162  & 0
163  \end{array}
164  \right]_n
165  \left\{
166  \begin{array}{c}
167  \delta \varepsilon^p \\
168  \delta \alpha \\
169  \delta \gamma
170  \end{array}
171  \right\}_{n+1}
172  =
173  \mathbf{R}_n
174  \f]
175  where
176  \f$\varepsilon^p_{n+1}=\varepsilon^p_n + \delta\varepsilon^p_{n+1}\f$,
177  \f$\alpha^p_{n+1}=\alpha_n + \delta\alpha_{n+1}\f$
178  and
179  \f$\Delta\gamma_{n+1} = \Delta\gamma_n+\delta\gamma_{n+1}\f$.
180  Usually it is assumed free energy split, i.e.
181  \f[
182  \psi(\varepsilon^e,\alpha) = \psi^{\varepsilon}(\varepsilon^e) + \psi^\alpha(\alpha)
183  \f]
184  in that case linearized residual takes the form
185  \f[
186  \left[
187  \begin{array}{ccc}
188  -1+\Delta\gamma
189  \frac{\partial\mathbf{n}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
190  &
191  \Delta\gamma
192  \frac{\partial\mathbf{n}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
193  &
194  \mathbf{n}
195  \\
196  \Delta\gamma
197  \frac{\partial\mathbf{h}}{\partial\sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
198  &
199  -1+\Delta\gamma
200  \frac{\partial\mathbf{h}}{\partial\mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
201  &
202  \mathbf{h}\\
203  \frac{\partial f}{\partial \sigma}\frac{\partial^2 \psi}{\partial \varepsilon^e \partial \varepsilon^p}
204  &
205  \frac{\partial f}{\partial \mathbf{q}}\frac{\partial^2 \psi}{\partial \alpha^2}
206  & 0
207  \end{array}
208  \right]_n
209  \left\{
210  \begin{array}{c}
211  \delta \varepsilon^p \\
212  \delta \alpha \\
213  \delta \gamma
214  \end{array}
215  \right\}_{n+1}
216  =
217  -\mathbf{R}_n
218  \f]
219 
220  \f[
221  \chi_{n+1} :=
222  \left\{
223  \begin{array}{c}
224  \delta \varepsilon^p \\
225  \delta \alpha \\
226  \delta \gamma
227  \end{array}
228  \right\}_{n+1}
229  =
230  -
231  \mathbf{A}^{-1}
232  \mathbf{R}_n
233  \f]
234 
235  \subsubsection small_strain_algebraic_tangent Algorithmic tangent matrix
236 
237  \f[
238  \delta \sigma := \mathbf{C}^e(\delta \varepsilon - \delta \varepsilon^p)
239  \f]
240  where \f$\mathbf{C}^e = \frac{\partial^2 \psi^e}{\partial {\varepsilon^e}^2}\f$ and
241  \f$\delta \varepsilon = \textrm{sym}[ \nabla \delta \mathbf{u}]\f$
242 
243  \f[
244  \left\{
245  \begin{array}{c}
246  \mathbf{E}^p \\
247  \mathbf{L}^p \\
248  \mathbf{Y}^p
249  \end{array}
250  \right\}
251  \delta \varepsilon
252  =
253  -
254  \mathbf{A}^{-1}
255  \left\{
256  \begin{array}{c}
257  \Delta\gamma \frac{\partial \mathbf{n}}{\partial \sigma}\mathbf{C}^e \\
258  \Delta\gamma \frac{\partial \mathbf{h}}{\partial \sigma}\mathbf{C}^e \\
259  \frac{\partial f}{\partial \sigma}\mathbf{C}^e
260  \end{array}
261  \right\}
262  \delta \varepsilon
263  \f]
264 
265  \f[
266  \delta \varepsilon^p = \mathbf{E}^p \delta \varepsilon
267  \f]
268 
269  \f[
270  \delta \alpha = \mathbf{L}^p \delta \varepsilon
271  \f]
272 
273  \f[
274  \delta \gamma = \mathbf{Y}^p \delta \varepsilon
275  \f]
276 
277  \f[
278  \delta \sigma = \mathbf{C}^e(\mathbf{1} - \mathbf{E}^p)\delta \varepsilon
279  \f]
280 
281  \f[
282  \delta \sigma = \mathbf{C}^{ep}\delta \varepsilon
283  \f]
284  where \f$\mathbf{C}^{ep} = \mathbf{C}^e-\mathbf{C}^e\mathbf{E}^p\f$
285 
286 For more details and easy access to all source files related to this see group
287 \ref small_strain_plasticity.
288 
289  */
290 
291 /*
292  * This file is part of MoFEM.
293  * MoFEM is free software: you can redistribute it and/or modify it under
294  * the terms of the GNU Lesser General Public License as published by the
295  * Free Software Foundation, either version 3 of the License, or (at your
296  * option) any later version.
297  *
298  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
299  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
300  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
301  * License for more details.
302  *
303  * You should have received a copy of the GNU Lesser General Public
304  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
305 
306 #ifndef __SMALL_STRAIN_PLASTICITY_HPP
307 #define __SMALL_STRAIN_PLASTICITY_HPP
308 
309 #ifndef WITH_ADOL_C
310  #error "MoFEM need to be compiled with ADOL-C"
311 #endif
312 
313 PetscErrorCode SmallStrainPlasticityfRes(SNES snes,Vec chi,Vec r,void *ctx);
314 PetscErrorCode SmallStrainPlasticityfJac(SNES,Vec,Mat,Mat,void *ctx);
315 
316 /** \brief Small strain finite element implementation
317 * \ingroup small_strain_plasticity
318 
319 */
321 
322 
323 
324 
326 
329 
330  /// \brief definition of volume element
332 
336  addToRule(0) { }
337 
338  /** \brief it is used to calculate nb. of Gauss integration points
339  *
340  * for more details pleas look
341  * Reference:
342  *
343  * Albert Nijenhuis, Herbert Wilf,
344  * Combinatorial Algorithms for Computers and Calculators,
345  * Second Edition,
346  * Academic Press, 1978,
347  * ISBN: 0-12-519260-6,
348  * LC: QA164.N54.
349  *
350  * More details about algorithm
351  * http://people.sc.fsu.edu/~jburkardt/cpp_src/gm_rule/gm_rule.html
352  **/
353  int getRule(int order) {
354  // cerr << "order " << order << endl;
355  return 2*(order-1)+addToRule;
356  }
357 
358  };
359 
360  MyVolumeFE feRhs; ///< calculate right hand side for tetrahedral elements
361  MyVolumeFE feLhs; //< calculate left hand side for tetrahedral elements
362  MyVolumeFE feUpdate; //< update history variables
363 
365  mField(m_field),
366  feRhs(m_field),
367  feLhs(m_field),
368  feUpdate(m_field) {
369  }
370 
371  PetscErrorCode createInternalVariablesTag();
372 
373  /** \brief common data used by volume elements
374  * \ingroup nonlinear_elastic_elem
375  */
376  struct CommonData {
377  map<string,vector<VectorDouble> > dataAtGaussPts;
378  map<string,vector<MatrixDouble> > gradAtGaussPts;
379  vector<VectorDouble> sTress;
380  vector<MatrixDouble> materialTangentOperator;
381  vector<VectorDouble> plasticStrain;
382  vector<VectorDouble> internalVariables;
383  vector<VectorDouble> internalFluxes;
384  vector<double> deltaGamma;
389  bool bBar;
391  bBar(true) {
392  }
393  };
395 
397  vector<VectorDouble > &valuesAtGaussPts;
398  vector<MatrixDouble > &gradientAtGaussPts;
399  const EntityType zeroAtType;
400 
401  OpGetDataAtGaussPts(const string field_name,
402  vector<VectorDouble > &values_at_gauss_pts,
403  vector<MatrixDouble > &gardient_at_gauss_pts
404  );
405 
406  /** \brief operator calculating field data and gradients
407  */
408  PetscErrorCode doWork(
409  int side,EntityType type,DataForcesAndSourcesCore::EntData &data
410  );
411  };
412 
414  OpGetCommonDataAtGaussPts(const string field_name,CommonData &common_data);
415  };
416 
417 
419 
421 
422  OpUpdate(
423  string field_name,
424  CommonData &common_data
425  );
426 
427  PetscErrorCode doWork(
428  int side,EntityType type,DataForcesAndSourcesCore::EntData &data
429  );
430 
431  };
432 
433  struct MakeB {
434  PetscErrorCode makeB(const MatrixAdaptor &diffN,MatrixDouble &B);
435  PetscErrorCode addVolumetricB(const MatrixAdaptor &diffN,MatrixDouble &volB,double alpha);
436  };
437 
439 
443 
444 
446  string field_name,
447  CommonData &common_data
448  );
449 
450  PetscErrorCode doWork(
451  int side,EntityType type,DataForcesAndSourcesCore::EntData &data
452  );
453 
454  };
455 
457 
462 
463 
465  string field_name,
466  CommonData &common_data
467  );
468 
469  PetscErrorCode doWork(
470  int row_side,int col_side,
471  EntityType row_type,EntityType col_type,
474  );
475 
476  };
477 
478  /** \brief Closest point projection algorithm
479  * \ingroup small_strain_plasticity
480 
481  */
483 
484 
490  double deltaGamma;
491  double tOl;
492  int gG;
493 
496  ublas::symmetric_matrix<double,ublas::lower> C,D;
502  ublas::symmetric_matrix<double,ublas::lower> partial2HSigma;
503  ublas::symmetric_matrix<double,ublas::lower> partial2HFlux;
505 
506  vector<int> tAgs;
507 
509 
511  virtual PetscErrorCode setActiveVariablesW();
512 
514  virtual PetscErrorCode setActiveVariablesYH();
515 
516  double w,y,h;
518  ublas::vector<adouble> a_sTress,a_internalFluxes;
520 
521  virtual PetscErrorCode rEcordW();
522 
523  virtual PetscErrorCode rEcordY();
524 
525  virtual PetscErrorCode rEcordH();
526 
529  virtual PetscErrorCode pLayW();
530  virtual PetscErrorCode pLayW_NoHessian();
531 
533  virtual PetscErrorCode pLayY();
534  virtual PetscErrorCode pLayY_NoGradient();
535 
538  virtual PetscErrorCode pLayH();
539 
540  Mat A;
541  ublas::matrix<double,ublas::column_major> dataA;
543  virtual PetscErrorCode createMatAVecR();
544 
545  virtual PetscErrorCode destroyMatAVecR();
546 
547  virtual PetscErrorCode evaluatePotentials();
548 
549  virtual PetscErrorCode cAlculateR(Vec R);
550 
551  virtual PetscErrorCode cAlculateA();
552 
553  friend PetscErrorCode SmallStrainPlasticityfRes(SNES,Vec,Vec,void *ctx);
554  friend PetscErrorCode SmallStrainPlasticityfJac(SNES,Vec,Mat,Mat,void *ctx);
555 
556  SNES sNes;
557  KSP kSp;
558  PC pC;
559  SNESLineSearch lineSearch;
560  PetscErrorCode snesCreate();
561 
562  PetscErrorCode snesDestroy();
563 
564  PetscErrorCode solveColasetProjection();
565 
568 
569  PetscErrorCode consistentTangent();
570 
571  virtual PetscErrorCode freeHemholtzEnergy() {
572  PetscFunctionBegin;
573  SETERRQ(
574  PETSC_COMM_SELF,
576  "Not implemented (overload class and implement this function)"
577  );
578  PetscFunctionReturn(0);
579  }
580 
581  virtual PetscErrorCode yieldFunction() {
582  PetscFunctionBegin;
583  SETERRQ(
584  PETSC_COMM_SELF,
586  "Not implemented (overload class and implement this function)"
587  );
588  PetscFunctionReturn(0);
589  }
590 
591  virtual PetscErrorCode flowPotential() {
592  PetscFunctionBegin;
593  SETERRQ(
594  PETSC_COMM_SELF,
596  "Not implemented (overload class and implement this function)"
597  );
598  PetscFunctionReturn(0);
599  }
600 
601  };
602 
604 
608 
609  bool initCp;
611 
612 
613 
615  MoFEM::Interface &m_field,
616  string field_name,
617  CommonData &common_data,
619  bool is_linear = true
620  ):
621  MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator(field_name,UserDataOperator::OPROW),
622  mField(m_field),
623  commonData(common_data),
624  cP(cp),
625  initCp(false),
626  isLinear(is_linear),
627  hessianWCalculated(false) {
628  }
629 
632  PetscErrorCode getTags();
633 
634  PetscErrorCode setTagsData(
635  const EntityHandle tet,const int nb_gauss_pts,const int nb_internal_variables
636  );
637 
638  PetscErrorCode doWork(
639  int side,EntityType type,DataForcesAndSourcesCore::EntData &data
640  );
641 
642  };
643 
644 };
645 
646 #endif //__SMALL_STRAIN_PLASTICITY_HPP
647 
648 /***************************************************************************//**
649  * \defgroup small_strain_plasticity Small Strain Plasticity
650  * \ingroup user_modules
651  ******************************************************************************/
652 
653 /***************************************************************************//**
654  * \defgroup user_modules User modules
655  ******************************************************************************/
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::OpGetCommonDataAtGaussPts
Definition: SmallStrainPlasticity.hpp:413
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::OpUpdate
Definition: SmallStrainPlasticity.hpp:418
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::ClosestPointProjection
Closest point projection algorithm.
Definition: SmallStrainPlasticity.hpp:482
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::MakeB
Definition: SmallStrainPlasticity.hpp:433
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
SmallStrainPlasticity::CommonData
common data used by volume elements
Definition: SmallStrainPlasticity.hpp:376
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::MyVolumeFE
definition of volume element
Definition: SmallStrainPlasticity.hpp:331
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::OpAssembleRhs
Definition: SmallStrainPlasticity.hpp:438
SmallStrainPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: SmallStrainPlasticity.hpp:519
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
SmallStrainPlasticity::OpGetDataAtGaussPts
Definition: SmallStrainPlasticity.hpp:396
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::OpAssembleLhs
Definition: SmallStrainPlasticity.hpp:456
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
SmallStrainPlasticity::OpCalculateStress
Definition: SmallStrainPlasticity.hpp:603