v0.15.4
Loading...
Searching...
No Matches
Classes | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes | List of all members
EshelbianPlasticity::HMHNeohookean Struct Reference
Inheritance diagram for EshelbianPlasticity::HMHNeohookean:
[legend]
Collaboration diagram for EshelbianPlasticity::HMHNeohookean:
[legend]

Classes

struct  BlockData
 
struct  OpGetScale
 
struct  OpJacobian
 
struct  OpSpatialPhysical
 
struct  OpSpatialPhysical_du_du
 
struct  OpSpatialPhysicalExternalStrain
 

Public Member Functions

 HMHNeohookean (MoFEM::Interface &m_field, const double c10, const double K)
 
auto getMaterialParameters (EntityHandle ent)
 
VolUserDataOperatorreturnOpSetScale (boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
 
UserDataOperatorreturnOpJacobian (const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
 
MoFEMErrorCode getOptions ()
 
MoFEMErrorCode extractBlockData (Sev sev)
 
MoFEMErrorCode extractBlockData (std::vector< const CubitMeshSets * > meshset_vec_ptr, Sev sev)
 
MoFEMErrorCode recordTape (const int tape, DTensor2Ptr *t_h_ptr)
 
virtual VolUserDataOperatorreturnOpSpatialPhysical (const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
 
virtual VolUserDataOperatorreturnOpSpatialPhysicalExternalStrain (const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< ExternalStrainVec > external_strain_vec_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
 
VolUserDataOperatorreturnOpSpatialPhysical_du_du (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
 
- Public Member Functions inherited from EshelbianPlasticity::PhysicalEquations
 PhysicalEquations ()=delete
 
 PhysicalEquations (const int size_active, const int size_dependent)
 
virtual ~PhysicalEquations ()=default
 
virtual VolUserDataOperatorreturnOpCalculateEnergy (boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< double > total_energy_ptr)
 
virtual VolUserDataOperatorreturnOpCalculateStretchFromStress (boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
 
virtual VolUserDataOperatorreturnOpCalculateVarStretchFromStress (boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
 
DTensor2Ptr get_P ()
 
DTensor3Ptr get_P_dh0 ()
 
DTensor3Ptr get_P_dh1 ()
 
DTensor3Ptr get_P_dh2 ()
 
DTensor2Ptr get_h ()
 

Static Public Member Functions

static double fun_neohookean (double c10, double v)
 Definition of Neo-hookean function.
 
static double fun_d_neohookean (double c10, double v)
 Definition of derivative of Neo-hookean function.
 
static double fun_neohookean_bulk (double K, double tr)
 Definition of axiator of Neo-hookean function.
 
static double fun_diff_neohookean_bulk (double K, double tr)
 Definition of derivative of axiator of Neo-hookean function.
 
- Static Public Member Functions inherited from EshelbianPlasticity::PhysicalEquations
template<int S>
static DTensor2Ptr get_VecTensor2 (std::vector< double > &v)
 
template<int S>
static DTensor0Ptr get_VecTensor0 (std::vector< double > &v)
 
template<int S0>
static DTensor3Ptr get_vecTensor3 (std::vector< double > &v, const int nba)
 

Static Public Attributes

static constexpr int numberOfActiveVariables = 9
 
static constexpr int numberOfDependentVariables = 9
 

Private Attributes

MoFEM::InterfacemField
 
double c10_default
 
double K_default
 
double alphaGradU
 
std::vector< BlockDatablockData
 
double eqScaling = 1.
 
PetscBool adolCOn = PETSC_FALSE
 
ATensor2 th
 
ATensor2 tF
 
adouble detF
 
ATensor2 tInvF
 
ATensor2 tCof
 
ATensor2 tP
 

Additional Inherited Members

- Public Types inherited from EshelbianPlasticity::PhysicalEquations
typedef FTensor::Tensor1< adouble, 3 > ATensor1
 
typedef FTensor::Tensor2< adouble, 3, 3 > ATensor2
 
typedef FTensor::Tensor3< adouble, 3, 3, 3 > ATensor3
 
typedef FTensor::Tensor1< double, 3 > DTensor1
 
typedef FTensor::Tensor2< double, 3, 3 > DTensor2
 
typedef FTensor::Tensor3< double, 3, 3, 3 > DTensor3
 
typedef FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > DTensor0Ptr
 
typedef FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, 3, 3 > DTensor2Ptr
 
typedef FTensor::Tensor3< FTensor::PackPtr< double *, 1 >, 3, 3, 3 > DTensor3Ptr
 
- Public Attributes inherited from EshelbianPlasticity::PhysicalEquations
std::vector< doubleactiveVariables
 
std::vector< doubledependentVariablesPiola
 
std::vector< doubledependentVariablesPiolaDirevatives
 

Detailed Description

Definition at line 13 of file HMHNeohookean.cpp.

Constructor & Destructor Documentation

◆ HMHNeohookean()

EshelbianPlasticity::HMHNeohookean::HMHNeohookean ( MoFEM::Interface m_field,
const double  c10,
const double  K 
)
inline

Definition at line 71 of file HMHNeohookean.cpp.

73 mField(m_field), c10_default(c10), K_default(K) {
74
75 CHK_THROW_MESSAGE(getOptions(), "get options failed");
77 "extract block data failed");
78
79#ifdef NEOHOOKEAN_SCALING
80 if (blockData.size()) {
81 double ave_K = 0;
82 for (auto b : blockData) {
83 ave_K += b.K;
84 }
85 eqScaling = ave_K / blockData.size();
86 }
87#endif
88
89 MOFEM_LOG("EP", Sev::inform) << "Neo-Hookean scaling = " << eqScaling;
90
91 if (adolCOn == PETSC_FALSE) {
94 "Stretch selector is not equal to LOG");
95 } else {
96 if (EshelbianCore::exponentBase != exp(1)) {
98 "Exponent base is not equal to exp(1)");
99 }
100 }
101 }
102 }
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MOFEM_LOG(channel, severity)
Log.
static enum StretchSelector stretchSelector
static double exponentBase
std::vector< BlockData > blockData
MoFEMErrorCode extractBlockData(Sev sev)
static constexpr int numberOfDependentVariables
static constexpr int numberOfActiveVariables

Member Function Documentation

◆ extractBlockData() [1/2]

MoFEMErrorCode EshelbianPlasticity::HMHNeohookean::extractBlockData ( Sev  sev)
inline

Definition at line 193 of file HMHNeohookean.cpp.

193 {
194 return extractBlockData(
195
196 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
197
198 (boost::format("%s(.*)") % "MAT_NEOHOOKEAN").str()
199
200 )),
201
202 sev);
203 }
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ extractBlockData() [2/2]

MoFEMErrorCode EshelbianPlasticity::HMHNeohookean::extractBlockData ( std::vector< const CubitMeshSets * >  meshset_vec_ptr,
Sev  sev 
)
inline

Definition at line 206 of file HMHNeohookean.cpp.

207 {
209
210 for (auto m : meshset_vec_ptr) {
211 MOFEM_LOG("EP", sev) << *m;
212 std::vector<double> block_data;
213 CHKERR m->getAttributes(block_data);
214 if (block_data.size() < 2) {
215 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
216 "Expected that block has atleast two attributes");
217 }
218 auto get_block_ents = [&]() {
219 Range ents;
220 CHKERR mField.get_moab().get_entities_by_handle(m->meshset, ents, true);
221 return ents;
222 };
223
224 double c10 = block_data[0];
225 double K = block_data[1];
226
227 blockData.push_back({c10, K, get_block_ents()});
228
229 MOFEM_LOG("EP", sev) << "MatBlock Neo-Hookean c10 = "
230 << blockData.back().c10
231 << " K = " << blockData.back().K << " nb ents. = "
232 << blockData.back().blockEnts.size();
233 }
235 }
#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()
#define CHKERR
Inline error check.
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0

◆ fun_d_neohookean()

static double EshelbianPlasticity::HMHNeohookean::fun_d_neohookean ( double  c10,
double  v 
)
inlinestatic

Definition of derivative of Neo-hookean function.

Parameters
c10
v
Returns
* double

Definition at line 38 of file HMHNeohookean.cpp.

38 {
39 return
40
41 (4.0 * c10) * EshelbianCore::d_f(2 * v);
42 }
const double v
phase velocity of light in medium (cm/ns)
static boost::function< double(const double)> d_f

◆ fun_diff_neohookean_bulk()

static double EshelbianPlasticity::HMHNeohookean::fun_diff_neohookean_bulk ( double  K,
double  tr 
)
inlinestatic

Definition of derivative of axiator of Neo-hookean function.

Parameters
K
tr
Returns
double

Definition at line 64 of file HMHNeohookean.cpp.

64 {
65 return K;//K * exp(tr) * (1. + tr); // K*2*ln(J)/J
66 }

◆ fun_neohookean()

static double EshelbianPlasticity::HMHNeohookean::fun_neohookean ( double  c10,
double  v 
)
inlinestatic

Definition of Neo-hookean function.

P(i, I) = 2. * c10_default * (tF(i, I) - tInvF(i, I)) + K_default * log(detF) * tInvF(i, I) tCof(i, I) = detF * tInvF(I, i) Psi = c10* (trace(tCof) + 2ln(det(F)) - 3) + (K/2.) * (log(detF))^2

