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

Public Member Functions

 HMHStVenantKirchhoff (const double lambda, const double mu, const double sigma_y)
 
MoFEMErrorCode getOptions ()
 
virtual OpJacobianreturnOpJacobian (const std::string &field_name, const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< PhysicalEquations > &physics_ptr)
 
MoFEMErrorCode recordTape (const int tape, DTensor2Ptr *t_h_ptr)
 
- Public Member Functions inherited from EshelbianPlasticity::PhysicalEquations
 PhysicalEquations ()=delete
 
 PhysicalEquations (const int size_active, const int size_dependent)
 
virtual ~PhysicalEquations ()
 
virtual MoFEMErrorCode recordTape (const int tag, DTensor2Ptr *t_h)=0
 
virtual OpJacobianreturnOpJacobian (const std::string &field_name, const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< PhysicalEquations > &physics_ptr)=0
 
DTensor2Ptr get_P ()
 
DTensor2Ptr get_Sigma ()
 
DTensor0Ptr get_Phi ()
 
doubleget_PhiRef ()
 
DTensor2Ptr get_Flow ()
 
DTensor3Ptr get_P_dh0 ()
 
DTensor3Ptr get_P_dH0 ()
 
DTensor3Ptr get_P_dh1 ()
 
DTensor3Ptr get_P_dH1 ()
 
DTensor3Ptr get_P_dh2 ()
 
DTensor3Ptr get_P_dH2 ()
 
DTensor3Ptr get_Sigma_dh0 ()
 
DTensor3Ptr get_Sigma_dH0 ()
 
DTensor3Ptr get_Sigma_dh1 ()
 
DTensor3Ptr get_Sigma_dH1 ()
 
DTensor3Ptr get_Sigma_dh2 ()
 
DTensor3Ptr get_Sigma_dH2 ()
 
DTensor2Ptr get_Phi_dh ()
 
DTensor2Ptr get_Phi_dH ()
 
DTensor3Ptr get_Flow_dh0 ()
 
DTensor3Ptr get_Flow_dH0 ()
 
DTensor3Ptr get_Flow_dh1 ()
 
DTensor3Ptr get_Flow_dH1 ()
 
DTensor3Ptr get_Flow_dh2 ()
 
DTensor3Ptr get_Flow_dH2 ()
 
DTensor2Ptr get_h ()
 
DTensor2Ptr get_H ()
 

Public Attributes

double lambda
 
double mu
 
double sigmaY
 
ATensor2 th
 
ATensor2 tH
 
ATensor2 tP
 
ATensor2 tSigma
 
adouble phi
 
adouble s2
 
adouble f
 
adouble energy
 
adouble detH
 
adouble detF
 
adouble trE
 
adouble meanCauchy
 
ATensor2 tInvH
 
ATensor2 tF
 
ATensor2 tInvF
 
ATensor2 tC
 
ATensor2 tE
 
ATensor2 tS
 
ATensor2 tCauchy
 
ATensor2 tDevCauchy
 
ATensor2 tPhi_dDevCauchy
 
ATensor2 tPhi_dCauchy
 
ATensor2 tPhi_dSigma
 
ATensor2 tPulledP
 
ATensor2 tPulledSigma
 
ATensor2 tPulledPhi_dSigma
 
- Public Attributes inherited from EshelbianPlasticity::PhysicalEquations
std::vector< doubleactiveVariables
 
std::vector< doubledependentVariables
 
std::vector< doubledependentVariablesDirevatives
 

Static Public Attributes

static constexpr int numberOfActiveVariables = 9 + 9
 
static constexpr int numberOfDependentVariables = 9 + 9 + 1 + 9
 

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::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
 
- 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, int S1, int S2>
static DTensor3Ptr get_vecTensor3 (std::vector< double > &v)
 

Detailed Description

Definition at line 231 of file EshelbianADOL-C.cpp.

Constructor & Destructor Documentation

◆ HMHStVenantKirchhoff()

EshelbianPlasticity::HMHStVenantKirchhoff::HMHStVenantKirchhoff ( const double  lambda,
const double  mu,
const double  sigma_y 
)
inline

Member Function Documentation

◆ getOptions()

MoFEMErrorCode EshelbianPlasticity::HMHStVenantKirchhoff::getOptions ( )
inline

Definition at line 241 of file EshelbianADOL-C.cpp.

241 {
243
244 double E = 1;
245 double nu = 0;
246
247 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "stvenant_", "", "none");
248
249 CHKERR PetscOptionsScalar("-young_modulus", "Young modulus", "", E, &E,
250 PETSC_NULL);
251 CHKERR PetscOptionsScalar("-poisson_ratio", "poisson ratio", "", nu, &nu,
252 PETSC_NULL);
253 CHKERR PetscOptionsScalar("-sigmaY", "plastic strain", "", sigmaY, &sigmaY,
254 PETSC_NULL);
255
256 ierr = PetscOptionsEnd();
257 CHKERRG(ierr);
258
259 lambda = LAMBDA(E, nu);
260 mu = MU(E, nu);
261
263 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MU(E, NU)
Definition: fem_tools.h:23
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76

◆ recordTape()

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

Implements EshelbianPlasticity::PhysicalEquations.

Definition at line 308 of file EshelbianADOL-C.cpp.

308 {
317
319
320 auto ih = get_h();
321 auto iH = get_H();
322 if (t_h_ptr)
323 ih(i, j) = (*t_h_ptr)(i, j);
324 else {
325 ih(i, j) = 0;
326 for (auto ii : {0, 1, 2})
327 ih(ii, ii) = 1;
328 }
329
330 iH(i, j) = 0;
331 for (auto ii : {0, 1, 2})
332 iH(ii, ii) = 1;
333
334 auto r_P = get_P();
335 auto r_Sigma = get_Sigma();
336 double &r_phi = get_PhiRef();
337 auto r_Flow = get_Flow();
338
339 enableMinMaxUsingAbs();
340 trace_on(tape);
341
342 // Set active variables to ADOL-C
343 th(i, j) <<= get_h()(i, j);
344 tH(i, j) <<= get_H()(i, j);
345
346 // Deformation gradient
349 tF(i, I) = th(i, J) * tInvH(J, I);
352
353 // Deformation and strain
354 tC(I, J) = tF(i, I) * tF(i, J);
355 tE(I, J) = tC(I, J);
356 tE(N0, N0) -= 1;
357 tE(N1, N1) -= 1;
358 tE(N2, N2) -= 1;
359 tE(I, J) *= 0.5;
360
361 // Energy
362 trE = tE(I, I);
363 energy = 0.5 * lambda * trE * trE + mu * (tE(I, J) * tE(I, J));
364
365 // Stress Piola II
366 tS(I, J) = tE(I, J);
367 tS(I, J) *= 2 * mu;
368 tS(N0, N0) += lambda * trE;
369 tS(N1, N1) += lambda * trE;
370 tS(N2, N2) += lambda * trE;
371 // Stress Piola I
372 tP(i, J) = tF(i, I) * tS(I, J);
373 // Stress Cauchy
374 tCauchy(i, j) = tP(i, J) * tF(j, J);
375 tCauchy(i, j) *= (1 / detF);
376
377 // Stress Eshelby
378 tSigma(I, J) = -tF(i, I) * tP(i, J);
379 tSigma(N0, N0) += energy;
380 tSigma(N1, N1) += energy;
381 tSigma(N2, N2) += energy;
382
383 // Deviator Cauchy Stress
384 meanCauchy = third * tCauchy(i, i);
385 tDevCauchy(i, j) = tCauchy(i, j);
389
390 // Plastic surface
391
392 s2 = tDevCauchy(i, j) * tDevCauchy(i, j); // s:s
393 f = sqrt(1.5 * s2);
394 phi = f - sigmaY;
395
396 // Flow
398 tPhi_dDevCauchy(i, j) *= 1.5 / f;
403 tPhi_dSigma(I, J) = tInvF(I, i) * (tPhi_dCauchy(i, j) * tF(j, J));
404 tPhi_dSigma(I, J) *= -(1 / detF);
405
406 // Pull back
407 tPulledP(i, J) = tP(i, I) * tInvH(J, I);
408 tPulledP(i, J) *= detH;
409 tPulledSigma(i, J) = tSigma(i, I) * tInvH(J, I);
410 tPulledSigma(i, J) *= detH;
413
414 // Set dependent variables to ADOL-C
415 tPulledP(i, j) >>= r_P(i, j);
416 tPulledSigma(i, j) >>= r_Sigma(i, j);
417 phi >>= r_phi;
418 tPulledPhi_dSigma(i, j) >>= r_Flow(i, j);
419
420 trace_off();
421
423 }
static Index< 'J', 3 > J
static Number< 2 > N2
static Number< 1 > N1
static Number< 0 > N0
constexpr double third
FTensor::Index< 'i', SPACE_DIM > i
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.
Definition: Templates.hpp:1363
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1352
constexpr IntegrationType I

◆ returnOpJacobian()

virtual OpJacobian * EshelbianPlasticity::HMHStVenantKirchhoff::returnOpJacobian ( const std::string &  field_name,
const int  tag,
const bool  eval_rhs,
const bool  eval_lhs,
boost::shared_ptr< DataAtIntegrationPts > &  data_ptr,
boost::shared_ptr< PhysicalEquations > &  physics_ptr 
)
inlinevirtual

