v0.9.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)
 
 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 ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
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
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
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 ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- 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 ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. 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< intdofsIndices
 
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
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- 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
TSContext ts_ctx
 
TS ts
 time solver 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...
 
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...
 

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::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
}
 
- 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
Vec snes_f
residual
MoFEM::Interface & mField
Definition: DirichletBC.hpp:45
Vec ts_F
residual vector
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:47
Vec snes_x
state vector
Vec ts_u
state vector
Mat snes_B
preconditioner of jacobian matrix
Mat ts_B
Preconditioner for ts_A.

◆ 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
Vec snes_f
residual
MoFEM::Interface & mField
Definition: DirichletBC.hpp:45
Vec ts_F
residual vector
double dIag
diagonal value set on zeroed column and rows
Definition: DirichletBC.hpp:47
Vec snes_x
state vector
Vec ts_u
state vector
Mat snes_B
preconditioner of jacobian matrix
Mat ts_B
Preconditioner for ts_A.

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  if (dof->getEntType() == MBVERTEX) {
105  if (dof->getDofCoeffIdx() == 0 && mydata.data.flag1) {
106  mapZeroRows[dof->getPetscGlobalDofIdx()] = scaled_values[0];
107  }
108  if (dof->getDofCoeffIdx() == 1 && mydata.data.flag2) {
109  mapZeroRows[dof->getPetscGlobalDofIdx()] = scaled_values[1];
110  }
111  if (dof->getDofCoeffIdx() == 2 && mydata.data.flag3) {
112  mapZeroRows[dof->getPetscGlobalDofIdx()] = scaled_values[2];
113  }
114  } else {
115  if (dof->getDofCoeffIdx() == 0 && mydata.data.flag1) {
116  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
117  }
118  if (dof->getDofCoeffIdx() == 1 && mydata.data.flag2) {
119  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
120  }
121  if (dof->getDofCoeffIdx() == 2 && mydata.data.flag3) {
122  mapZeroRows[dof->getPetscGlobalDofIdx()] = 0;
123  }
124  }
126  };
128  for_each_dof);
129  }
130  }
131  dofsIndices.resize(mapZeroRows.size());
132  dofsValues.resize(mapZeroRows.size());
133  int ii = 0;
134  auto mit = mapZeroRows.begin();
135  for (; mit != mapZeroRows.end(); mit++, ii++) {
136  dofsIndices[ii] = mit->first;
137  dofsValues[ii] = mit->second;
138  }
139  }
141 }
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:501
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
MoFEM::Interface & mField
Definition: DirichletBC.hpp:45
const Problem * problemPtr
raw pointer to problem
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:596
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:407
static MoFEMErrorCode applyScale(const FEMethod *fe, boost::ptr_vector< MethodForForceScaling > &methods_op, VectorDouble &nf)
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:72
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 176 of file DirichletBC.cpp.

