v0.15.5
Loading...
Searching...
No Matches
AdolCElastic.hpp
Go to the documentation of this file.
1/**
2 * @file AdolCElastic.hpp
3 * @brief
4 * @version 0.1
5 * @date 2026-04-02
6 *
7 * @copyright Copyright (c) 2026
8 *
9 */
10
11#ifndef ADOLC_OPS_ELEASTIC_HPP
12#define ADOLC_OPS_ELEASTIC_HPP
13
14namespace AdolCOps {
15
16struct ELASTICITY {
17 ELASTICITY() = delete;
18 // Behaviors
19 struct META {};
20 struct NEOHOOKEAN {};
22};
23
24/** \name Elasticity Operator Factory */
25/**@{*/
26template <int MODEL_TYPE> struct OpMaterialFactory<ELASTICITY, MODEL_TYPE> {
28
29 template <AssemblyType A, IntegrationType I, typename DomainEleOp>
30 static MoFEMErrorCode
32 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
33 std::string field_name,
34 boost::shared_ptr<PhysicalEquations> physical_equations_ptr) {
36 constexpr int DIM = (MODEL_TYPE == MODEL_3D) ? 3 : 2;
37 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
38 A>::template LinearForm<I>;
39 using OpInternalForce = typename B::template OpGradTimesTensor<1, DIM, DIM>;
40 auto m_grad =
41 physical_equations_ptr->adolcDataPtr->getCommonDataPtr("grad");
42 auto m_P = physical_equations_ptr->adolcDataPtr->getCommonDataPtr("P");
43 pip.push_back(
44 new OpCalculateVectorFieldGradient<DIM, DIM>(field_name, m_grad));
45 pip.push_back(
46 physical_equations_ptr->createOp(physical_equations_ptr, true, false));
47 pip.push_back(new OpInternalForce(field_name, m_P));
49 }
50
51 template <AssemblyType A, IntegrationType I, typename DomainEleOp>
52 static MoFEMErrorCode
54 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
55 std::string field_name,
56 boost::shared_ptr<PhysicalEquations> physical_equations_ptr) {
58 constexpr int DIM = (MODEL_TYPE == MODEL_3D) ? 3 : 2;
59 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
60 A>::template BiLinearForm<I>;
61 using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
62 auto m_grad =
63 physical_equations_ptr->adolcDataPtr->getCommonDataPtr("grad");
64 auto m_P_dF =
65 physical_equations_ptr->adolcDataPtr->getCommonDataPtr("P_dF");
66 pip.push_back(
67 new OpCalculateVectorFieldGradient<DIM, DIM>(field_name, m_grad));
68 pip.push_back(
69 physical_equations_ptr->createOp(physical_equations_ptr, false, true));
70 pip.push_back(new OpKPiola(field_name, field_name, m_P_dF));
72 }
73
74 static MoFEMErrorCode opPostProcFactory(
75 MoFEM::Interface &m_field,
76 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
77 std::string field_name,
78 boost::shared_ptr<PhysicalEquations> physical_equations_ptr) {
80 constexpr int DIM = (MODEL_TYPE == MODEL_3D) ? 3 : 2;
81 auto m_grad =
82 physical_equations_ptr->adolcDataPtr->getCommonDataPtr("grad");
83 auto m_P = physical_equations_ptr->adolcDataPtr->getCommonDataPtr("P");
84 pip.push_back(
85 new OpCalculateVectorFieldGradient<DIM, DIM>(field_name, m_grad));
86 pip.push_back(
87 physical_equations_ptr->createOp(physical_equations_ptr, true, false));
89 }
90};
91/**@}*/
92
93/** \name Elasticity Material Models Creators */
94/**@{*/
95
96template <>
97boost::shared_ptr<PhysicalEquations>
99 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
100
101template <>
102boost::shared_ptr<PhysicalEquations>
104 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
105
106template <>
107boost::shared_ptr<PhysicalEquations>
109 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
110
111template <>
112boost::shared_ptr<PhysicalEquations>
114 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
115
123
124template <>
125boost::shared_ptr<PhysicalEquations>
127 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
128
129/**@}*/
130
131} // namespace AdolCOps
132
133#endif // ADOLC_OPS_ELEASTIC_HPP
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::NEOHOOKEAN, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::VOLUMELENGTHQUALITY, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::META, MODEL_2D_PLANE_STRAIN >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
@ LASTOP_VOLUMELENGTHQUALITYTYPE
@ BARRIER_AND_CHANGE_QUALITY
@ BARRIER_AND_CHANGE_QUALITY_SCALED_BY_VOLUME
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::META, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::NEOHOOKEAN, MODEL_2D_PLANE_STRAIN >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
constexpr AssemblyType A
constexpr auto field_name
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition seepage.cpp:63
static MoFEMErrorCode opLhsFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr)
static MoFEMErrorCode opPostProcFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr)
static MoFEMErrorCode opRhsFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr)
Deprecated interface functions.