Implements EshelbianPlasticity::PhysicalEquations.

Definition at line 266 of file EshelbianADOL-C.cpp.

269 {
270 return (
271 new OpHMHH(field_name, tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
272 }
constexpr auto field_name

Member Data Documentation

◆ detF

adouble EshelbianPlasticity::HMHStVenantKirchhoff::detF

Definition at line 288 of file EshelbianADOL-C.cpp.

◆ detH

adouble EshelbianPlasticity::HMHStVenantKirchhoff::detH

Definition at line 287 of file EshelbianADOL-C.cpp.

◆ energy

adouble EshelbianPlasticity::HMHStVenantKirchhoff::energy

Definition at line 285 of file EshelbianADOL-C.cpp.

◆ f

adouble EshelbianPlasticity::HMHStVenantKirchhoff::f

Definition at line 284 of file EshelbianADOL-C.cpp.

◆ lambda

double EshelbianPlasticity::HMHStVenantKirchhoff::lambda

Definition at line 274 of file EshelbianADOL-C.cpp.

◆ meanCauchy

adouble EshelbianPlasticity::HMHStVenantKirchhoff::meanCauchy

Definition at line 290 of file EshelbianADOL-C.cpp.

◆ mu

double EshelbianPlasticity::HMHStVenantKirchhoff::mu

Definition at line 275 of file EshelbianADOL-C.cpp.

◆ numberOfActiveVariables

constexpr int EshelbianPlasticity::HMHStVenantKirchhoff::numberOfActiveVariables = 9 + 9
staticconstexpr

Definition at line 233 of file EshelbianADOL-C.cpp.

◆ numberOfDependentVariables

constexpr int EshelbianPlasticity::HMHStVenantKirchhoff::numberOfDependentVariables = 9 + 9 + 1 + 9
staticconstexpr

Definition at line 234 of file EshelbianADOL-C.cpp.

◆ phi

adouble EshelbianPlasticity::HMHStVenantKirchhoff::phi

Definition at line 282 of file EshelbianADOL-C.cpp.

◆ s2

adouble EshelbianPlasticity::HMHStVenantKirchhoff::s2

Definition at line 283 of file EshelbianADOL-C.cpp.

◆ sigmaY

double EshelbianPlasticity::HMHStVenantKirchhoff::sigmaY

Definition at line 276 of file EshelbianADOL-C.cpp.

◆ tC

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tC

Definition at line 294 of file EshelbianADOL-C.cpp.

◆ tCauchy

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tCauchy

Definition at line 297 of file EshelbianADOL-C.cpp.

◆ tDevCauchy

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tDevCauchy

Definition at line 298 of file EshelbianADOL-C.cpp.

◆ tE

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tE

Definition at line 295 of file EshelbianADOL-C.cpp.

◆ tF

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tF

Definition at line 292 of file EshelbianADOL-C.cpp.

◆ th

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::th

Definition at line 278 of file EshelbianADOL-C.cpp.

◆ tH

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tH

Definition at line 279 of file EshelbianADOL-C.cpp.

◆ tInvF

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tInvF

Definition at line 293 of file EshelbianADOL-C.cpp.

◆ tInvH

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tInvH

Definition at line 291 of file EshelbianADOL-C.cpp.

◆ tP

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tP

Definition at line 280 of file EshelbianADOL-C.cpp.

◆ tPhi_dCauchy

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tPhi_dCauchy

Definition at line 301 of file EshelbianADOL-C.cpp.

◆ tPhi_dDevCauchy

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tPhi_dDevCauchy

Definition at line 300 of file EshelbianADOL-C.cpp.

◆ tPhi_dSigma

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tPhi_dSigma

Definition at line 302 of file EshelbianADOL-C.cpp.

◆ tPulledP

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tPulledP

Definition at line 304 of file EshelbianADOL-C.cpp.

◆ tPulledPhi_dSigma

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tPulledPhi_dSigma

Definition at line 306 of file EshelbianADOL-C.cpp.

◆ tPulledSigma

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tPulledSigma

Definition at line 305 of file EshelbianADOL-C.cpp.

◆ trE

adouble EshelbianPlasticity::HMHStVenantKirchhoff::trE

Definition at line 289 of file EshelbianADOL-C.cpp.

◆ tS

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tS

Definition at line 296 of file EshelbianADOL-C.cpp.

◆ tSigma

ATensor2 EshelbianPlasticity::HMHStVenantKirchhoff::tSigma

Definition at line 281 of file EshelbianADOL-C.cpp.


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