176  {
178 
179  switch (ts_ctx) {
180  case CTX_TSSETIFUNCTION: {
182  snes_x = ts_u;
183  snes_f = ts_F;
184  break;
185  }
186  case CTX_TSSETIJACOBIAN: {
188  snes_B = ts_B;
189  break;
190  }
191  default:
192  break;
193  }
194 
195  if (snes_ctx == CTX_SNESNONE && ts_ctx == CTX_TSNONE) {
196  if (snes_B) {
197  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
198  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
199  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
200  dofsIndices.empty() ? PETSC_NULL
201  : &dofsIndices[0],
202  dIag, PETSC_NULL, PETSC_NULL);
203  }
204  if (snes_f) {
205  CHKERR VecAssemblyBegin(snes_f);
206  CHKERR VecAssemblyEnd(snes_f);
207  for (std::vector<int>::iterator vit = dofsIndices.begin();
208  vit != dofsIndices.end(); vit++) {
209  CHKERR VecSetValue(snes_f, *vit, 0, INSERT_VALUES);
210  }
211  CHKERR VecAssemblyBegin(snes_f);
212  CHKERR VecAssemblyEnd(snes_f);
213  }
214  }
215 
216  switch (snes_ctx) {
217  case CTX_SNESNONE:
218  break;
219  case CTX_SNESSETFUNCTION: {
220  if (!dofsIndices.empty()) {
221  dofsXValues.resize(dofsIndices.size());
222  const double *a_snes_x;
223  CHKERR VecGetArrayRead(snes_x, &a_snes_x);
224  auto &dofs_by_glob_idx =
226  int idx = 0;
227  for (auto git : dofsIndices) {
228  auto dof_it = dofs_by_glob_idx.find(git);
229  if (dof_it != dofs_by_glob_idx.end()) {
230  dofsXValues[idx] = a_snes_x[dof_it->get()->getPetscLocalDofIdx()];
231  ++idx;
232  } else
233  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
234  "Dof with global %d id not found", git);
235  }
236  CHKERR VecRestoreArrayRead(snes_x, &a_snes_x);
237  }
238  CHKERR VecAssemblyBegin(snes_f);
239  CHKERR VecAssemblyEnd(snes_f);
240 
241  if (!dofsIndices.empty()) {
242  int ii = 0;
243  for (std::vector<int>::iterator vit = dofsIndices.begin();
244  vit != dofsIndices.end(); vit++, ii++) {
245  double val = 0;
246  if (!dofsXValues.empty()) {
247  val += dofsXValues[ii];
248  val += -mapZeroRows[*vit]; // in snes it is on the left hand side,
249  // that way -1
250  dofsXValues[ii] = val;
251  }
252  }
254  snes_f, dofsIndices.size(),
255  dofsIndices.empty() ? PETSC_NULL : &*dofsIndices.begin(),
256  dofsXValues.empty() ? PETSC_NULL : &*dofsXValues.begin(),
257  INSERT_VALUES);
258  }
259  CHKERR VecAssemblyBegin(snes_f);
260  CHKERR VecAssemblyEnd(snes_f);
261  } break;
262  case CTX_SNESSETJACOBIAN: {
263 
264  CHKERR MatAssemblyBegin(snes_B, MAT_FINAL_ASSEMBLY);
265  CHKERR MatAssemblyEnd(snes_B, MAT_FINAL_ASSEMBLY);
266  CHKERR MatZeroRowsColumns(snes_B, dofsIndices.size(),
267  dofsIndices.empty() ? PETSC_NULL
268  : &*dofsIndices.begin(),
269  dIag, PETSC_NULL, PETSC_NULL);
270 
271  } break;
272  default:
273  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown snes stage");
274  }
275 
277 }
Vec snes_f
residual
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
const boost::shared_ptr< NumeredDofEntity_multiIndex > & getNumeredDofsCols() const
get access to numeredDofsCols storing DOFs on cols
const Problem * problemPtr
raw pointer to problem
Vec ts_F
residual vector
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:596
Vec ts_u
state vector
Mat snes_B
preconditioner of jacobian matrix
SNESContext snes_ctx
Mat ts_B
Preconditioner for ts_A.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
TSContext ts_ctx
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 143 of file DirichletBC.cpp.

143  {
145 
146  switch (ts_ctx) {
147  case CTX_TSSETIFUNCTION: {
149  snes_x = ts_u;
150  snes_f = ts_F;
151  break;
152  }
153  case CTX_TSSETIJACOBIAN: {
155  snes_B = ts_B;
156  break;
157  }
158  default:
159  break;
160  }
161 
162  CHKERR iNitalize();
163 
164  if (snes_ctx == CTX_SNESNONE && ts_ctx == CTX_TSNONE) {
165  if (!dofsIndices.empty()) {
166  CHKERR VecSetValues(snes_x, dofsIndices.size(), &*dofsIndices.begin(),
167  &*dofsValues.begin(), INSERT_VALUES);
168  }
169  CHKERR VecAssemblyBegin(snes_x);
170  CHKERR VecAssemblyEnd(snes_x);
171  }
172 
174 }
Vec snes_f
residual
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
virtual MoFEMErrorCode iNitalize()
Definition: DirichletBC.cpp:73
Vec ts_F
residual vector
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:596
Vec ts_u
state vector
Mat snes_B
preconditioner of jacobian matrix
SNESContext snes_ctx
Mat ts_B
Preconditioner for ts_A.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
TSContext ts_ctx

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: