v0.15.0
Loading...
Searching...
No Matches
DataAtIntegrationPts Struct Reference

#include "tutorials/cor-7/src/ElasticityMixedFormulation.hpp"

Collaboration diagram for DataAtIntegrationPts:
[legend]

Public Member Functions

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

Public Attributes

boost::shared_ptr< MatrixDoublegradDispPtr
 
boost::shared_ptr< VectorDoublepPtr
 
FTensor::Ddg< double, 3, 3 > tD
 
double pOisson
 
double yOung
 
double lAmbda
 
double mU
 
std::map< int, BlockDatasetOfBlocksData
 
SmartPetscObj< Vec > petscVec
 
double blockPermittivity
 
double blockChrgDens
 
boost::shared_ptr< MatrixDoublesmallStrainMat
 
boost::shared_ptr< MatrixDoublehMat
 
boost::shared_ptr< MatrixDoubleFMat
 
boost::shared_ptr< MatrixDoubleHMat
 
boost::shared_ptr< VectorDoubledetHVec
 
boost::shared_ptr< MatrixDoubleinvHMat
 
boost::shared_ptr< MatrixDoublecauchyStressMat
 
boost::shared_ptr< MatrixDoublestiffnessMat
 
boost::shared_ptr< VectorDoubleenergyVec
 
boost::shared_ptr< MatrixDoubleeshelbyStressMat
 
boost::shared_ptr< MatrixDoubleeshelbyStress_dx
 
Range forcesOnlyOnEntitiesRow
 
Range forcesOnlyOnEntitiesCol
 

Private Attributes

MoFEM::InterfacemField
 

Detailed Description

Constructor & Destructor Documentation

◆ DataAtIntegrationPts() [1/3]

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 MoFEMErrorCodeGeneric< PetscErrorCode > ierr
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
UBlasVector< double > VectorDouble
Definition Types.hpp:68
boost::shared_ptr< MatrixDouble > gradDispPtr
boost::shared_ptr< VectorDouble > pPtr

◆ DataAtIntegrationPts() [2/3]

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

Definition at line 62 of file electrostatics.hpp.

62 {
64 blockChrgDens = 0.0;
65 PetscInt ghosts[2] = {0, 1};
66 if (!m_field.get_comm_rank())
67 petscVec = createGhostVector(m_field.get_comm(), 2, 2, 0, ghosts);
68 else
69 petscVec = createGhostVector(m_field.get_comm(), 0, 2, 2, ghosts);
70 }
auto createGhostVector(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[])
Create smart ghost vector.
SmartPetscObj< Vec > petscVec
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0

◆ DataAtIntegrationPts() [3/3]

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 }
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

Member Function Documentation

◆ getBlockData()

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

Definition at line 51 of file ElasticityMixedFormulation.hpp.

51 {
53
54 yOung = data.yOung;
55 pOisson = data.pOisson;
56 lAmbda = (yOung * pOisson) / ((1. + pOisson) * (1. - 2. * pOisson));
57 mU = yOung / (2. * (1. + pOisson));
58
59 FTensor::Index<'i', 3> i;
60 FTensor::Index<'j', 3> j;
61 FTensor::Index<'k', 3> k;
62 FTensor::Index<'l', 3> l;
63
64 tD(i, j, k, l) = 0.;
65
66 tD(0, 0, 0, 0) = 2 * mU;
67 tD(0, 1, 0, 1) = mU;
68 tD(0, 1, 1, 0) = mU;
69 tD(0, 2, 0, 2) = mU;
70 tD(0, 2, 2, 0) = mU;
71 tD(1, 0, 0, 1) = mU;
72 tD(1, 0, 1, 0) = mU;
73 tD(1, 1, 1, 1) = 2 * mU;
74 tD(1, 2, 1, 2) = mU;
75 tD(1, 2, 2, 1) = mU;
76 tD(2, 0, 0, 2) = mU;
77 tD(2, 0, 2, 0) = mU;
78 tD(2, 1, 1, 2) = mU;
79 tD(2, 1, 2, 1) = mU;
80 tD(2, 2, 2, 2) = 2 * mU;
81
83 }
#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()
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 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Problem", "none");
46
47 PetscOptionsEnd();
49 }

◆ setBlocks()

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

Definition at line 85 of file ElasticityMixedFormulation.hpp.

85 {
89 Mat_Elastic mydata;
90 CHKERR it->getAttributeDataStructure(mydata);
91 int id = it->getMeshsetId();
92 EntityHandle meshset = it->getMeshset();
93 CHKERR mField.get_moab().get_entities_by_type(
94 meshset, MBTET, setOfBlocksData[id].tEts, true);
95 setOfBlocksData[id].iD = id;
96 setOfBlocksData[id].yOung = mydata.data.Young;
97 setOfBlocksData[id].pOisson = mydata.data.Poisson;
98 }
100 }
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
@ BLOCKSET
#define CHKERR
Inline error check.
#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

◆ blockChrgDens

double DataAtIntegrationPts::blockChrgDens

Definition at line 61 of file electrostatics.hpp.

◆ blockPermittivity

double DataAtIntegrationPts::blockPermittivity

Definition at line 60 of file electrostatics.hpp.

◆ 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 85 of file HookeElement.hpp.

◆ hMat

boost::shared_ptr<MatrixDouble> DataAtIntegrationPts::hMat

Definition at line 82 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

◆ petscVec

SmartPetscObj<Vec> DataAtIntegrationPts::petscVec

Definition at line 59 of file electrostatics.hpp.

◆ 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: