v0.15.0
Loading...
Searching...
No Matches
EshelbianPlasticity.hpp
Go to the documentation of this file.
1/**
2 * \file EshelbianPlasticity.hpp
3 * \brief Eshelbian plasticity interface
4 *
5 * \brief Problem implementation for mix element for large-strain elasticity
6 *
7 * For reference on mixed formulation see: \cite gopalakrishnan2012second and
8 * \cite cockburn2010new
9 *
10 * \todo Implementation of plasticity
11 */
12
13constexpr int SPACE_DIM = 3;
14// constexpr auto A = AssemblyType::BLOCK_SCHUR;
15constexpr auto A = AssemblyType::BLOCK_MAT;
16
17#ifndef __ESHELBIAN_PLASTICITY_HPP__
18 #define __ESHELBIAN_PLASTICITY_HPP__
19
20 #ifdef ENABLE_PYTHON_BINDING
21 #include <boost/python.hpp>
22 #include <boost/python/def.hpp>
23 #include <boost/python/numpy.hpp>
24namespace bp = boost::python;
25namespace np = boost::python::numpy;
26 #endif
27namespace EshelbianPlasticity {
28
30
31 using CachePhi = boost::tuple<int, int, MatrixDouble>;
32
33 CGGUserPolynomialBase(boost::shared_ptr<CachePhi> cache_phi = nullptr);
35
36 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
37 BaseFunctionUnknownInterface **iface) const;
38 MoFEMErrorCode getValue(MatrixDouble &pts,
39 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
40
41private:
42 MatrixDouble shapeFun;
43 boost::shared_ptr<CachePhi> cachePhiPtr;
44
45 MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts);
46};
47
48struct ContactTree;
49
54
55using MatrixPtr = boost::shared_ptr<MatrixDouble>;
56using VectorPtr = boost::shared_ptr<VectorDouble>;
57
58using EntData = EntitiesFieldData::EntData;
59using UserDataOperator = ForcesAndSourcesCore::UserDataOperator;
60using VolUserDataOperator = VolumeElementForcesAndSourcesCore::UserDataOperator;
61using FaceUserDataOperator = FaceElementForcesAndSourcesCore::UserDataOperator;
62using EleOnSide = PipelineManager::ElementsAndOpsByDim<SPACE_DIM>::FaceSideEle;
63using SideEleOp = EleOnSide::UserDataOperator;
64
66
67struct AnalyticalExprPython;
68
70 : public boost::enable_shared_from_this<DataAtIntegrationPts> {
71
72 MatrixDouble approxPAtPts;
73 MatrixDouble approxSigmaAtPts;
74 MatrixDouble divPAtPts;
75 MatrixDouble divSigmaAtPts;
76 MatrixDouble wL2AtPts;
77 MatrixDouble wL2DotAtPts;
78 MatrixDouble wL2DotDotAtPts;
80 MatrixDouble stretchTensorAtPts;
81 VectorDouble jacobianAtPts;
82 MatrixDouble tractionAtPts;
83
89 MatrixDouble rotAxisAtPts;
90 MatrixDouble rotAxisDotAtPts;
91 MatrixDouble rotAxisGradDotAtPts;
92 MatrixDouble rotAxis0AtPts;
93 MatrixDouble WAtPts;
94 MatrixDouble W0AtPts;
95 MatrixDouble GAtPts;
96 MatrixDouble G0AtPts;
97 MatrixDouble wH1AtPts;
98 MatrixDouble XH1AtPts;
99 MatrixDouble contactL2AtPts;
100 MatrixDouble wGradH1AtPts;
101 MatrixDouble logStretch2H1AtPts;
103
104 MatrixDouble hAtPts;
105 MatrixDouble hdOmegaAtPts;
106 MatrixDouble hdLogStretchAtPts;
107 MatrixDouble leviKirchhoffAtPts;
112 MatrixDouble adjointPdUAtPts;
114 MatrixDouble adjointPdUdPAtPts;
115 MatrixDouble rotMatAtPts;
116 MatrixDouble PAtPts;
117 MatrixDouble SigmaAtPts;
118 VectorDouble energyAtPts; //< this is density of energy at integration points
119 MatrixDouble flowL2AtPts;
120
121 MatrixDouble varRotAxis;
122 MatrixDouble varLogStreach;
123 MatrixDouble varPiola;
124 MatrixDouble varDivPiola;
125 MatrixDouble varWL2;
126
127 MatrixDouble P_du;
128
129 MatrixDouble eigenVals;
130 MatrixDouble eigenVecs;
131 VectorInt nbUniq;
132
133 MatrixDouble matD;
134 MatrixDouble matAxiatorD;
135 MatrixDouble matDeviatorD;
136 MatrixDouble matInvD;
137
140
141 MatrixDouble facePiolaAtPts;
145
146 double mu;
147 double lambda;
148 double piolaScale = 1.;
149 double faceEnergy = 0.;
150
151 inline auto getPiolaScalePtr() {
152 return boost::shared_ptr<double>(shared_from_this(), &piolaScale);
153 }
154
156 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
158 }
160 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &approxPAtPts);
161 }
162
164 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &divPAtPts);
165 }
166
168 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &divSigmaAtPts);
169 }
170
172 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wL2AtPts);
173 }
174
176 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wL2DotAtPts);
177 }
178
180 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wL2DotDotAtPts);
181 }
182
184 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
186 }
187
189 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
191 }
192
194 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
196 }
197
199 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
201 }
202
204 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &rotAxisAtPts);
205 }
206
208 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &rotAxis0AtPts);
209 }
210
212 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
214 }
215
217 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
219 }
220
222 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &GAtPts);
223 }
224
226 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &GAtPts);
227 }
228
230 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matD);
231 }
232
234 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matInvD);
235 }
236
238 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matAxiatorD);
239 }
240
242 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matDeviatorD);
243 }
244
246 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wH1AtPts);
247 }
248
250 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &XH1AtPts);
251 }
252
254 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &contactL2AtPts);
255 }
256
258 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wGradH1AtPts);
259 }
260
262 return boost::shared_ptr<VectorDouble>(shared_from_this(), &jacobianAtPts);
263 }
264
266 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
268 };
269
271 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varRotAxis);
272 }
273
275 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varLogStreach);
276 }
277
279 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varPiola);
280 }
281
283 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varDivPiola);
284 }
285
287 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varWL2);
288 }
289
291 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
293 }
294
296 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &facePiolaAtPts);
297 }
298
300 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
302 }
303
304 boost::shared_ptr<PhysicalEquations> physicsPtr;
305};
306
307struct OpJacobian;
308
310
317
319
322
324 PhysicalEquations(const int size_active, const int size_dependent)
325 : activeVariables(size_active, 0),
326 dependentVariablesPiola(size_dependent, 0),
327 dependentVariablesPiolaDirevatives(size_dependent * size_active, 0) {}
328 virtual ~PhysicalEquations() = default;
329
330 virtual MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h) = 0;
331
332 virtual OpJacobian *
333 returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
334 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
335 boost::shared_ptr<PhysicalEquations> physics_ptr) = 0;
336
337 virtual VolUserDataOperator *
338 returnOpSpatialPhysical(const std::string &field_name,
339 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
340 const double alpha_u);
341
343 std::string row_field, std::string col_field,
344 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha);
345
346 virtual VolUserDataOperator *
347 returnOpCalculateEnergy(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
348 boost::shared_ptr<double> total_energy_ptr);
349
351 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
352 boost::shared_ptr<PhysicalEquations> physics_ptr);
353
355 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
356 boost::shared_ptr<PhysicalEquations> physics_ptr);
357
358 virtual VolUserDataOperator *
359 returnOpSetScale(boost::shared_ptr<double> scale_ptr,
360 boost::shared_ptr<PhysicalEquations> physics_ptr);
361
362 std::vector<double> activeVariables;
363 std::vector<double> dependentVariablesPiola;
365
366 /** \name Active variables */
367
368 /**@{*/
369
370 template <int S>
371 inline static DTensor2Ptr get_VecTensor2(std::vector<double> &v) {
372 return DTensor2Ptr(&v[S + 0], &v[S + 1], &v[S + 2], &v[S + 3], &v[S + 4],
373 &v[S + 5], &v[S + 6], &v[S + 7], &v[S + 8]);
374 }
375
376 template <int S>
377 inline static DTensor0Ptr get_VecTensor0(std::vector<double> &v) {
378 return DTensor0Ptr(&v[S + 0]);
379 }
380
381 template <int S0>
382 inline static DTensor3Ptr get_vecTensor3(std::vector<double> &v,
383 const int nba) {
384
385 const int A00 = nba * 0 + S0;
386 const int A01 = nba * 1 + S0;
387 const int A02 = nba * 2 + S0;
388 const int A10 = nba * 3 + S0;
389 const int A11 = nba * 4 + S0;
390 const int A12 = nba * 5 + S0;
391 const int A20 = nba * 6 + S0;
392 const int A21 = nba * 7 + S0;
393 const int A22 = nba * 8 + S0;
394
395 return DTensor3Ptr(
396
397 &v[A00 + 0], &v[A00 + 1], &v[A00 + 2], &v[A01 + 0], &v[A01 + 1],
398 &v[A01 + 2], &v[A02 + 0], &v[A02 + 1], &v[A02 + 2],
399
400 &v[A10 + 0], &v[A10 + 1], &v[A10 + 2], &v[A11 + 0], &v[A11 + 1],
401 &v[A11 + 2], &v[A12 + 0], &v[A12 + 1], &v[A12 + 2],
402
403 &v[A20 + 0], &v[A20 + 1], &v[A20 + 2], &v[A21 + 0], &v[A21 + 1],
404 &v[A21 + 2], &v[A22 + 0], &v[A22 + 1], &v[A22 + 2]
405
406 );
407 }
408
409 /**@}*/
410
411 /** \name Active variables */
412
413 /**@{*/
414
415 inline DTensor2Ptr get_h() { return get_VecTensor2<0>(activeVariables); }
416
417 /**@}*/
418
419 /** \name Dependent variables */
420
421 /**@{*/
422
424 return get_VecTensor2<0>(dependentVariablesPiola);
425 }
426
427 /**@}*/
428
429 /** \name Derivatives of dependent variables */
430
431 /**@{*/
432
434 return get_vecTensor3<0>(dependentVariablesPiolaDirevatives,
435 activeVariables.size());
436 }
438 return get_vecTensor3<3>(dependentVariablesPiolaDirevatives,
439 activeVariables.size());
440 }
442 return get_vecTensor3<6>(dependentVariablesPiolaDirevatives,
443 activeVariables.size());
444 }
445
446 /**@}*/
447};
448
449struct BcDisp {
450 BcDisp(std::string name, std::vector<double> attr, Range faces);
451 std::string blockName;
453 VectorDouble3 vals;
454 VectorInt3 flags;
455};
456using BcDispVec = std::vector<BcDisp>;
457
458struct BcRot {
459 BcRot(std::string name, std::vector<double> attr, Range faces);
460 std::string blockName;
462 VectorDouble vals;
463 double theta;
464};
465using BcRotVec = std::vector<BcRot>;
466
467typedef std::vector<Range> TractionFreeBc;
468
470 TractionBc(std::string name, std::vector<double> attr, Range faces);
471 std::string blockName;
473 VectorDouble3 vals;
474 VectorInt3 flags;
475};
476using TractionBcVec = std::vector<TractionBc>;
477
479 NormalDisplacementBc(std::string name, std::vector<double> attr, Range faces);
480 std::string blockName;
482 double val;
483};
484using NormalDisplacementBcVec = std::vector<NormalDisplacementBc>;
485
487 AnalyticalDisplacementBc(std::string name, std::vector<double> attr, Range faces);
488 std::string blockName;
490 VectorInt3 flags;
491};
492using AnalyticalDisplacementBcVec = std::vector<AnalyticalDisplacementBc>;
493
495 AnalyticalTractionBc(std::string name, std::vector<double> attr, Range faces);
496 std::string blockName;
498 VectorInt3 flags;
499};
500using AnalyticalTractionBcVec = std::vector<AnalyticalTractionBc>;
501
503 PressureBc(std::string name, std::vector<double> attr, Range faces);
504 std::string blockName;
506 double val;
507};
508using PressureBcVec = std::vector<PressureBc>;
509
511 ExternalStrain(std::string name, std::vector<double> attr, Range ents);
512 std::string blockName;
514 double val;
515};
516using ExternalStrainVec = std::vector<ExternalStrain>;
517
518#ifdef ENABLE_PYTHON_BINDING
519struct AnalyticalExprPython {
520 AnalyticalExprPython() = default;
521 virtual ~AnalyticalExprPython() = default;
522
523 MoFEMErrorCode analyticalExprInit(const std::string py_file);
524 MoFEMErrorCode evalAnalyticalDisp(double delta_t, double t, np::ndarray x,
525 np::ndarray y, np::ndarray z, const std::string& block_name,
526 np::ndarray &analytical_expr);
527 MoFEMErrorCode evalAnalyticalTraction(double delta_t, double t, np::ndarray x,
528 np::ndarray y, np::ndarray z, const std::string& block_name,
529 np::ndarray &analytical_expr);
530
531 template <typename T>
532 inline std::vector<T>
533 py_list_to_std_vector(const boost::python::object &iterable) {
534 return std::vector<T>(boost::python::stl_input_iterator<T>(iterable),
535 boost::python::stl_input_iterator<T>());
536 }
537
538private:
539 bp::object mainNamespace;
540 bp::object analyticalDispFun;
541 bp::object analyticalTractionFun;
542};
543
544extern boost::weak_ptr<AnalyticalExprPython> AnalyticalExprPythonWeakPtr;
545#endif
546
547 #include "EshelbianCore.hpp"
548 #include "EshelbianOperators.hpp"
549
550} // namespace EshelbianPlasticity
551
552#endif //__ESHELBIAN_PLASTICITY_HPP__
constexpr int SPACE_DIM
constexpr auto A
const double v
phase velocity of light in medium (cm/ns)
boost::shared_ptr< MatrixDouble > MatrixPtr
std::vector< AnalyticalTractionBc > AnalyticalTractionBcVec
std::vector< TractionBc > TractionBcVec
std::vector< AnalyticalDisplacementBc > AnalyticalDisplacementBcVec
std::vector< PressureBc > PressureBcVec
std::vector< Range > TractionFreeBc
std::vector< ExternalStrain > ExternalStrainVec
ForcesAndSourcesCore::UserDataOperator UserDataOperator
std::vector< BcRot > BcRotVec
std::vector< NormalDisplacementBc > NormalDisplacementBcVec
FaceElementForcesAndSourcesCore::UserDataOperator FaceUserDataOperator
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
std::vector< BcDisp > BcDispVec
boost::shared_ptr< VectorDouble > VectorPtr
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
boost::tuple< int, int, MatrixDouble > CachePhi
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, BaseFunctionUnknownInterface **iface) const
boost::shared_ptr< PhysicalEquations > physicsPtr
virtual VolUserDataOperator * returnOpCalculateVarStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpCalculateEnergy(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< double > total_energy_ptr)
FTensor::Tensor3< double, 3, 3, 3 > DTensor3
FTensor::Tensor3< FTensor::PackPtr< double *, 1 >, 3, 3, 3 > DTensor3Ptr
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, 3, 3 > DTensor2Ptr
virtual VolUserDataOperator * returnOpCalculateStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
static DTensor3Ptr get_vecTensor3(std::vector< double > &v, const int nba)
FTensor::Tensor2< adouble, 3, 3 > ATensor2
FTensor::Tensor2< double, 3, 3 > DTensor2
PhysicalEquations(const int size_active, const int size_dependent)
virtual VolUserDataOperator * returnOpSetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
static DTensor2Ptr get_VecTensor2(std::vector< double > &v)
FTensor::Tensor1< adouble, 3 > ATensor1
FTensor::Tensor3< adouble, 3, 3, 3 > ATensor3
virtual MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h)=0
virtual 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)=0
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > DTensor0Ptr
virtual VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
static DTensor0Ptr get_VecTensor0(std::vector< double > &v)
Data on single entity (This is passed as argument to DataOperator::doWork)