v0.10.0
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 iNitalize ()
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop 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 - is the rotation angle 2,3,4 are x,y,z coords of the center of rotation 5,6,7 are x,y,z coords of the normal of rotation. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
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
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. More...
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size 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 operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of 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
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
DEPRECATED MoFEMErrorCode setSnesCtx (SNESContext ctx)
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
DEPRECATED MoFEMErrorCode setTsCtx (TSContext ctx)
 

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
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process 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
 
- 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 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)
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

Set Dirichlet boundary conditions on displacements.

Examples
bone_adaptation.cpp, cell_forces.cpp, and elasticity_mixed_formulation.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" 
)

◆ DirichletDisplacementBc() [2/2]

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

Member Function Documentation

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

74  {
75 
77 
78  // sets kinetic boundary conditions by blockset.
79  bool flag_cubit_disp = false;
81  mField, NODESET | DISPLACEMENTSET, it)) {
82  flag_cubit_disp = true;
83  }
84  // Loop over meshsets with Dirichlet boundary condition on displacements
85  if (flag_cubit_disp) {
87  mField, NODESET | DISPLACEMENTSET, it)) {
88  bc_data.push_back(DataFromBc());
90  CHKERR bc_data.back().getBcData(mydata, &(*it));
91  CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
92  }
93  } else {
95  if (it->getName().compare(0, blocksetName.length(), blocksetName) == 0) {
96  bc_data.push_back(DataFromBc());
97  std::vector<double> mydata;
98  CHKERR bc_data.back().getBcData(mydata, &(*it));
99  CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
100  }
101  }
102  }
103 
105 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
MoFEM::Interface & mField
Definition: DirichletBC.hpp:71
Data from Cubit blocksets.
Definition: DirichletBC.hpp:43
Definition of the displacement bc data structure.
Definition: BCData.hpp:87
#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.
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
const std::string blocksetName
Definition: DirichletBC.hpp:86

◆ getRotationBcFromBlock()

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

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

Parameters
bc_data
Returns
MoFEMErrorCode

Definition at line 153 of file DirichletBC.cpp.

154  {
155 
157 
159  if (it->getName().compare(0, 8, "ROTATION") == 0) {
160  bc_data.push_back(DataFromBc());
161  std::vector<double> mydata;
162  CHKERR bc_data.back().getEntitiesFromBc(mField, &(*it));
163  CHKERR it->getAttributes(mydata);
164  if (mydata.size() < 6) {
165  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
166  "7 attributes are required for Rotation (angle + 3 center "
167  "coords + 3 normal coords)");
168  }
169  for (int ii : {0, 1, 2}) {
170  bc_data.back().bc_flags[ii] = 1;
171  bc_data.back().t_centr(ii) = mydata[ii + 1];
172  bc_data.back().t_normal(ii) = mydata[ii + 4];
173 
174  }
175  bc_data.back().scaled_values[0] = mydata[0];
176  bc_data.back().is_rotation = true;
177  }
178  }
179 
181 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
MoFEM::Interface & mField
Definition: DirichletBC.hpp:71
Data from Cubit blocksets.
Definition: DirichletBC.hpp:43
#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.
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ iNitalize()

MoFEMErrorCode DirichletDisplacementBc::iNitalize ( )
virtual

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

Definition at line 183 of file DirichletBC.cpp.

183  {
185  if (mapZeroRows.empty() || !methodsOp.empty()) {
186 
187  std::vector<DataFromBc> bcData;
190 
191  for (auto &bc_it : bcData) {
192 
194  bc_it.scaled_values);
195 
196  bc_it.theta = bc_it.scaled_values[0]; // for rotation only
197  auto apply_rotation = [&](auto &dof) {
199  if (bc_it.is_rotation) {
200  double coords[3];
201  auto ent = dof->getEnt();
202  CHKERR mField.get_moab().get_coords(&ent, 1, coords);
203  bc_it.scaled_values = get_displacement(coords, bc_it.t_centr,
204  bc_it.t_normal, bc_it.theta);
205  }
207  };
208 
209  for (int dim = 0; dim < 3; dim++) {
210 
211  auto for_each_dof = [&](auto &dof) {
213  if (dof->getEntType() == MBVERTEX) {
214  CHKERR apply_rotation(dof);
215  if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
216  // set boundary values to field data
217  mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[0];
218  dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
219  }
220  if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
221  mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[1];
222  dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
223  }
224  if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
225  mapZeroRows[dof->getPetscGlobalDofIdx()] = bc_it.scaled_values[2];
226  dof->getFieldData() = mapZeroRows[dof->getPetscGlobalDofIdx()];
227  }
228 
229  } else {
230  if (dof->getDofCoeffIdx() == 0 && bc_it.bc_flags[0]) {
231  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
232  dof->getFieldData() = 0;
233  }
234  if (dof->getDofCoeffIdx() == 1 && bc_it.bc_flags[1]) {
235  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
236  dof->getFieldData() = 0;
237  }
238  if (dof->getDofCoeffIdx() == 2 && bc_it.bc_flags[2]) {
239  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
240  dof->getFieldData() = 0;
241  }
242  }
244  };
245 
248  bc_it.bc_ents[dim], for_each_dof);
249  }
250  }
251  dofsIndices.resize(mapZeroRows.size());
252  dofsValues.resize(mapZeroRows.size());
253  int ii = 0;
254  auto mit = mapZeroRows.begin();
255  for (; mit != mapZeroRows.end(); mit++, ii++) {
256  dofsIndices[ii] = mit->first;
257  dofsValues[ii] = mit->second;
258  }
259  }
261 }
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:72
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
MoFEM::Interface & mField
Definition: DirichletBC.hpp:71
const Problem * problemPtr
raw pointer to problem
const int dim
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:73
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:83
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
std::vector< double > dofsValues
Definition: DirichletBC.hpp:84
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
boost::ptr_vector< MethodForForceScaling > methodsOp
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
static MoFEMErrorCode applyScale(const FEMethod *fe, boost::ptr_vector< MethodForForceScaling > &methods_op, VectorDouble &nf)
auto get_displacement(double *coords, FTensor1 t_centr, FTensor1 t_normal, double theta)
MoFEMErrorCode getRotationBcFromBlock(std::vector< DataFromBc > &bc_data)
Get the Rotation Bc From Block object Use ROTATION blockset name with 7 atributes: 1 - is the rotatio...
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:82

◆ 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.

Definition at line 280 of file DirichletBC.cpp.

280  {
282 
283  if (snes_ctx == CTX_SNESNONE) {
284  if (snes_B) {
285  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
286  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
287  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
288  dofsIndices.empty() ? PETSC_NULL
289  : &dofsIndices[0],
290  dIag, PETSC_NULL, PETSC_NULL);
291  }
292  if (snes_f) {
293  CHKERR VecAssemblyBegin(snes_f);
294  CHKERR VecAssemblyEnd(snes_f);
295  for (std::vector<int>::iterator vit = dofsIndices.begin();
296  vit != dofsIndices.end(); vit++) {
297  CHKERR VecSetValue(snes_f, *vit, 0, INSERT_VALUES);
298  }
299  CHKERR VecAssemblyBegin(snes_f);
300  CHKERR VecAssemblyEnd(snes_f);
301  }
302  }
303 
304  switch (snes_ctx) {
305  case CTX_SNESNONE:
306  break;
307  case CTX_SNESSETFUNCTION: {
308  if (!dofsIndices.empty()) {
309  dofsXValues.resize(dofsIndices.size());
310  const double *a_snes_x;
311  CHKERR VecGetArrayRead(snes_x, &a_snes_x);
312  auto &dofs_by_glob_idx =
314  int idx = 0;
315  for (auto git : dofsIndices) {
316  auto dof_it = dofs_by_glob_idx.find(git);
317  if (dof_it != dofs_by_glob_idx.end()) {
318  dofsXValues[idx] = a_snes_x[dof_it->get()->getPetscLocalDofIdx()];
319  ++idx;
320  } else
321  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
322  "Dof with global %d id not found", git);
323  }
324  CHKERR VecRestoreArrayRead(snes_x, &a_snes_x);
325  }
326  CHKERR VecAssemblyBegin(snes_f);
327  CHKERR VecAssemblyEnd(snes_f);
328 
329  if (!dofsIndices.empty()) {
330  int ii = 0;
331  for (std::vector<int>::iterator vit = dofsIndices.begin();
332  vit != dofsIndices.end(); vit++, ii++) {
333  double val = 0;
334  if (!dofsXValues.empty()) {
335  val += dofsXValues[ii];
336  val += -mapZeroRows[*vit]; // in snes it is on the left hand side,
337  // that way -1
338  dofsXValues[ii] = val;
339  }
340  }
342  snes_f, dofsIndices.size(),
343  dofsIndices.empty() ? PETSC_NULL : &*dofsIndices.begin(),
344  dofsXValues.empty() ? PETSC_NULL : &*dofsXValues.begin(),
345  INSERT_VALUES);
346  }
347  CHKERR VecAssemblyBegin(snes_f);
348  CHKERR VecAssemblyEnd(snes_f);
349  } break;
350  case CTX_SNESSETJACOBIAN: {
351 
352  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
353  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
354  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
355  dofsIndices.empty() ? PETSC_NULL
356  : &*dofsIndices.begin(),
357  dIag, PETSC_NULL, PETSC_NULL);
358 
359  } break;
360  default:
361  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown snes stage");
362  }
363 
365 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
const Problem * problemPtr
raw pointer to problem
Vec & snes_f
residual
Mat & snes_B
preconditioner of jacobian matrix
std::vector< double > dofsXValues
Definition: DirichletBC.hpp:85
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:83
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:73
Vec & snes_x
state vector
#define CHKERR
Inline error check.
Definition: definitions.h:604
SNESContext snes_ctx
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
boost::shared_ptr< NumeredDofEntity_multiIndex > & getNumeredColDofsPtr() const
get access to numeredColDofsPtr storing DOFs on cols
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:82

◆ 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.

Definition at line 263 of file DirichletBC.cpp.

263  {
265 
266  CHKERR iNitalize();
267 
268  if (snes_ctx == CTX_SNESNONE) {
269  if (!dofsIndices.empty()) {
270  CHKERR VecSetValues(snes_x, dofsIndices.size(), &*dofsIndices.begin(),
271  &*dofsValues.begin(), INSERT_VALUES);
272  }
273  CHKERR VecAssemblyBegin(snes_x);
274  CHKERR VecAssemblyEnd(snes_x);
275  }
276 
278 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
virtual MoFEMErrorCode iNitalize()
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:83
std::vector< double > dofsValues
Definition: DirichletBC.hpp:84
Vec & snes_x
state vector
#define CHKERR
Inline error check.
Definition: definitions.h:604
SNESContext snes_ctx
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

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

Definition at line 113 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: