v0.13.2
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 ()
 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 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 displacements.

Examples
bone_adaptation.cpp, cell_forces.cpp, elasticity_mixed_formulation.cpp, navier_stokes.cpp, and nonlinear_dynamics.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
Definition: DirichletBC.hpp:74
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:60
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:61
MoFEM::Interface & mField
Definition: DirichletBC.hpp:59
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 51 of file DirichletBC.cpp.

54 : mField(m_field), fieldName(field_name), blocksetName(blockset_name),
55 dIag(1) {
56 snes_B = PETSC_NULL;
57 snes_x = PETSC_NULL;
58 snes_f = PETSC_NULL;
59 ts_B = PETSC_NULL;
60 ts_u = PETSC_NULL;
61 ts_F = PETSC_NULL;
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 ...
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
VectorDouble scaled_values
Definition: DirichletBC.hpp:32
VectorDouble initial_values
Definition: DirichletBC.hpp:33
double theta
Definition: DirichletBC.hpp:41
boost::ptr_vector< MethodForForceScaling > methodsOp
Definition: DirichletBC.hpp:76
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
Definition: DirichletBC.hpp:40
FTensor::Tensor1< double, 3 > t_normal
Definition: DirichletBC.hpp:39
bool is_rotation
Definition: DirichletBC.hpp:38

◆ 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
Definition: definitions.h:146
@ DISPLACEMENTSET
Definition: definitions.h:150
@ BLOCKSET
Definition: definitions.h:148
#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:31
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)
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 MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
const int dim
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:70
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: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)
std::vector< double > dofsValues
Definition: DirichletBC.hpp:72
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() ( )
inlinevirtual

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 80 of file DirichletBC.hpp.

80{ 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 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_NULL
290 : &dofsIndices[0],
291 dIag, PETSC_NULL, PETSC_NULL);
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 SETERRQ1(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_NULL : &*dofsIndices.begin(),
345 dofsXValues.empty() ? PETSC_NULL : &*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_NULL
357 : &*dofsIndices.begin(),
358 dIag, PETSC_NULL, PETSC_NULL);
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
Definition: DirichletBC.hpp:73
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 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
Examples
bone_adaptation.cpp.

Definition at line 76 of file DirichletBC.hpp.

◆ mField

MoFEM::Interface& DirichletDisplacementBc::mField

Definition at line 59 of file DirichletBC.hpp.


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