v0.14.0
Classes | Public Member Functions | Private Attributes | List of all members
EshelbianPlasticity::HMHHencky Struct Reference
Inheritance diagram for EshelbianPlasticity::HMHHencky:
[legend]
Collaboration diagram for EshelbianPlasticity::HMHHencky:
[legend]

Classes

struct  BlockData
 
struct  OpHenckyJacobian
 

Public Member Functions

 HMHHencky (MoFEM::Interface &m_field, const double E, const double nu)
 
MoFEMErrorCode recordTape (const int tag, DTensor2Ptr *t_h)
 
virtual OpJacobianreturnOpJacobian (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 (boost::shared_ptr< DataAtIntegrationPts > data_ptr)
 
MoFEMErrorCode extractBlockData (Sev sev)
 
MoFEMErrorCode extractBlockData (std::vector< const CubitMeshSets * > meshset_vec_ptr, Sev sev)
 
MoFEMErrorCode getMatDPtr (boost::shared_ptr< MatrixDouble > mat_D_ptr, boost::shared_ptr< MatrixDouble > mat_axiator_D_ptr, boost::shared_ptr< MatrixDouble > mat_deviator_D_ptr, double bulk_modulus_K, double shear_modulus_G)
 
- Public Member Functions inherited from EshelbianPlasticity::PhysicalEquations
 PhysicalEquations ()=delete
 
 PhysicalEquations (const int size_active, const int size_dependent)
 
virtual ~PhysicalEquations ()=default
 
DTensor2Ptr get_P ()
 
DTensor3Ptr get_P_dh0 ()
 
DTensor3Ptr get_P_dh1 ()
 
DTensor3Ptr get_P_dh2 ()
 
DTensor2Ptr get_h ()
 

Private Attributes

MoFEM::InterfacemField
 
std::vector< BlockDatablockData
 
double E
 
double nu
 

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
 
- 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)
 
- Public Attributes inherited from EshelbianPlasticity::PhysicalEquations
std::vector< doubleactiveVariables
 
std::vector< doubledependentVariablesPiola
 
std::vector< doubledependentVariablesPiolaDirevatives
 

Detailed Description

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

Constructor & Destructor Documentation

◆ HMHHencky()

EshelbianPlasticity::HMHHencky::HMHHencky ( MoFEM::Interface m_field,
const double  E,
const double  nu 
)
inline

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

468  : PhysicalEquations(0, 0), mField(m_field), E(E), nu(nu) {}

Member Function Documentation

◆ extractBlockData() [1/2]

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

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

547  {
548  return extractBlockData(
549 
551 
552  (boost::format("%s(.*)") % "MAT_ELASTIC").str()
553 
554  )),
555 
556  sev);
557  }

◆ extractBlockData() [2/2]

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

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

561  {
563 
564  for (auto m : meshset_vec_ptr) {
565  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock") << *m;
566  std::vector<double> block_data;
567  CHKERR m->getAttributes(block_data);
568  if (block_data.size() != 2) {
569  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
570  "Expected that block has two attribute");
571  }
572  auto get_block_ents = [&]() {
573  Range ents;
574  CHKERR mField.get_moab().get_entities_by_handle(m->meshset, ents, true);
575  return ents;
576  };
577 
578  double young_modulus = block_data[0];
579  double poisson_ratio = block_data[1];
580  double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
581  double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio));
582 
583  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
584  << "E = " << young_modulus << " nu = " << poisson_ratio;
585 
586  blockData.push_back({bulk_modulus_K, shear_modulus_G, get_block_ents()});
587  }
588  MOFEM_LOG_CHANNEL("WORLD");
590  }

◆ getMatDPtr()

MoFEMErrorCode EshelbianPlasticity::HMHHencky::getMatDPtr ( boost::shared_ptr< MatrixDouble >  mat_D_ptr,
boost::shared_ptr< MatrixDouble >  mat_axiator_D_ptr,
boost::shared_ptr< MatrixDouble >  mat_deviator_D_ptr,
double  bulk_modulus_K,
double  shear_modulus_G 
)
inline

[Calculate elasticity tensor]

[Calculate elasticity tensor]

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

595  {
597  //! [Calculate elasticity tensor]
598  auto set_material_stiffness = [&]() {
604  double A = (SPACE_DIM == 2)
605  ? 2 * shear_modulus_G /
606  (bulk_modulus_K + (4. / 3.) * shear_modulus_G)
607  : 1;
608  auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mat_D_ptr);
609  auto t_axiator_D =
610  getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mat_axiator_D_ptr);
611  auto t_deviator_D =
612  getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*mat_deviator_D_ptr);
613  t_axiator_D(i, j, k, l) = A *
614  (bulk_modulus_K - (2. / 3.) * shear_modulus_G) *
615  t_kd(i, j) * t_kd(k, l);
616  t_deviator_D(i, j, k, l) =
617  2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.);
618  t_D(i, j, k, l) = t_axiator_D(i, j, k, l) + t_deviator_D(i, j, k, l);
619  };
620  //! [Calculate elasticity tensor]
621  constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
622  mat_D_ptr->resize(size_symm * size_symm, 1);
623  mat_axiator_D_ptr->resize(size_symm * size_symm, 1);
624  mat_deviator_D_ptr->resize(size_symm * size_symm, 1);
625  set_material_stiffness();
627  }

◆ getOptions()

MoFEMErrorCode EshelbianPlasticity::HMHHencky::getOptions ( boost::shared_ptr< DataAtIntegrationPts data_ptr)
inline

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

532  {
534  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "hencky_", "", "none");
535 
536  CHKERR PetscOptionsScalar("-young_modulus", "Young modulus", "", E, &E,
537  PETSC_NULL);
538  CHKERR PetscOptionsScalar("-poisson_ratio", "poisson ratio", "", nu, &nu,
539  PETSC_NULL);
540 
541  ierr = PetscOptionsEnd();
542  CHKERRG(ierr);
543 
545  }

◆ recordTape()

MoFEMErrorCode EshelbianPlasticity::HMHHencky::recordTape ( const int  tag,
DTensor2Ptr t_h 
)
inlinevirtual

Implements EshelbianPlasticity::PhysicalEquations.

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

470 { return 0; }

◆ returnOpJacobian()

virtual OpJacobian* EshelbianPlasticity::HMHHencky::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

Implements EshelbianPlasticity::PhysicalEquations.

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

527  {
528  return (new OpHenckyJacobian(
529  data_ptr, boost::dynamic_pointer_cast<HMHHencky>(physics_ptr)));
530  }

Member Data Documentation

◆ blockData

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

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

◆ E

double EshelbianPlasticity::HMHHencky::E
private

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

◆ mField

MoFEM::Interface& EshelbianPlasticity::HMHHencky::mField
private

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

◆ nu

double EshelbianPlasticity::HMHHencky::nu
private

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


The documentation for this struct was generated from the following file:
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
EshelbianPlasticity::size_symm
constexpr static auto size_symm
Definition: EshelbianOperators.cpp:47
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:172
E
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
bulk_modulus_K
double bulk_modulus_K
Definition: dynamic_first_order_con_law.cpp:96
EshelbianPlasticity::HMHHencky::blockData
std::vector< BlockData > blockData
Definition: EshelbianADOL-C.cpp:637
EshelbianPlasticity::HMHHencky::nu
double nu
Definition: EshelbianADOL-C.cpp:640
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:173
EshelbianPlasticity::HMHHencky::E
double E
Definition: EshelbianADOL-C.cpp:639
EshelbianPlasticity::HMHHencky::mField
MoFEM::Interface & mField
Definition: EshelbianADOL-C.cpp:630
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
EshelbianPlasticity::HMHHencky::extractBlockData
MoFEMErrorCode extractBlockData(Sev sev)
Definition: EshelbianADOL-C.cpp:547
FTensor::Index< 'i', SPACE_DIM >
Range
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
shear_modulus_G
double shear_modulus_G
Definition: dynamic_first_order_con_law.cpp:97
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
MoFEM::MeshsetsManager::getCubitMeshsetPtr
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
Definition: MeshsetsManager.cpp:575
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
EshelbianPlasticity::PhysicalEquations::PhysicalEquations
PhysicalEquations()=delete
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21