v0.15.0
Loading...
Searching...
No Matches
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 ()
 Pre-processing function executed at loop initialization.
 
MoFEMErrorCode operator() ()
 Main operator function executed for each loop iteration.
 
MoFEMErrorCode postProcess ()
 Post-processing function executed at loop completion.
 
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)
 
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.
 
MoFEMErrorCode calculateRotationForDof (VectorDouble3 &coords, DataFromBc &bc_data)
 Calculate displacements from rotation for particular dof.
 
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
 Query interface for type casting.
 
 FEMethod ()=default
 Default constructor.
 
auto getFEName () const
 Get the name of the current finite element.
 
auto getDataDofsPtr () const
 Get pointer to DOF data for the current finite element.
 
auto getDataVectorDofsPtr () const
 Get pointer to vector DOF data for the current finite element.
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 Get reference to data field entities for the current finite element.
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 Get shared pointer to data field entities for the current finite element.
 
auto & getRowFieldEnts () const
 Get reference to row field entities for the current finite element.
 
auto & getRowFieldEntsPtr () const
 Get shared pointer to row field entities for the current finite element.
 
auto & getColFieldEnts () const
 Get reference to column field entities for the current finite element.
 
auto & getColFieldEntsPtr () const
 Get shared pointer to column field entities for the current finite element.
 
auto getRowDofsPtr () const
 Get pointer to row DOFs for the current finite element.
 
auto getColDofsPtr () const
 Get pointer to column DOFs for the current finite element.
 
auto getNumberOfNodes () const
 Get the number of nodes in the current finite element.
 
EntityHandle getFEEntityHandle () const
 Get the entity handle of the current finite element.
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 Get nodal data for a specific field.
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 Default constructor.
 
virtual ~BasicMethod ()=default
 Virtual destructor.
 
int getNinTheLoop () const
 Get current loop iteration index.
 
int getLoopSize () const
 Get total loop size.
 
auto getLoHiFERank () const
 Get processor rank range for finite element iteration.
 
auto getLoFERank () const
 Get lower processor rank for finite element iteration.
 
auto getHiFERank () const
 Get upper processor rank for finite element iteration.
 
unsigned int getFieldBitNumber (std::string field_name) const
 Get bit number for a specific field by name.
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from another BasicMethod instance.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak pointer object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 KspMethod ()
 Default constructor.
 
virtual ~KspMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 Copy data from another KSP method.
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 Default constructor.
 
virtual ~PetscData ()=default
 Virtual destructor.
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 Copy PETSc data from another instance.
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 SnesMethod ()
 Default constructor.
 
virtual ~SnesMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy SNES data from another instance.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TSMethod ()
 Default constructor.
 
virtual ~TSMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data from another instance.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Query interface for type casting.
 
 TaoMethod ()
 Default constructor.
 
virtual ~TaoMethod ()=default
 Virtual destructor.
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data from another instance.
 

Public Attributes

MoFEM::InterfacemField
 
const std::string fieldName
 field name to set Dirichlet BC
 
double dIag
 diagonal value set on zeroed column and rows
 
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 the finite element being processed.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 Shared pointer to finite element database structure.
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Test function to determine if element should be skipped.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 Current index of processed method in the loop.
 
int loopSize
 Total number of methods to process in the loop.
 
std::pair< int, int > loHiFERank
 Processor rank range for distributed finite element iteration.
 
int rAnk
 Current processor rank in MPI communicator.
 
int sIze
 Total number of processors in MPI communicator.
 
const RefEntity_multiIndexrefinedEntitiesPtr
 Pointer to container of refined MoFEM DOF entities.
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 Pointer to container of refined finite element entities.
 
const ProblemproblemPtr
 Raw pointer to current MoFEM problem instance.
 
const Field_multiIndexfieldsPtr
 Raw pointer to fields multi-index container.
 
const FieldEntity_multiIndexentitiesPtr
 Raw pointer to container of field entities.
 
const DofEntity_multiIndexdofsPtr
 Raw pointer to container of degree of freedom entities.
 
const FiniteElement_multiIndexfiniteElementsPtr
 Raw pointer to container of finite elements.
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 Raw pointer to container of finite element entities.
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 Raw pointer to container of adjacencies between DOFs and finite elements.
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing operations.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing operations.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for main operator execution.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 Switch for vector assembly operations.
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 Switch for matrix assembly operations.
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 Weak pointer to cached entity data.
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Current KSP computation context.
 
KSP ksp
 PETSc KSP linear solver object.
 
Vec & ksp_f
 Reference to residual vector in KSP context.
 
Mat & ksp_A
 Reference to system matrix in KSP context.
 
Mat & ksp_B
 Reference to preconditioner matrix in KSP context.
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 Current data context switches.
 
Vec f
 PETSc residual vector.
 
Mat A
 PETSc Jacobian matrix.
 
Mat B
 PETSc preconditioner matrix.
 
Vec x
 PETSc solution vector.
 
Vec dx
 PETSc solution increment vector.
 
Vec x_t
 PETSc first time derivative vector.
 
Vec x_tt
 PETSc second time derivative vector.
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 Current SNES computation context.
 
SNES snes
 PETSc SNES nonlinear solver object.
 
Vec & snes_x
 Reference to current solution state vector.
 
Vec & snes_dx
 Reference to solution update/increment vector.
 
Vec & snes_f
 Reference to residual vector.
 
Mat & snes_A
 Reference to Jacobian matrix.
 
Mat & snes_B
 Reference to preconditioner of Jacobian matrix.
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 PETSc time stepping solver object.
 
TSContext ts_ctx
 Current TS computation context.
 
PetscInt ts_step
 Current time step number.
 
PetscReal ts_a
 Shift parameter for U_t (see PETSc Time Solver documentation)
 
PetscReal ts_aa
 Shift parameter for U_tt (second time derivative)
 
PetscReal ts_t
 Current time value.
 
PetscReal ts_dt
 Current time step size.
 
Vec & ts_u
 Reference to current state vector U(t)
 
Vec & ts_u_t
 Reference to first time derivative of state vector dU/dt.
 
Vec & ts_u_tt
 Reference to second time derivative of state vector d²U/dt²
 
Vec & ts_F
 Reference to residual vector F(t,U,U_t,U_tt)
 
Mat & ts_A
 Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
 
Mat & ts_B
 Reference to preconditioner matrix for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 Current TAO computation context.
 
Tao tao
 PETSc TAO optimization solver object.
 
Vec & tao_x
 Reference to optimization variables vector.
 
Vec & tao_f
 Reference to gradient vector.
 
Mat & tao_A
 Reference to Hessian matrix.
 
Mat & tao_B
 Reference to preconditioner matrix for Hessian.
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 Context enumeration for KSP solver phases. 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_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 Enumeration for data context flags. More...
 
using Switches = std::bitset< 8 >
 Bitset type for context switches.
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 Context enumeration for SNES solver phases. More...
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 Context enumeration for TS solver phases. More...
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 Context enumeration for TAO solver phases. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 No data switch.
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 Residual vector switch.
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 Jacobian matrix switch.
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 Preconditioner matrix switch.
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 Solution vector switch.
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 Solution increment switch.
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 First time derivative switch.
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 Second time derivative switch.
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 Time value switch.
 

Detailed Description

Set Dirichlet boundary conditions on displacements.

Examples
mofem/tutorials/cor-10/navier_stokes.cpp, mofem/tutorials/cor-7/elasticity_mixed_formulation.cpp, mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp, and mofem/users_modules/bone_remodelling/bone_adaptation.cpp.

Definition at line 57 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 37 of file DirichletBC.cpp.

41 : mField(m_field), fieldName(field_name), blocksetName(blockset_name),
42 dIag(1) {
43 snes_B = Aij;
44 snes_x = X;
45 snes_f = F;
46 ts_B = Aij;
47 ts_u = X;
48 ts_F = F;
49};
@ F
constexpr auto field_name
const std::string blocksetName
const std::string fieldName
field name to set Dirichlet BC
double dIag
diagonal value set on zeroed column and rows
MoFEM::Interface & mField
Vec & snes_f
Reference to residual vector.
Vec & snes_x
Reference to current solution state vector.
Mat & snes_B
Reference to preconditioner of Jacobian matrix.
Vec & ts_F
Reference to residual vector F(t,U,U_t,U_tt)
Vec & ts_u
Reference to current state vector U(t)
Mat & ts_B
Reference to preconditioner matrix for ts_A.

◆ DirichletDisplacementBc() [2/2]

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

Definition at line 51 of file DirichletBC.cpp.

54 : mField(m_field), fieldName(field_name), blocksetName(blockset_name),
55 dIag(1) {
56 snes_B = PETSC_NULLPTR;
57 snes_x = PETSC_NULLPTR;
58 snes_f = PETSC_NULLPTR;
59 ts_B = PETSC_NULLPTR;
60 ts_u = PETSC_NULLPTR;
61 ts_F = PETSC_NULLPTR;
62};

Member Function Documentation

◆ applyScaleBcData()

MoFEMErrorCode DirichletDisplacementBc::applyScaleBcData ( DataFromBc bc_data)

Definition at line 187 of file DirichletBC.cpp.

187 {
189
190 bc_data.scaled_values = bc_data.initial_values;
192 bc_data.scaled_values);
193 // for rotation
194 bc_data.theta = bc_data.scaled_values[0];
196}
#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()
#define CHKERR
Inline error check.
VectorDouble scaled_values
VectorDouble initial_values
boost::ptr_vector< MethodForForceScaling > methodsOp
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 178 of file DirichletBC.cpp.

179 {
181 VectorDouble3 coords(3);
182 CHKERR mField.get_moab().get_coords(&ent, 1, &*coords.begin());
183 CHKERR calculateRotationForDof(coords, bc_data);
185}
VectorBoundedArray< double, 3 > VectorDouble3
Definition Types.hpp:92
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 168 of file DirichletBC.cpp.

169 {
171 if (bc_data.is_rotation)
172 bc_data.scaled_values = get_displacement(coords, bc_data.t_centr,
173 bc_data.t_normal, bc_data.theta);
175}
auto get_displacement(VectorDouble3 &coords, FTensor1 t_centr, FTensor1 t_normal, double theta)
FTensor::Tensor1< double, 3 > t_centr
FTensor::Tensor1< double, 3 > t_normal

◆ 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 64 of file DirichletBC.cpp.

65 {
66
68
69 // Loop over meshsets with Dirichlet boundary condition on displacements
72 bc_data.push_back(DataFromBc());
74 CHKERR bc_data.back().getBcData(mydata, &(*it));
75 CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
76 }
77 // Loop over blocksets with DISPLACEMENT (default) boundary condition
79 if (it->getName().compare(0, blocksetName.length(), blocksetName) == 0) {
80 bc_data.push_back(DataFromBc());
81 std::vector<double> mydata;
82 CHKERR bc_data.back().getBcData(mydata, &(*it));
83 CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
84 }
85 }
86
88}
@ NODESET
@ DISPLACEMENTSET
@ BLOCKSET
#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 of the displacement bc data structure.
Definition BCData.hpp:72

◆ 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 133 of file DirichletBC.cpp.

134 {
135
137
139 if (it->getName().compare(0, 8, "ROTATION") == 0) {
140 bc_data.push_back(DataFromBc());
141 std::vector<double> mydata;
142 CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
143 CHKERR it->getAttributes(mydata);
144 if (mydata.size() < 6) {
145 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
146 "6 attributes are required for Rotation (3 center coords + 3 "
147 "angles, (+ 3 optional) flags for xyz)");
148 }
149 for (int ii : {0, 1, 2}) {
150 bc_data.back().bc_flags[ii] = 1;
151 bc_data.back().t_centr(ii) = mydata[ii + 0];
152 bc_data.back().t_normal(ii) = mydata[ii + 3];
153 }
154 if (mydata.size() > 8)
155 for (int ii : {0, 1, 2})
156 bc_data.back().bc_flags[ii] = mydata[6 + ii];
157
158 bc_data.back().scaled_values[0] = bc_data.back().initial_values[0] =
159 bc_data.back().t_normal.l2();
160 bc_data.back().is_rotation = true;
161 }
162 }
163
165}
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31

◆ iNitialize()

MoFEMErrorCode DirichletDisplacementBc::iNitialize ( )
virtual

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

Definition at line 198 of file DirichletBC.cpp.

198 {
200 if (mapZeroRows.empty() || !methodsOp.empty()) {
201
202 std::vector<DataFromBc> bcData;
205
206 for (auto &bc_it : bcData) {
207
209
210 for (int dim = 0; dim < 3; dim++) {
211
212 auto for_each_dof = [&](auto &dof) {
214 if (dof->getEntType() == MBVERTEX) {
215 CHKERR calculateRotationForDof(dof->getEnt(), bc_it);
216 if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
217 // set boundary values to field data
218 mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[0];
219 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
220 }
221 if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
222 mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[1];
223 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
224 }
225 if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
226 mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[2];
227 dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
228 }
229
230 } else {
231 if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
232 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
233 dof->getFieldData() = 0;
234 }
235 if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
236 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
237 dof->getFieldData() = 0;
238 }
239 if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
240 mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
241 dof->getFieldData() = 0;
242 }
243 }
245 };
246
249 bc_it.bc_ents[dim], for_each_dof);
250 }
251 }
252 dofsIndices.resize(mapZeroRows.size());
253 dofsValues.resize(mapZeroRows.size());
254 int ii = 0;
255 auto mit = mapZeroRows.begin();
256 for (; mit != mapZeroRows.end(); mit++, ii++) {
257 dofsIndices[ii] = mit->first;
258 dofsValues[ii] = mit->second;
259 }
260 }
262}
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)
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
std::map< DofIdx, FieldData > mapZeroRows
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
MoFEMErrorCode getBcDataFromSetsAndBlocks(std::vector< DataFromBc > &bc_data)
Get the Bc Data From Sets And Blocks object Use DISPLACEMENT blockset name (default) with 6 atributes...
MoFEMErrorCode applyScaleBcData(DataFromBc &bc_data)
std::vector< double > dofsValues
const Problem * problemPtr
Raw pointer to current MoFEM problem instance.
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number

◆ operator()()

MoFEMErrorCode DirichletDisplacementBc::operator() ( )
inlinevirtual

Main operator function executed for each loop iteration.

This virtual function is called for every item (finite element, entity, etc.) in the loop. It is the core computation function typically used for:

  • Calculating element local matrices and vectors
  • Assembling contributions to global system
  • Element-level post-processing operations
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Reimplemented in DirichletDisplacementRemoveDofsBc.

Definition at line 80 of file DirichletBC.hpp.

80{ return 0; }

◆ postProcess()

MoFEMErrorCode DirichletDisplacementBc::postProcess ( )
virtual

Post-processing function executed at loop completion.

This virtual function is called once at the end of the loop. It is typically used for:

  • Assembling matrices and vectors to global system
  • Computing global variables (e.g., total internal energy)
  • Finalizing computation results
  • Cleanup operations

Example of iterating over DOFs:

for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) {
// Process each DOF
}
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Reimplemented in DirichletFixFieldAtEntitiesBc, and DirichletDisplacementRemoveDofsBc.

Definition at line 281 of file DirichletBC.cpp.

