v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Private Attributes | List of all members
DataAtIntegrationPts Struct Reference

#include <users_modules/basic_finite_elements/src/HookeElement.hpp>

Collaboration diagram for DataAtIntegrationPts:
[legend]

Public Member Functions

 DataAtIntegrationPts ()
 
 DataAtIntegrationPts (MoFEM::Interface &m_field)
 
MoFEMErrorCode getParameters ()
 
MoFEMErrorCode getBlockData (BlockData &data)
 
MoFEMErrorCode setBlocks ()
 

Public Attributes

boost::shared_ptr< MatrixDouble > smallStrainMat
 
boost::shared_ptr< MatrixDouble > hMat
 
boost::shared_ptr< MatrixDouble > FMat
 
boost::shared_ptr< MatrixDouble > HMat
 
boost::shared_ptr< VectorDouble > detHVec
 
boost::shared_ptr< MatrixDouble > invHMat
 
boost::shared_ptr< MatrixDouble > cauchyStressMat
 
boost::shared_ptr< MatrixDouble > stiffnessMat
 
boost::shared_ptr< VectorDouble > energyVec
 
boost::shared_ptr< MatrixDouble > eshelbyStressMat
 
boost::shared_ptr< MatrixDouble > eshelbyStress_dx
 
Range forcesOnlyOnEntitiesRow
 
Range forcesOnlyOnEntitiesCol
 
boost::shared_ptr< MatrixDouble > gradDispPtr
 
boost::shared_ptr< VectorDouble > pPtr
 
FTensor::Ddg< double, 3, 3 > tD
 
double pOisson
 
double yOung
 
double lAmbda
 
double mU
 
std::map< int, BlockDatasetOfBlocksData
 

Private Attributes

MoFEM::InterfacemField
 

Detailed Description

Examples
ElasticityMixedFormulation.hpp, EshelbianPlasticity.cpp, HookeElement.cpp, HookeElement.hpp, HookeInternalStressElement.hpp, and elasticity_mixed_formulation.cpp.

Definition at line 79 of file HookeElement.hpp.

Constructor & Destructor Documentation

◆ DataAtIntegrationPts() [1/2]

DataAtIntegrationPts::DataAtIntegrationPts ( )
inline
Examples
ElasticityMixedFormulation.hpp.

Definition at line 96 of file HookeElement.hpp.

96 {
97
98 smallStrainMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
99 hMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
100 FMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
101
102 HMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
103 detHVec = boost::shared_ptr<VectorDouble>(new VectorDouble());
104 invHMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
105
106 cauchyStressMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
107 stiffnessMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
108 energyVec = boost::shared_ptr<VectorDouble>(new VectorDouble());
109 eshelbyStressMat = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
110
111 eshelbyStress_dx = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
112 }
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
boost::shared_ptr< MatrixDouble > HMat
boost::shared_ptr< MatrixDouble > invHMat
boost::shared_ptr< MatrixDouble > smallStrainMat
boost::shared_ptr< MatrixDouble > cauchyStressMat
boost::shared_ptr< MatrixDouble > eshelbyStress_dx
boost::shared_ptr< MatrixDouble > FMat
boost::shared_ptr< MatrixDouble > eshelbyStressMat
boost::shared_ptr< MatrixDouble > stiffnessMat
boost::shared_ptr< MatrixDouble > hMat
boost::shared_ptr< VectorDouble > energyVec
boost::shared_ptr< VectorDouble > detHVec

◆ DataAtIntegrationPts() [2/2]

DataAtIntegrationPts::DataAtIntegrationPts ( MoFEM::Interface m_field)
inline

Definition at line 33 of file ElasticityMixedFormulation.hpp.

33 : mField(m_field) {
34
35 // Setting default values for coeffcients
36 gradDispPtr = boost::shared_ptr<MatrixDouble>(new MatrixDouble());
37 pPtr = boost::shared_ptr<VectorDouble>(new VectorDouble());
38
39 ierr = setBlocks();
40 CHKERRABORT(PETSC_COMM_WORLD, ierr);
41 }
static PetscErrorCode ierr
boost::shared_ptr< MatrixDouble > gradDispPtr
boost::shared_ptr< VectorDouble > pPtr

Member Function Documentation

◆ getBlockData()

MoFEMErrorCode DataAtIntegrationPts::getBlockData ( BlockData data)
inline
Examples
ElasticityMixedFormulation.hpp.

Definition at line 52 of file ElasticityMixedFormulation.hpp.

52 {
54
55 yOung = data.yOung;
56 pOisson = data.pOisson;
57 lAmbda = (yOung * pOisson) / ((1. + pOisson) * (1. - 2. * pOisson));
58 mU = yOung / (2. * (1. + pOisson));
59
64
65 tD(i, j, k, l) = 0.;
66
67 tD(0, 0, 0, 0) = 2 * mU;
68 tD(0, 1, 0, 1) = mU;
69 tD(0, 1, 1, 0) = mU;
70 tD(0, 2, 0, 2) = mU;
71 tD(0, 2, 2, 0) = mU;
72 tD(1, 0, 0, 1) = mU;
73 tD(1, 0, 1, 0) = mU;
74 tD(1, 1, 1, 1) = 2 * mU;
75 tD(1, 2, 1, 2) = mU;
76 tD(1, 2, 2, 1) = mU;
77 tD(2, 0, 0, 2) = mU;
78 tD(2, 0, 2, 0) = mU;
79 tD(2, 1, 1, 2) = mU;
80 tD(2, 1, 2, 1) = mU;
81 tD(2, 2, 2, 2) = 2 * mU;
82
84 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
FTensor::Ddg< double, 3, 3 > tD

◆ getParameters()

MoFEMErrorCode DataAtIntegrationPts::getParameters ( )
inline
Examples
ElasticityMixedFormulation.hpp, and elasticity_mixed_formulation.cpp.

Definition at line 43 of file ElasticityMixedFormulation.hpp.

43 {
44 MoFEMFunctionBegin; // They will be overwriten by BlockData
45 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Problem", "none");
46
47 ierr = PetscOptionsEnd();
48 CHKERRQ(ierr);
50 }
#define CHKERR
Inline error check.
Definition: definitions.h:535

◆ setBlocks()

MoFEMErrorCode DataAtIntegrationPts::setBlocks ( )
inline
Examples
ElasticityMixedFormulation.hpp.

Definition at line 86 of file ElasticityMixedFormulation.hpp.

86 {
90 Mat_Elastic mydata;
91 CHKERR it->getAttributeDataStructure(mydata);
92 int id = it->getMeshsetId();
93 EntityHandle meshset = it->getMeshset();
94 CHKERR mField.get_moab().get_entities_by_type(
95 meshset, MBTET, setOfBlocksData[id].tEts, true);
96 setOfBlocksData[id].iD = id;
97 setOfBlocksData[id].yOung = mydata.data.Young;
98 setOfBlocksData[id].pOisson = mydata.data.Poisson;
99 }
101 }
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
Definition: definitions.h:159
@ BLOCKSET
Definition: definitions.h:148
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
std::map< int, BlockData > setOfBlocksData
virtual moab::Interface & get_moab()=0

Member Data Documentation

◆ cauchyStressMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::cauchyStressMat

Definition at line 89 of file HookeElement.hpp.

◆ detHVec

boost::shared_ptr<VectorDouble> DataAtIntegrationPts::detHVec

Definition at line 86 of file HookeElement.hpp.

◆ energyVec

boost::shared_ptr<VectorDouble> DataAtIntegrationPts::energyVec

Definition at line 91 of file HookeElement.hpp.

◆ eshelbyStress_dx

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::eshelbyStress_dx

Definition at line 94 of file HookeElement.hpp.

◆ eshelbyStressMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::eshelbyStressMat

Definition at line 92 of file HookeElement.hpp.

◆ FMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::FMat

Definition at line 83 of file HookeElement.hpp.

◆ forcesOnlyOnEntitiesCol

Range DataAtIntegrationPts::forcesOnlyOnEntitiesCol

Definition at line 115 of file HookeElement.hpp.

◆ forcesOnlyOnEntitiesRow

Range DataAtIntegrationPts::forcesOnlyOnEntitiesRow

Definition at line 114 of file HookeElement.hpp.

◆ gradDispPtr

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::gradDispPtr

◆ hMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::hMat

Definition at line 82 of file HookeElement.hpp.

◆ HMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::HMat

Definition at line 85 of file HookeElement.hpp.

◆ invHMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::invHMat

Definition at line 87 of file HookeElement.hpp.

◆ lAmbda

double DataAtIntegrationPts::lAmbda

◆ mField

MoFEM::Interface& DataAtIntegrationPts::mField
private

◆ mU

double DataAtIntegrationPts::mU

◆ pOisson

double DataAtIntegrationPts::pOisson

◆ pPtr

boost::shared_ptr<VectorDouble> DataAtIntegrationPts::pPtr

◆ setOfBlocksData

std::map<int, BlockData> DataAtIntegrationPts::setOfBlocksData

◆ smallStrainMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::smallStrainMat

Definition at line 81 of file HookeElement.hpp.

◆ stiffnessMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::stiffnessMat

Definition at line 90 of file HookeElement.hpp.

◆ tD

FTensor::Ddg<double, 3, 3> DataAtIntegrationPts::tD

◆ yOung

double DataAtIntegrationPts::yOung

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