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

Set Dirichlet boundary conditions on displacements. More...

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

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

Public Member Functions

 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name, Mat Aij, Vec X, Vec F)
 
 DirichletDisplacementBc (MoFEM::Interface &m_field, const std::string &field_name)
 
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...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
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...
 
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
 
boost::ptr_vector< MethodForForceScalingmethodsOp
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- 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_tt (see PETSc Time Solver) More...
 
PetscReal ts_v
 shift for U_t shift for U_t 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, elasticity.cpp, and elasticity_mixed_formulation.cpp.

Definition at line 43 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 
)

Definition at line 50 of file DirichletBC.cpp.

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

◆ DirichletDisplacementBc() [2/2]

DirichletDisplacementBc::DirichletDisplacementBc ( MoFEM::Interface m_field,
const std::string &  field_name 
)

Definition at line 62 of file DirichletBC.cpp.

64  : mField(m_field), fieldName(field_name), dIag(1) {
65  snes_B = PETSC_NULL;
66  snes_x = PETSC_NULL;
67  snes_f = PETSC_NULL;
68  ts_B = PETSC_NULL;
69  ts_u = PETSC_NULL;
70  ts_F = PETSC_NULL;
71 };
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:46
Mat & ts_B
Preconditioner for ts_A.
MoFEM::Interface & mField
Definition: DirichletBC.hpp:45
Vec & snes_f
residual
Mat & snes_B
preconditioner of jacobian matrix
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:47
Vec & snes_x
state vector
Vec & ts_F
residual vector
Vec & ts_u
state vector

Member Function Documentation

◆ iNitalize()

MoFEMErrorCode DirichletDisplacementBc::iNitalize ( )
virtual

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

Definition at line 73 of file DirichletBC.cpp.

73  {
75  if (mapZeroRows.empty() || !methodsOp.empty()) {
77  mField, NODESET | DISPLACEMENTSET, it)) {
79  CHKERR it->getBcDataStructure(mydata);
80  VectorDouble scaled_values(3);
81  scaled_values[0] = mydata.data.value1;
82  scaled_values[1] = mydata.data.value2;
83  scaled_values[2] = mydata.data.value3;
85  for (int dim = 0; dim < 3; dim++) {
86  Range ents;
87  CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), dim, ents,
88  true);
89  if (dim > 1) {
90  Range _edges;
91  CHKERR mField.get_moab().get_adjacencies(ents, 1, false, _edges,
92  moab::Interface::UNION);
93  ents.insert(_edges.begin(), _edges.end());
94  }
95  if (dim > 0) {
96  Range _nodes;
97  CHKERR mField.get_moab().get_connectivity(ents, _nodes, true);
98  ents.insert(_nodes.begin(), _nodes.end());
99  }
100 
101  auto for_each_dof = [&](auto &dof) {
103 
104 
105  if (dof->getEntType() == MBVERTEX) {
106  if (dof->getDofCoeffIdx() == 0 && mydata.data.flag1) {
107  mapZeroRows[dof->getPetscGlobalDofIdx()] = scaled_values[0];
108  }
109  if (dof->getDofCoeffIdx() == 1 && mydata.data.flag2) {
110  mapZeroRows[dof->getPetscGlobalDofIdx()] = scaled_values[1];
111  }
112  if (dof->getDofCoeffIdx() == 2 && mydata.data.flag3) {
113  mapZeroRows[dof->getPetscGlobalDofIdx()] = scaled_values[2];
114  }
115  } else {
116  if (dof->getDofCoeffIdx() == 0 && mydata.data.flag1) {
117  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
118  }
119  if (dof->getDofCoeffIdx() == 1 && mydata.data.flag2) {
120  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
121  }
122  if (dof->getDofCoeffIdx() == 2 && mydata.data.flag3) {
123  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
124  }
125  }
127  };
129  for_each_dof);
130  }
131  }
132  dofsIndices.resize(mapZeroRows.size());
133  dofsValues.resize(mapZeroRows.size());
134  int ii = 0;
135  auto mit = mapZeroRows.begin();
136  for (; mit != mapZeroRows.end(); mit++, ii++) {
137  dofsIndices[ii] = mit->first;
138  dofsValues[ii] = mit->second;
139  }
140  }
142 }
const std::string fieldName
field name to set Dirichlet BC
Definition: DirichletBC.hpp:46
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:507
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
MoFEM::Interface & mField
Definition: DirichletBC.hpp:45
const Problem * problemPtr
raw pointer to problem
const int dim
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.
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:55
std::vector< double > dofsValues
Definition: DirichletBC.hpp:56
boost::ptr_vector< MethodForForceScaling > methodsOp
Definition: DirichletBC.hpp:63
#define CHKERR
Inline error check.
Definition: definitions.h:602
static MoFEMErrorCode set_numered_dofs_on_ents(const Problem *problem_ptr, const string &field_name, Range &ents, boost::function< MoFEMErrorCode(const boost::shared_ptr< MoFEM::NumeredDofEntity > &dof)> for_each_dof)
Definition: DirichletBC.cpp:23
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
static MoFEMErrorCode applyScale(const FEMethod *fe, boost::ptr_vector< MethodForForceScaling > &methods_op, VectorDouble &nf)
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:73
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:54

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