281 {
283
284 if (snes_ctx == CTX_SNESNONE) {
285 if (snes_B) {
286 CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
287 CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
288 CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
289 dofsIndices.empty() ? PETSC_NULLPTR
290 : &dofsIndices[0],
291 dIag, PETSC_NULLPTR, PETSC_NULLPTR);
292 }
293 if (snes_f) {
294 CHKERR VecAssemblyBegin(snes_f);
295 CHKERR VecAssemblyEnd(snes_f);
296 for (std::vector<int>::iterator vit = dofsIndices.begin();
297 vit != dofsIndices.end(); vit++) {
298 CHKERR VecSetValue(snes_f, *vit, 0, INSERT_VALUES);
299 }
300 CHKERR VecAssemblyBegin(snes_f);
301 CHKERR VecAssemblyEnd(snes_f);
302 }
303 }
304
305 switch (snes_ctx) {
306 case CTX_SNESNONE:
307 break;
308 case CTX_SNESSETFUNCTION: {
309 if (!dofsIndices.empty()) {
310 dofsXValues.resize(dofsIndices.size());
311 const double *a_snes_x;
312 CHKERR VecGetArrayRead(snes_x, &a_snes_x);
313 auto &dofs_by_glob_idx =
315 int idx = 0;
316 for (auto git : dofsIndices) {
317 auto dof_it = dofs_by_glob_idx.find(git);
318 if (dof_it != dofs_by_glob_idx.end()) {
319 dofsXValues[idx] = a_snes_x[dof_it->get()->getPetscLocalDofIdx()];
320 ++idx;
321 } else
322 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
323 "Dof with global %d id not found", git);
324 }
325 CHKERR VecRestoreArrayRead(snes_x, &a_snes_x);
326 }
327 CHKERR VecAssemblyBegin(snes_f);
328 CHKERR VecAssemblyEnd(snes_f);
329
330 if (!dofsIndices.empty()) {
331 int ii = 0;
332 for (std::vector<int>::iterator vit = dofsIndices.begin();
333 vit != dofsIndices.end(); vit++, ii++) {
334 double val = 0;
335 if (!dofsXValues.empty()) {
336 val += dofsXValues[ii];
337 val += -mapZeroRows[*vit]; // in snes it is on the left hand side,
338 // that way -1
339 dofsXValues[ii] = val;
340 }
341 }
343 snes_f, dofsIndices.size(),
344 dofsIndices.empty() ? PETSC_NULLPTR : &*dofsIndices.begin(),
345 dofsXValues.empty() ? PETSC_NULLPTR : &*dofsXValues.begin(),
346 INSERT_VALUES);
347 }
348 CHKERR VecAssemblyBegin(snes_f);
349 CHKERR VecAssemblyEnd(snes_f);
350 } break;
351 case CTX_SNESSETJACOBIAN: {
352
353 CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
354 CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
355 CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
356 dofsIndices.empty() ? PETSC_NULLPTR
357 : &*dofsIndices.begin(),
358 dIag, PETSC_NULLPTR, PETSC_NULLPTR);
359
360 } break;
361 default:
362 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown snes stage");
363 }
364
366}
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
std::vector< double > dofsXValues
auto & getNumeredColDofsPtr() const
get access to numeredColDofsPtr storing DOFs on cols
@ CTX_SNESSETJACOBIAN
Setting up Jacobian matrix computation.
@ CTX_SNESSETFUNCTION
Setting up nonlinear function evaluation.
@ CTX_SNESNONE
No specific SNES context.
SNESContext snes_ctx
Current SNES computation context.

◆ preProcess()

MoFEMErrorCode DirichletDisplacementBc::preProcess ( )
virtual

Pre-processing function executed at loop initialization.

This virtual function is called once at the beginning of the loop. It is typically used for:

  • Zeroing matrices and vectors
  • Computing shape functions on reference elements
  • Preprocessing boundary conditions
  • Setting up initial computation state
Returns
Error code indicating success or failure

Reimplemented from MoFEM::BasicMethod.

Reimplemented in DirichletFixFieldAtEntitiesBc, and DirichletDisplacementRemoveDofsBc.

Definition at line 264 of file DirichletBC.cpp.

264 {
266
268
269 if (snes_ctx == CTX_SNESNONE) {
270 if (!dofsIndices.empty()) {
272 &*dofsValues.begin(), INSERT_VALUES);
273 }
274 CHKERR VecAssemblyBegin(snes_x);
275 CHKERR VecAssemblyEnd(snes_x);
276 }
277
279}
virtual MoFEMErrorCode iNitialize()

Member Data Documentation

◆ blocksetName

const std::string DirichletDisplacementBc::blocksetName

Definition at line 74 of file DirichletBC.hpp.

◆ dIag

double DirichletDisplacementBc::dIag

diagonal value set on zeroed column and rows

Definition at line 61 of file DirichletBC.hpp.

◆ dofsIndices

std::vector<int> DirichletDisplacementBc::dofsIndices

Definition at line 71 of file DirichletBC.hpp.

◆ dofsValues

std::vector<double> DirichletDisplacementBc::dofsValues

Definition at line 72 of file DirichletBC.hpp.

◆ dofsXValues

std::vector<double> DirichletDisplacementBc::dofsXValues

Definition at line 73 of file DirichletBC.hpp.

◆ fieldName

const std::string DirichletDisplacementBc::fieldName

field name to set Dirichlet BC

Definition at line 60 of file DirichletBC.hpp.

◆ mapZeroRows

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

Definition at line 70 of file DirichletBC.hpp.

◆ methodsOp

boost::ptr_vector<MethodForForceScaling> DirichletDisplacementBc::methodsOp

◆ mField

MoFEM::Interface& DirichletDisplacementBc::mField

Definition at line 59 of file DirichletBC.hpp.


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