Parameters
c10
v
Returns
double

Definition at line 27 of file HMHNeohookean.cpp.

27 {
28 return (2.0 * c10) * (-1. + EshelbianCore::f(2 * v));
29 }
static boost::function< double(const double)> f

◆ fun_neohookean_bulk()

static double EshelbianPlasticity::HMHNeohookean::fun_neohookean_bulk ( double  K,
double  tr 
)
inlinestatic

Definition of axiator of Neo-hookean function.

Psi = (K/2.) * (log(J))^2 = (K/2.) * (log(exp tr H))^2 = (K/2.) * (tr H)^2

Parameters
c10
v
Returns
double

Definition at line 53 of file HMHNeohookean.cpp.

53 {
54 return K * tr;
55 }

◆ getMaterialParameters()

auto EshelbianPlasticity::HMHNeohookean::getMaterialParameters ( EntityHandle  ent)
inline

Definition at line 104 of file HMHNeohookean.cpp.

104 {
105 for (auto &b : blockData) {
106 if (b.blockEnts.find(ent) != b.blockEnts.end()) {
107 return std::make_pair(b.c10, b.K);
108 }
109 }
110 if (blockData.size() != 0)
112 "Block not found for entity handle. If you mat set "
113 "block, set it to all elements");
114 return std::make_pair(c10_default, K_default);
115 }

◆ getOptions()

MoFEMErrorCode EshelbianPlasticity::HMHNeohookean::getOptions ( )
inline

Definition at line 162 of file HMHNeohookean.cpp.

162 {
164 PetscOptionsBegin(PETSC_COMM_WORLD, "neo_hookean_", "", "none");
165
166 c10_default = 1;
167 CHKERR PetscOptionsScalar("-c10", "C10", "", c10_default, &c10_default,
168 PETSC_NULLPTR);
169 K_default = 1;
170 CHKERR PetscOptionsScalar("-K", "Bulk modulus K", "", K_default, &K_default,
171 PETSC_NULLPTR);
172
173 alphaGradU = 0;
174 CHKERR PetscOptionsScalar("-viscosity_alpha_grad_u", "viscosity", "",
175 alphaGradU, &alphaGradU, PETSC_NULLPTR);
176
177 CHKERR PetscOptionsBool("-adolc_on", "adolc ON", "", adolCOn, &adolCOn,
178 PETSC_NULLPTR);
179
180 PetscOptionsEnd();
181
182 MOFEM_LOG_CHANNEL("WORLD");
183 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "MatBlock Neo-Hookean (default)")
184 << "c10 = " << c10_default << " K = " << K_default
185 << " grad alpha u = " << alphaGradU;
186 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "MatBlock Neo-Hookean (default)")
187 << "ADOL-C On " << (adolCOn)
188 ? "Yes"
189 : "No";
191 }
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.

◆ recordTape()

MoFEMErrorCode EshelbianPlasticity::HMHNeohookean::recordTape ( const int  tape,
DTensor2Ptr t_h_ptr 
)
inlinevirtual

Implements EshelbianPlasticity::PhysicalEquations.

Definition at line 237 of file HMHNeohookean.cpp.

237 {
239
240 FTENSOR_INDEX(3, i);
241 FTENSOR_INDEX(3, j);
242 FTENSOR_INDEX(3, I);
243 FTENSOR_INDEX(3, J);
244
245 auto large = [&]() {
247
248 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
249
250 auto ih = get_h();
251 if (t_h_ptr)
252 ih(i, j) = (*t_h_ptr)(i, j);
253 else {
254 ih(i, j) = t_kd(i, j);
255 }
256
257 auto r_P = get_P();
258
259 enableMinMaxUsingAbs();
260 trace_on(tape);
261
262 // Set active variables to ADOL-C
263 th(i, j) <<= get_h()(i, j);
264
265 tF(i, I) = th(i, I);
268
269 // tCof(i, I) = detF * tInvF(I, i);
270
271 // Psi = c10* (trace(tCof) + 2ln(det(F)) - 3) + (K/2.) * (log(detF))^2;
272
273 tP(i, I) = 2. * c10_default * (tF(i, I) - tInvF(i, I)) +
274 K_default * log(detF) * tInvF(i, I);
275
276 // Set dependent variables to ADOL-C
277 tP(i, j) >>= r_P(i, j);
278
279 trace_off();
280
282 };
283
286 case LARGE_ROT:
287 case MODERATE_ROT:
288 CHKERR large();
289 break;
290 default:
291 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
292 "gradApproximator not handled");
293 };
294
296 }
#define FTENSOR_INDEX(DIM, I)
Kronecker Delta class.
constexpr auto t_kd
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
MoFEMErrorCode determinantTensor3by3(T1 &t, T2 &det)
Calculate determinant 3 by 3.
constexpr IntegrationType I
static enum RotSelector gradApproximator

◆ returnOpJacobian()

UserDataOperator * EshelbianPlasticity::HMHNeohookean::returnOpJacobian ( const int  tag,
const bool  eval_rhs,
const bool  eval_lhs,
boost::shared_ptr< DataAtIntegrationPts data_ptr,
boost::shared_ptr< PhysicalEquations physics_ptr 
)
inlinevirtual

Reimplemented from EshelbianPlasticity::PhysicalEquations.

Definition at line 151 of file HMHNeohookean.cpp.

153 {
154 if (adolCOn) {
155 return PhysicalEquations::returnOpJacobian(tag, eval_rhs, eval_lhs,
156 data_ptr, physics_ptr);
157 } else {
158 return (new OpJacobian(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
159 }
160 }
virtual UserDataOperator * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)

◆ returnOpSetScale()

VolUserDataOperator * EshelbianPlasticity::HMHNeohookean::returnOpSetScale ( boost::shared_ptr< double scale_ptr,
boost::shared_ptr< PhysicalEquations physics_ptr 
)
inlinevirtual

Reimplemented from EshelbianPlasticity::PhysicalEquations.

Definition at line 138 of file HMHNeohookean.cpp.

139 {
140
141 return new OpGetScale(scale_ptr, physics_ptr);
142 };

◆ returnOpSpatialPhysical()

virtual VolUserDataOperator * EshelbianPlasticity::HMHNeohookean::returnOpSpatialPhysical ( const std::string &  field_name,
boost::shared_ptr< DataAtIntegrationPts data_ptr,
const double  alpha_u 
)
inlinevirtual

Reimplemented from EshelbianPlasticity::PhysicalEquations.

Definition at line 311 of file HMHNeohookean.cpp.

313 {
314 if (adolCOn) {
316 alpha_u);
317 } else {
318 return new OpSpatialPhysical(field_name, data_ptr, alpha_u);
319 }
320 }
constexpr auto field_name
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)

