v0.13.1
Public Member Functions | Public Attributes | List of all members
DirichletDisplacementBc Struct Reference

Set Dirichlet boundary conditions on displacements. More...

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

Inheritance diagram for DirichletDisplacementBc:
[legend]
Collaboration diagram for DirichletDisplacementBc:
[legend]

Public Member Functions

 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name, Mat Aij, Vec X, Vec F, string blockset_name="DISPLACEMENT")
 
 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name, string blockset_name="DISPLACEMENT")
 
virtual MoFEMErrorCode iNitialize ()
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
MoFEMErrorCode getBcDataFromSetsAndBlocks (std::vector< DataFromBc > &bc_data)
 Get the Bc Data From Sets And Blocks object Use DISPLACEMENT blockset name (default) with 6 atributes: 1,2,3 are values of displacements x,y,z 4,5,6 are flags for x,y,z (0 or 1) More...
 
MoFEMErrorCode getRotationBcFromBlock (std::vector< DataFromBc > &bc_data)
 Get the Rotation Bc From Block object Use ROTATION blockset name with 7 atributes: 1,2,3 are x,y,z coords of the center of rotation 4,5,6 are are angular velocities in x,y,z. More...
 
MoFEMErrorCode calculateRotationForDof (VectorDouble3 &coords, DataFromBc &bc_data)
 Calculate displacements from rotation for particular dof. More...
 
MoFEMErrorCode calculateRotationForDof (EntityHandle ent, DataFromBc &bc_data)
 
MoFEMErrorCode applyScaleBcData (DataFromBc &bc_data)
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Public Attributes

MoFEM::InterfacemField
 
const std::string fieldName
 field name to set Dirichlet BC More...
 
double dIag
 diagonal value set on zeroed column and rows More...
 
std::map< DofIdx, FieldData > mapZeroRows
 
std::vector< int > dofsIndices
 
std::vector< doubledofsValues
 
std::vector< doubledofsXValues
 
const std::string blocksetName
 
boost::ptr_vector< MethodForForceScalingmethodsOp
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_X_T = 1 << 4 , CTX_SET_X_TT = 1 << 6 , CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

Set Dirichlet boundary conditions on displacements.

Examples
bone_adaptation.cpp, cell_forces.cpp, elasticity_mixed_formulation.cpp, mortar_contact.cpp, navier_stokes.cpp, and nonlinear_dynamics.cpp.

Definition at line 69 of file DirichletBC.hpp.

Constructor & Destructor Documentation

◆ DirichletDisplacementBc() [1/2]

DirichletDisplacementBc::DirichletDisplacementBc ( MoFEM::Interface m_field,
const std::string &  field_name,
Mat  Aij,
Vec  X,
Vec  F,
string  blockset_name = "DISPLACEMENT" 
)

Definition at line 49 of file DirichletBC.cpp.

53 : mField(m_field), fieldName(field_name), blocksetName(blockset_name),
54 dIag(1) {
55 snes_B = Aij;
56 snes_x = X;
57 snes_f = F;
58 ts_B = Aij;
59 ts_u = X;
60 ts_F = F;
61};
constexpr auto field_name
const std::string blocksetName
Definition: DirichletBC.hpp:86
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:72
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:73
MoFEM::Interface & mField
Definition: DirichletBC.hpp:71
Vec & snes_f
residual
Vec & snes_x
state vector
Mat & snes_B
preconditioner of jacobian matrix
Vec & ts_F
residual vector
Vec & ts_u
state vector
Mat & ts_B
Preconditioner for ts_A.

◆ DirichletDisplacementBc() [2/2]

DirichletDisplacementBc::DirichletDisplacementBc ( MoFEM::Interface m_field,
const std::string &  field_name,
string  blockset_name = "DISPLACEMENT" 
)

Definition at line 63 of file DirichletBC.cpp.

66 : mField(m_field), fieldName(field_name), blocksetName(blockset_name),
67 dIag(1) {
68 snes_B = PETSC_NULL;
69 snes_x = PETSC_NULL;
70 snes_f = PETSC_NULL;
71 ts_B = PETSC_NULL;
72 ts_u = PETSC_NULL;
73 ts_F = PETSC_NULL;
74};

Member Function Documentation

◆ applyScaleBcData()

MoFEMErrorCode DirichletDisplacementBc::applyScaleBcData ( DataFromBc bc_data)

Definition at line 199 of file DirichletBC.cpp.

199 {
201
202 bc_data.scaled_values = bc_data.initial_values;
204 bc_data.scaled_values);
205 // for rotation
206 bc_data.theta = bc_data.scaled_values[0];
208}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
VectorDouble scaled_values
Definition: DirichletBC.hpp:44
VectorDouble initial_values
Definition: DirichletBC.hpp:45
double theta
Definition: DirichletBC.hpp:53
boost::ptr_vector< MethodForForceScaling > methodsOp
Definition: DirichletBC.hpp:88
static MoFEMErrorCode applyScale(const FEMethod *fe, boost::ptr_vector< MethodForForceScaling > &methods_op, VectorDouble &nf)

◆ calculateRotationForDof() [1/2]

MoFEMErrorCode DirichletDisplacementBc::calculateRotationForDof ( EntityHandle  ent,
DataFromBc bc_data 
)

Definition at line 190 of file DirichletBC.cpp.

191 {
193 VectorDouble3 coords(3);
194 CHKERR mField.get_moab().get_coords(&ent, 1, &*coords.begin());
195 CHKERR calculateRotationForDof(coords, bc_data);
197}
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:103
MoFEMErrorCode calculateRotationForDof(VectorDouble3 &coords, DataFromBc &bc_data)
Calculate displacements from rotation for particular dof.
virtual moab::Interface & get_moab()=0

◆ calculateRotationForDof() [2/2]

MoFEMErrorCode DirichletDisplacementBc::calculateRotationForDof ( VectorDouble3 &  coords,
DataFromBc bc_data 
)

Calculate displacements from rotation for particular dof.

Parameters
dof
bc_data
Returns
MoFEMErrorCode

Definition at line 180 of file DirichletBC.cpp.

181 {
183 if (bc_data.is_rotation)
184 bc_data.scaled_values = get_displacement(coords, bc_data.t_centr,
185 bc_data.t_normal, bc_data.theta);
187}
auto get_displacement(VectorDouble3 &coords, FTensor1 t_centr, FTensor1 t_normal, double theta)
FTensor::Tensor1< double, 3 > t_centr
Definition: DirichletBC.hpp:52
FTensor::Tensor1< double, 3 > t_normal
Definition: DirichletBC.hpp:51
bool is_rotation
Definition: DirichletBC.hpp:50

◆ getBcDataFromSetsAndBlocks()

MoFEMErrorCode DirichletDisplacementBc::getBcDataFromSetsAndBlocks ( std::vector< DataFromBc > &  bc_data)

Get the Bc Data From Sets And Blocks object Use DISPLACEMENT blockset name (default) with 6 atributes: 1,2,3 are values of displacements x,y,z 4,5,6 are flags for x,y,z (0 or 1)

Parameters
bc_data
Returns
MoFEMErrorCode

Definition at line 76 of file DirichletBC.cpp.

77 {
78
80
81 // Loop over meshsets with Dirichlet boundary condition on displacements
84 bc_data.push_back(DataFromBc());
86 CHKERR bc_data.back().getBcData(mydata, &(*it));
87 CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
88 }
89 // Loop over blocksets with DISPLACEMENT (default) boundary condition
91 if (it->getName().compare(0, blocksetName.length(), blocksetName) == 0) {
92 bc_data.push_back(DataFromBc());
93 std::vector<double> mydata;
94 CHKERR bc_data.back().getBcData(mydata, &(*it));
95 CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
96 }
97 }
98
100}
@ NODESET
Definition: definitions.h:159
@ DISPLACEMENTSET
Definition: definitions.h:163
@ BLOCKSET
Definition: definitions.h:161
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
#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.
Data from Cubit blocksets.
Definition: DirichletBC.hpp:43
Definition of the displacement bc data structure.
Definition: BCData.hpp:87

◆ getRotationBcFromBlock()

MoFEMErrorCode DirichletDisplacementBc::getRotationBcFromBlock ( std::vector< DataFromBc > &  bc_data)

Get the Rotation Bc From Block object Use ROTATION blockset name with 7 atributes: 1,2,3 are x,y,z coords of the center of rotation 4,5,6 are are angular velocities in x,y,z.

Parameters
bc_data
Returns
MoFEMErrorCode

Definition at line 145 of file DirichletBC.cpp.

146 {
147
149
151 if (it->getName().compare(0, 8, "ROTATION") == 0) {
152 bc_data.push_back(DataFromBc());
153 std::vector<double> mydata;
154 CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
155 CHKERR it->getAttributes(mydata);
156 if (mydata.size() < 6) {
157 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
158 "6 attributes are required for Rotation (3 center coords + 3 "
159 "angles, (+ 3 optional) flags for xyz)");
160 }
161 for (int ii : {0, 1, 2}) {
162 bc_data.back().bc_flags[ii] = 1;
163 bc_data.back().t_centr(ii) = mydata[ii + 0];
164 bc_data.back().t_normal(ii) = mydata[ii + 3];
165 }
166 if (mydata.size() > 8)
167 for (int ii : {0, 1, 2})
168 bc_data.back().bc_flags[ii] = mydata[6 + ii];
169
170 bc_data.back().scaled_values[0] = bc_data.back().initial_values[0] =
171 bc_data.back().t_normal.l2();
172 bc_data.back().is_rotation = true;
173 }
174 }
175
177}
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:44

◆ iNitialize()

MoFEMErrorCode DirichletDisplacementBc::iNitialize ( )
virtual

Reimplemented in AnalyticalDirichletBC::DirichletBC, DirichletSpatialPositionsBc, DirichletTemperatureBc, DirichletFixFieldAtEntitiesBc, and DirichletDisplacementRemoveDofsBc.

Definition at line 210 of file DirichletBC.cpp.

210 {
212 if (mapZeroRows.empty() || !methodsOp.empty()) {
213
214 std::vector<DataFromBc> bcData;
217
218 for (auto &bc_it : bcData) {
219
221
222 for (int dim = 0; dim < 3; dim++) {
223
224 auto for_each_dof = [&](auto &dof) {
226 if (dof->getEntType() == MBVERTEX) {
227 CHKERR calculateRotationForDof(dof->getEnt(), bc_it);
228 if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
229 // set boundary values to field data
230 mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[0];
231 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
232 }
233 if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
234 mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[1];
235 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
236 }
237 if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
238 mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[2];
239 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
240 }
241
242 } else {
243 if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
244 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
245 dof->getFieldData() = 0;
246 }
247 if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
248 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
249 dof->getFieldData() = 0;
250 }
251 if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
252 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
253 dof->getFieldData() = 0;
254 }
255 }
257 };
258
261 bc_it.bc_ents[dim], for_each_dof);
262 }
263 }
264 dofsIndices.resize(mapZeroRows.size());
265 dofsValues.resize(mapZeroRows.size());
266 int ii = 0;
267 auto mit = mapZeroRows.begin();
268 for (; mit != mapZeroRows.end(); mit++, ii++) {
269 dofsIndices[ii] = mit->first;
270 dofsValues[ii] = mit->second;
271 }
272 }
274}
static MoFEMErrorCode set_numered_dofs_on_ents(const Problem *problem_ptr, const FieldBitNumber bit_number, Range &ents, boost::function< MoFEMErrorCode(const boost::shared_ptr< MoFEM::NumeredDofEntity > &dof)> for_each_dof)
Definition: DirichletBC.cpp:23
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
const int dim
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:82
MoFEMErrorCode getRotationBcFromBlock(std::vector< DataFromBc > &bc_data)
Get the Rotation Bc From Block object Use ROTATION blockset name with 7 atributes: 1,...
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:83
MoFEMErrorCode getBcDataFromSetsAndBlocks(std::vector< DataFromBc > &bc_data)
Get the Bc Data From Sets And Blocks object Use DISPLACEMENT blockset name (default) with 6 atributes...
Definition: DirichletBC.cpp:76
MoFEMErrorCode applyScaleBcData(DataFromBc &bc_data)
std::vector< double > dofsValues
Definition: DirichletBC.hpp:84
const Problem * problemPtr
raw pointer to problem
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number

◆ operator()()

MoFEMErrorCode DirichletDisplacementBc::operator() ( )
virtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from MoFEM::BasicMethod.

Reimplemented in DirichletDisplacementRemoveDofsBc.

Definition at line 92 of file DirichletBC.hpp.

92{ return 0; }

◆ postProcess()

MoFEMErrorCode DirichletDisplacementBc::postProcess ( )
virtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::BasicMethod.

Reimplemented in DirichletFixFieldAtEntitiesBc, DirichletDisplacementRemoveDofsBc, and DirichletDensityBc.

Definition at line 293 of file DirichletBC.cpp.

293 {
295
296 if (snes_ctx == CTX_SNESNONE) {
297 if (snes_B) {
298 CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
299 CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
300 CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
301 dofsIndices.empty() ? PETSC_NULL
302 : &dofsIndices[0],
303 dIag, PETSC_NULL, PETSC_NULL);
304 }
305 if (snes_f) {
306 CHKERR VecAssemblyBegin(snes_f);
307 CHKERR VecAssemblyEnd(snes_f);
308 for (std::vector<int>::iterator vit = dofsIndices.begin();
309 vit != dofsIndices.end(); vit++) {
310 CHKERR VecSetValue(snes_f, *vit, 0, INSERT_VALUES);
311 }
312 CHKERR VecAssemblyBegin(snes_f);
313 CHKERR VecAssemblyEnd(snes_f);
314 }
315 }
316
317 switch (snes_ctx) {
318 case CTX_SNESNONE:
319 break;
320 case CTX_SNESSETFUNCTION: {
321 if (!dofsIndices.empty()) {
322 dofsXValues.resize(dofsIndices.size());
323 const double *a_snes_x;
324 CHKERR VecGetArrayRead(snes_x, &a_snes_x);
325 auto &dofs_by_glob_idx =
327 int idx = 0;
328 for (auto git : dofsIndices) {
329 auto dof_it = dofs_by_glob_idx.find(git);
330 if (dof_it != dofs_by_glob_idx.end()) {
331 dofsXValues[idx] = a_snes_x[dof_it->get()->getPetscLocalDofIdx()];
332 ++idx;
333 } else
334 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
335 "Dof with global %d id not found", git);
336 }
337 CHKERR VecRestoreArrayRead(snes_x, &a_snes_x);
338 }
339 CHKERR VecAssemblyBegin(snes_f);
340 CHKERR VecAssemblyEnd(snes_f);
341
342 if (!dofsIndices.empty()) {
343 int ii = 0;
344 for (std::vector<int>::iterator vit = dofsIndices.begin();
345 vit != dofsIndices.end(); vit++, ii++) {
346 double val = 0;
347 if (!dofsXValues.empty()) {
348 val += dofsXValues[ii];
349 val += -mapZeroRows[*vit]; // in snes it is on the left hand side,
350 // that way -1
351 dofsXValues[ii] = val;
352 }
353 }
355 snes_f, dofsIndices.size(),
356 dofsIndices.empty() ? PETSC_NULL : &*dofsIndices.begin(),
357 dofsXValues.empty() ? PETSC_NULL : &*dofsXValues.begin(),
358 INSERT_VALUES);
359 }
360 CHKERR VecAssemblyBegin(snes_f);
361 CHKERR VecAssemblyEnd(snes_f);
362 } break;
363 case CTX_SNESSETJACOBIAN: {
364
365 CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
366 CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
367 CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
368 dofsIndices.empty() ? PETSC_NULL
369 : &*dofsIndices.begin(),
370 dIag, PETSC_NULL, PETSC_NULL);
371
372 } break;
373 default:
374 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown snes stage");
375 }
376
378}
std::vector< double > dofsXValues
Definition: DirichletBC.hpp:85
auto & getNumeredColDofsPtr() const
get access to numeredColDofsPtr storing DOFs on cols
SNESContext snes_ctx

◆ preProcess()

MoFEMErrorCode DirichletDisplacementBc::preProcess ( )
virtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from MoFEM::BasicMethod.

Reimplemented in DirichletFixFieldAtEntitiesBc, and DirichletDisplacementRemoveDofsBc.

Definition at line 276 of file DirichletBC.cpp.

276 {
278
280
281 if (snes_ctx == CTX_SNESNONE) {
282 if (!dofsIndices.empty()) {
284 &*dofsValues.begin(), INSERT_VALUES);
285 }
286 CHKERR VecAssemblyBegin(snes_x);
287 CHKERR VecAssemblyEnd(snes_x);
288 }
289
291}
virtual MoFEMErrorCode iNitialize()

Member Data Documentation

◆ blocksetName

const std::string DirichletDisplacementBc::blocksetName

Definition at line 86 of file DirichletBC.hpp.

◆ dIag

double DirichletDisplacementBc::dIag

diagonal value set on zeroed column and rows

Definition at line 73 of file DirichletBC.hpp.

◆ dofsIndices

std::vector<int> DirichletDisplacementBc::dofsIndices

Definition at line 83 of file DirichletBC.hpp.

◆ dofsValues

std::vector<double> DirichletDisplacementBc::dofsValues

Definition at line 84 of file DirichletBC.hpp.

◆ dofsXValues

std::vector<double> DirichletDisplacementBc::dofsXValues

Definition at line 85 of file DirichletBC.hpp.

◆ fieldName

const std::string DirichletDisplacementBc::fieldName

field name to set Dirichlet BC

Definition at line 72 of file DirichletBC.hpp.

◆ mapZeroRows

std::map<DofIdx, FieldData> DirichletDisplacementBc::mapZeroRows

Definition at line 82 of file DirichletBC.hpp.

◆ methodsOp

boost::ptr_vector<MethodForForceScaling> DirichletDisplacementBc::methodsOp
Examples
bone_adaptation.cpp.

Definition at line 88 of file DirichletBC.hpp.

◆ mField

MoFEM::Interface& DirichletDisplacementBc::mField

Definition at line 71 of file DirichletBC.hpp.


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