v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
DirichletSpatialPositionsBc Struct Reference

Set Dirichlet boundary conditions on spatial displacements. More...

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

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

Public Member Functions

 DirichletSpatialPositionsBc (MoFEM::Interface &m_field, const std::string &field_name, Mat aij, Vec x, Vec f, const std::string material_positions="MESH_NODE_POSITIONS", const std::string blockset_name="DISPLACEMENT")
 
 DirichletSpatialPositionsBc (MoFEM::Interface &m_field, const std::string &field_name, const std::string material_positions="MESH_NODE_POSITIONS", const std::string blockset_name="DISPLACEMENT")
 
MoFEMErrorCode iNitialize ()
 
- Public Member Functions inherited from DirichletDisplacementBc
 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

std::string materialPositions
 
std::vector< std::string > fixFields
 
VectorDouble cOords
 
- Public Attributes inherited from DirichletDisplacementBc
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 number 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...
 
PetscReal ts_dt
 time step size 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 spatial displacements.

Examples
mortar_contact_thermal.cpp, simple_contact.cpp, and simple_contact_thermal.cpp.

Definition at line 211 of file DirichletBC.hpp.

Constructor & Destructor Documentation

◆ DirichletSpatialPositionsBc() [1/2]

DirichletSpatialPositionsBc::DirichletSpatialPositionsBc ( MoFEM::Interface m_field,
const std::string &  field_name,
Mat  aij,
Vec  x,
Vec  f,
const std::string  material_positions = "MESH_NODE_POSITIONS",
const std::string  blockset_name = "DISPLACEMENT" 
)
inline

Definition at line 213 of file DirichletBC.hpp.

217 : DirichletDisplacementBc(m_field, field_name, aij, x, f, blockset_name),
218 materialPositions(material_positions) {}
constexpr auto field_name
Set Dirichlet boundary conditions on displacements.
Definition: DirichletBC.hpp:57

◆ DirichletSpatialPositionsBc() [2/2]

DirichletSpatialPositionsBc::DirichletSpatialPositionsBc ( MoFEM::Interface m_field,
const std::string &  field_name,
const std::string  material_positions = "MESH_NODE_POSITIONS",
const std::string  blockset_name = "DISPLACEMENT" 
)
inline

Definition at line 220 of file DirichletBC.hpp.

224 : DirichletDisplacementBc(m_field, field_name, blockset_name),
225 materialPositions(material_positions) {}

Member Function Documentation

◆ iNitialize()

MoFEMErrorCode DirichletSpatialPositionsBc::iNitialize ( )
virtual

Reimplemented from DirichletDisplacementBc.

Definition at line 368 of file DirichletBC.cpp.

368 {
370
371 if (mapZeroRows.empty() || !methodsOp.empty()) {
372
373 std::vector<DataFromBc> bcData;
375
376 const FieldEntity_multiIndex *field_ents;
377 CHKERR mField.get_field_ents(&field_ents);
378 auto &field_ents_by_uid = field_ents->get<Unique_mi_tag>();
379
380 VectorDouble3 coords(3);
381
382 for (auto &bc_it : bcData) {
383
385
386 for (int dim = 0; dim < 3; dim++) {
387
388 auto for_each_dof = [&](auto &dof) {
390
391 if (!dim) {
392
393 EntityHandle node = dof->getEnt();
394 if (!dof->getDofCoeffIdx()) {
395 auto eit =
396 field_ents_by_uid.find(FieldEntity::getLocalUniqueIdCalculate(
398 if (eit != field_ents_by_uid.end())
399 noalias(coords) = (*eit)->getEntFieldData();
400 else
401 CHKERR mField.get_moab().get_coords(&node, 1,
402 &*coords.data().begin());
403 }
404 CHKERR calculateRotationForDof(coords, bc_it);
405 if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
406 mapZeroRows[dof->getPetscGlobalDofIdx()] =
407 coords[0] + bc_it.scaled_values[0];
408 // set boundary values to field data
409 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
410 }
411 if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
412 mapZeroRows[dof->getPetscGlobalDofIdx()] =
413 coords[1] + bc_it.scaled_values[1];
414 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
415 }
416 if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
417 mapZeroRows[dof->getPetscGlobalDofIdx()] =
418 coords[2] + bc_it.scaled_values[2];
419 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
420 }
421 } else {
422 if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
423 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
424 dof->getFieldData() = 0;
425 }
426 if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
427 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
428 dof->getFieldData() = 0;
429 }
430 if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
431 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
432 dof->getFieldData() = 0;
433 }
434 }
435
437 };
438
441 bc_it.bc_ents[dim], for_each_dof);
442
443 auto fix_field_dof = [&](auto &dof) {
445 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
447 };
448
449 for (auto &fix_field : fixFields) {
452 bc_it.bc_ents[dim], fix_field_dof);
453 }
454 }
455 }
456
457 // set vector of values and indices
458 dofsIndices.resize(mapZeroRows.size());
459 dofsValues.resize(mapZeroRows.size());
460 auto mit = mapZeroRows.begin();
461 for (int ii = 0; mit != mapZeroRows.end(); mit++, ii++) {
462 dofsIndices[ii] = mit->first;
463 dofsValues[ii] = mit->second;
464 }
465 }
467}
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:11
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#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
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
const int dim
multi_index_container< boost::shared_ptr< FieldEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, member< FieldEntity, UId, &FieldEntity::localUId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FieldEntity::interface_type_RefEntity, EntityHandle, &FieldEntity::getEnt > > > > FieldEntity_multiIndex
MultiIndex container keeps FieldEntity.
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:70
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:60
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:71
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:64
MoFEMErrorCode applyScaleBcData(DataFromBc &bc_data)
boost::ptr_vector< MethodForForceScaling > methodsOp
Definition: DirichletBC.hpp:76
std::vector< double > dofsValues
Definition: DirichletBC.hpp:72
MoFEM::Interface & mField
Definition: DirichletBC.hpp:59
MoFEMErrorCode calculateRotationForDof(VectorDouble3 &coords, DataFromBc &bc_data)
Calculate displacements from rotation for particular dof.
std::vector< std::string > fixFields
const Problem * problemPtr
raw pointer to problem
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.

Member Data Documentation

◆ cOords

VectorDouble DirichletSpatialPositionsBc::cOords

Definition at line 231 of file DirichletBC.hpp.

◆ fixFields

std::vector<std::string> DirichletSpatialPositionsBc::fixFields

Definition at line 229 of file DirichletBC.hpp.

◆ materialPositions

std::string DirichletSpatialPositionsBc::materialPositions

name of the field with reference material positions

Definition at line 227 of file DirichletBC.hpp.


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