◆ returnOpSpatialPhysical_du_du()

VolUserDataOperator * EshelbianPlasticity::HMHNeohookean::returnOpSpatialPhysical_du_du ( std::string  row_field,
std::string  col_field,
boost::shared_ptr< DataAtIntegrationPts data_ptr,
const double  alpha 
)
inlinevirtual

Reimplemented from EshelbianPlasticity::PhysicalEquations.

Definition at line 356 of file HMHNeohookean.cpp.

358 {
359 if (adolCOn) {
361 row_field, col_field, data_ptr, alpha);
362 } else {
363 return new OpSpatialPhysical_du_du(row_field, col_field, data_ptr, alpha);
364 }
365 }
virtual VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)

◆ returnOpSpatialPhysicalExternalStrain()

virtual VolUserDataOperator * EshelbianPlasticity::HMHNeohookean::returnOpSpatialPhysicalExternalStrain ( const std::string &  field_name,
boost::shared_ptr< DataAtIntegrationPts data_ptr,
boost::shared_ptr< ExternalStrainVec external_strain_vec_ptr,
std::map< std::string, boost::shared_ptr< ScalingMethod > >  smv 
)
inlinevirtual

Reimplemented from EshelbianPlasticity::PhysicalEquations.

Definition at line 337 of file HMHNeohookean.cpp.

341 {
342 return new OpSpatialPhysicalExternalStrain(field_name, data_ptr,
343 external_strain_vec_ptr, smv);
344 }

Member Data Documentation

◆ adolCOn

PetscBool EshelbianPlasticity::HMHNeohookean::adolCOn = PETSC_FALSE
private

Definition at line 383 of file HMHNeohookean.cpp.

◆ alphaGradU

double EshelbianPlasticity::HMHNeohookean::alphaGradU
private

Definition at line 372 of file HMHNeohookean.cpp.

◆ blockData

std::vector<BlockData> EshelbianPlasticity::HMHNeohookean::blockData
private

Definition at line 378 of file HMHNeohookean.cpp.

◆ c10_default

double EshelbianPlasticity::HMHNeohookean::c10_default
private

Definition at line 370 of file HMHNeohookean.cpp.

◆ detF

adouble EshelbianPlasticity::HMHNeohookean::detF
private

Definition at line 387 of file HMHNeohookean.cpp.

◆ eqScaling

double EshelbianPlasticity::HMHNeohookean::eqScaling = 1.
private

Definition at line 380 of file HMHNeohookean.cpp.

◆ K_default

double EshelbianPlasticity::HMHNeohookean::K_default
private

Definition at line 371 of file HMHNeohookean.cpp.

◆ mField

MoFEM::Interface& EshelbianPlasticity::HMHNeohookean::mField
private

Definition at line 368 of file HMHNeohookean.cpp.

◆ numberOfActiveVariables

constexpr int EshelbianPlasticity::HMHNeohookean::numberOfActiveVariables = 9
staticconstexpr

Definition at line 68 of file HMHNeohookean.cpp.

◆ numberOfDependentVariables

constexpr int EshelbianPlasticity::HMHNeohookean::numberOfDependentVariables = 9
staticconstexpr

Definition at line 69 of file HMHNeohookean.cpp.

◆ tCof

ATensor2 EshelbianPlasticity::HMHNeohookean::tCof
private

Definition at line 389 of file HMHNeohookean.cpp.

◆ tF

ATensor2 EshelbianPlasticity::HMHNeohookean::tF
private

Definition at line 386 of file HMHNeohookean.cpp.

◆ th

ATensor2 EshelbianPlasticity::HMHNeohookean::th
private

Definition at line 385 of file HMHNeohookean.cpp.

◆ tInvF

ATensor2 EshelbianPlasticity::HMHNeohookean::tInvF
private

Definition at line 388 of file HMHNeohookean.cpp.

◆ tP

ATensor2 EshelbianPlasticity::HMHNeohookean::tP
private

Definition at line 390 of file HMHNeohookean.cpp.


The documentation for this struct was generated from the following file: