v0.15.0
Loading...
Searching...
No Matches
MetaSpringBC::DataAtIntegrationPtsSprings Struct Reference

#include "users_modules/basic_finite_elements/src/SpringElement.hpp"

Inheritance diagram for MetaSpringBC::DataAtIntegrationPtsSprings:
[legend]
Collaboration diagram for MetaSpringBC::DataAtIntegrationPtsSprings:
[legend]

Public Member Functions

 DataAtIntegrationPtsSprings (MoFEM::Interface &m_field, double scale_stiffness=1.)
 
MoFEMErrorCode getParameters ()
 
MoFEMErrorCode getBlockData (BlockOptionDataSprings &data)
 
MoFEMErrorCode setBlocks ()
 

Public Attributes

boost::shared_ptr< MatrixDouble > gradDispPtr
 
boost::shared_ptr< MatrixDouble > xAtPts = boost::make_shared<MatrixDouble>()
 
boost::shared_ptr< MatrixDouble > xInitAtPts
 
boost::shared_ptr< MatrixDouble > hMat = boost::make_shared<MatrixDouble>()
 
boost::shared_ptr< MatrixDouble > FMat = boost::make_shared<MatrixDouble>()
 
boost::shared_ptr< MatrixDouble > HMat = boost::make_shared<MatrixDouble>()
 
boost::shared_ptr< MatrixDouble > invHMat
 
boost::shared_ptr< VectorDouble > detHVec
 
MatrixDouble tangent1
 
MatrixDouble tangent2
 
MatrixDouble normalVector
 
double springStiffnessNormal
 
double springStiffnessTangent
 
double scaleStiffness
 
Range forcesOnlyOnEntitiesRow
 
EntitiesFieldData::EntData * faceRowData
 
std::map< int, BlockOptionDataSpringsmapSpring
 

Private Attributes

MoFEM::InterfacemField
 

Detailed Description

Definition at line 32 of file SpringElement.hpp.

Constructor & Destructor Documentation

◆ DataAtIntegrationPtsSprings()

MetaSpringBC::DataAtIntegrationPtsSprings::DataAtIntegrationPtsSprings ( MoFEM::Interface & m_field,
double scale_stiffness = 1. )
inline

Definition at line 63 of file SpringElement.hpp.

64 : mField(m_field), faceRowData(nullptr), scaleStiffness(scale_stiffness) {
65
66 ierr = setBlocks();
67 CHKERRABORT(PETSC_COMM_WORLD, ierr);
68 }
static PetscErrorCode ierr
EntitiesFieldData::EntData * faceRowData

Member Function Documentation

◆ getBlockData()

MoFEMErrorCode MetaSpringBC::DataAtIntegrationPtsSprings::getBlockData ( BlockOptionDataSprings & data)
inline

Definition at line 78 of file SpringElement.hpp.

78 {
80
81 springStiffnessNormal = data.springStiffnessNormal * scaleStiffness;
82 springStiffnessTangent = data.springStiffnessTangent * scaleStiffness;
83
85 }
#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()

◆ getParameters()

MoFEMErrorCode MetaSpringBC::DataAtIntegrationPtsSprings::getParameters ( )
inline

Definition at line 70 of file SpringElement.hpp.

70 {
71 MoFEMFunctionBegin; // They will be overwritten by BlockData
72 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Problem", "none");
73
74 PetscOptionsEnd();
76 }

◆ setBlocks()

MoFEMErrorCode MetaSpringBC::DataAtIntegrationPtsSprings::setBlocks ( )
inline

Definition at line 87 of file SpringElement.hpp.