161  {
163 
164  if (snes_ctx == CTX_SNESNONE) {
165  if (snes_B) {
166  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
167  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
168  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
169  dofsIndices.empty() ? PETSC_NULL
170  : &dofsIndices[0],
171  dIag, PETSC_NULL, PETSC_NULL);
172  }
173  if (snes_f) {
174  CHKERR VecAssemblyBegin(snes_f);
175  CHKERR VecAssemblyEnd(snes_f);
176  for (std::vector<int>::iterator vit = dofsIndices.begin();
177  vit != dofsIndices.end(); vit++) {
178  CHKERR VecSetValue(snes_f, *vit, 0, INSERT_VALUES);
179  }
180  CHKERR VecAssemblyBegin(snes_f);
181  CHKERR VecAssemblyEnd(snes_f);
182  }
183  }
184 
185  switch (snes_ctx) {
186  case CTX_SNESNONE:
187  break;
188  case CTX_SNESSETFUNCTION: {
189  if (!dofsIndices.empty()) {
190  dofsXValues.resize(dofsIndices.size());
191  const double *a_snes_x;
192  CHKERR VecGetArrayRead(snes_x, &a_snes_x);
193  auto &dofs_by_glob_idx =
195  int idx = 0;
196  for (auto git : dofsIndices) {
197  auto dof_it = dofs_by_glob_idx.find(git);
198  if (dof_it != dofs_by_glob_idx.end()) {
199  dofsXValues[idx] = a_snes_x[dof_it->get()->getPetscLocalDofIdx()];
200  ++idx;
201  } else
202  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
203  "Dof with global %d id not found", git);
204  }
205  CHKERR VecRestoreArrayRead(snes_x, &a_snes_x);
206  }
207  CHKERR VecAssemblyBegin(snes_f);
208  CHKERR VecAssemblyEnd(snes_f);
209 
210  if (!dofsIndices.empty()) {
211  int ii = 0;
212  for (std::vector<int>::iterator vit = dofsIndices.begin();
213  vit != dofsIndices.end(); vit++, ii++) {
214  double val = 0;
215  if (!dofsXValues.empty()) {
216  val += dofsXValues[ii];
217  val += -mapZeroRows[*vit]; // in snes it is on the left hand side,
218  // that way -1
219  dofsXValues[ii] = val;
220  }
221  }
223  snes_f, dofsIndices.size(),
224  dofsIndices.empty() ? PETSC_NULL : &*dofsIndices.begin(),
225  dofsXValues.empty() ? PETSC_NULL : &*dofsXValues.begin(),
226  INSERT_VALUES);
227  }
228  CHKERR VecAssemblyBegin(snes_f);
229  CHKERR VecAssemblyEnd(snes_f);
230  } break;
231  case CTX_SNESSETJACOBIAN: {
232 
233  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
234  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
235  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
236  dofsIndices.empty() ? PETSC_NULL
237  : &*dofsIndices.begin(),
238  dIag, PETSC_NULL, PETSC_NULL);
239 
240  } break;
241  default:
242  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown snes stage");
243  }
244 
246 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
const boost::shared_ptr< NumeredDofEntity_multiIndex > & getNumeredDofsCols() const
get access to numeredDofsCols storing DOFs on cols
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:57
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:55
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:47
Vec & snes_x
state vector
MoFEMErrorCode VecSetValues(Vec V, const DataForcesAndSourcesCore::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
#define CHKERR
Inline error check.
Definition: definitions.h:602
SNESContext snes_ctx
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
std::map< DofIdx, FieldData > mapZeroRows
Definition: DirichletBC.hpp:54

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

144  {
146 
147  CHKERR iNitalize();
148 
149  if (snes_ctx == CTX_SNESNONE) {
150  if (!dofsIndices.empty()) {
151  CHKERR VecSetValues(snes_x, dofsIndices.size(), &*dofsIndices.begin(),
152  &*dofsValues.begin(), INSERT_VALUES);
153  }
154  CHKERR VecAssemblyBegin(snes_x);
155  CHKERR VecAssemblyEnd(snes_x);
156  }
157 
159 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
virtual MoFEMErrorCode iNitalize()
Definition: DirichletBC.cpp:73
std::vector< int > dofsIndices
Definition: DirichletBC.hpp:55
std::vector< double > dofsValues
Definition: DirichletBC.hpp:56
Vec & snes_x
state vector
MoFEMErrorCode VecSetValues(Vec V, const DataForcesAndSourcesCore::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
#define CHKERR
Inline error check.
Definition: definitions.h:602
SNESContext snes_ctx
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413

Member Data Documentation

◆ dIag

double DirichletDisplacementBc::dIag

diagonal value set on zeroed column and rows

Definition at line 47 of file DirichletBC.hpp.

◆ dofsIndices

std::vector<int> DirichletDisplacementBc::dofsIndices

Definition at line 55 of file DirichletBC.hpp.

◆ dofsValues

std::vector<double> DirichletDisplacementBc::dofsValues

Definition at line 56 of file DirichletBC.hpp.

◆ dofsXValues

std::vector<double> DirichletDisplacementBc::dofsXValues

Definition at line 57 of file DirichletBC.hpp.

◆ fieldName

const std::string DirichletDisplacementBc::fieldName

field name to set Dirichlet BC

Definition at line 46 of file DirichletBC.hpp.

◆ mapZeroRows

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

Definition at line 54 of file DirichletBC.hpp.

◆ methodsOp

boost::ptr_vector<MethodForForceScaling> DirichletDisplacementBc::methodsOp

Definition at line 63 of file DirichletBC.hpp.

◆ mField

MoFEM::Interface& DirichletDisplacementBc::mField

Definition at line 45 of file DirichletBC.hpp.


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