87 {
89
91 if (bit->getName().compare(0, 9, "SPRING_BC") == 0) {
92
93 const int id = bit->getMeshsetId();
94 mapSpring[id].tRis.clear();
95 CHKERR mField.get_moab().get_entities_by_dimension(
96 bit->getMeshset(), 2, mapSpring[id].tRis, true);
97
98 std::vector<double> attributes;
99 bit->getAttributes(attributes);
100 if (attributes.size() < 2) {
101 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
102 "Springs should have 2 attributes but there is %ld",
103 attributes.size());
104 }
105 mapSpring[id].iD = id;
106 mapSpring[id].springStiffnessNormal = attributes[0];
107 mapSpring[id].springStiffnessTangent = attributes[1];
108
109 // Print spring blocks after being read
110 MOFEM_LOG_C("WORLD", Sev::verbose, "\nSpring block %d\n", id);
111 MOFEM_LOG_C("WORLD", Sev::verbose, "\tNormal stiffness %3.4g\n",
112 attributes[0]);
113 MOFEM_LOG_C("WORLD", Sev::verbose, "\tTangent stiffness %3.4g\n",
114 attributes[1]);
115 }
116 }
117
119 }
#define MOFEM_LOG_C(channel, severity, format,...)
@ BLOCKSET
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define CHKERR
Inline error check.
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
auto bit
set bit
std::map< int, BlockOptionDataSprings > mapSpring
virtual moab::Interface & get_moab()=0

Member Data Documentation

◆ detHVec

boost::shared_ptr<VectorDouble> MetaSpringBC::DataAtIntegrationPtsSprings::detHVec
Initial value:
=
boost::make_shared<VectorDouble>()

Definition at line 46 of file SpringElement.hpp.

◆ faceRowData

EntitiesFieldData::EntData* MetaSpringBC::DataAtIntegrationPtsSprings::faceRowData

Definition at line 59 of file SpringElement.hpp.

◆ FMat

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::FMat = boost::make_shared<MatrixDouble>()

Definition at line 42 of file SpringElement.hpp.

◆ forcesOnlyOnEntitiesRow

Range MetaSpringBC::DataAtIntegrationPtsSprings::forcesOnlyOnEntitiesRow

Definition at line 57 of file SpringElement.hpp.

◆ gradDispPtr

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::gradDispPtr
Initial value:
=
boost::make_shared<MatrixDouble>()

Definition at line 35 of file SpringElement.hpp.

◆ HMat

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::HMat = boost::make_shared<MatrixDouble>()

Definition at line 43 of file SpringElement.hpp.

◆ hMat

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::hMat = boost::make_shared<MatrixDouble>()

Definition at line 41 of file SpringElement.hpp.

◆ invHMat

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::invHMat
Initial value:
=
boost::make_shared<MatrixDouble>()

Definition at line 44 of file SpringElement.hpp.

◆ mapSpring

std::map<int, BlockOptionDataSprings> MetaSpringBC::DataAtIntegrationPtsSprings::mapSpring

Definition at line 61 of file SpringElement.hpp.

◆ mField

MoFEM::Interface& MetaSpringBC::DataAtIntegrationPtsSprings::mField
private

Definition at line 122 of file SpringElement.hpp.

◆ normalVector

MatrixDouble MetaSpringBC::DataAtIntegrationPtsSprings::normalVector

Definition at line 51 of file SpringElement.hpp.

◆ scaleStiffness

double MetaSpringBC::DataAtIntegrationPtsSprings::scaleStiffness

Definition at line 55 of file SpringElement.hpp.

◆ springStiffnessNormal

double MetaSpringBC::DataAtIntegrationPtsSprings::springStiffnessNormal

Definition at line 53 of file SpringElement.hpp.

◆ springStiffnessTangent

double MetaSpringBC::DataAtIntegrationPtsSprings::springStiffnessTangent

Definition at line 54 of file SpringElement.hpp.

◆ tangent1

MatrixDouble MetaSpringBC::DataAtIntegrationPtsSprings::tangent1

Definition at line 49 of file SpringElement.hpp.

◆ tangent2

MatrixDouble MetaSpringBC::DataAtIntegrationPtsSprings::tangent2

Definition at line 50 of file SpringElement.hpp.

◆ xAtPts

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::xAtPts = boost::make_shared<MatrixDouble>()

Definition at line 37 of file SpringElement.hpp.

◆ xInitAtPts

boost::shared_ptr<MatrixDouble> MetaSpringBC::DataAtIntegrationPtsSprings::xInitAtPts
Initial value:
=
boost::make_shared<MatrixDouble>()

Definition at line 38 of file SpringElement.hpp.


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