v0.14.0
LoopMethods.hpp
Go to the documentation of this file.
1 /** \file LoopMethods.hpp
2  * \brief MoFEM interface
3  *
4  * Data structures for making loops over finite elements and entities in the
5  * problem or MoFEM database.
6  *
7  */
8 
9 #ifndef __LOOPMETHODS_HPP__
10 #define __LOOPMETHODS_HPP__
11 
12 namespace MoFEM {
13 struct PetscData : public UnknownInterface {
14 
15  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
16  UnknownInterface **iface) const;
17 
18  PetscData();
19 
20  virtual ~PetscData() = default;
21 
22  enum DataContext {
24  CTX_SET_F = 1 << 0,
25  CTX_SET_A = 1 << 1,
26  CTX_SET_B = 1 << 2,
27  CTX_SET_X = 1 << 3,
28  CTX_SET_X_T = 1 << 4,
29  CTX_SET_X_TT = 1 << 6,
30  CTX_SET_TIME = 1 << 7
31  };
32 
33  using Switches = std::bitset<8>;
34 
43 
45 
46  MoFEMErrorCode copyPetscData(const PetscData &petsc_data);
47 
48  Vec f;
49  Mat A;
50  Mat B;
51  Vec x;
54 };
55 
56 /**
57  * \brief data structure for ksp (linear solver) context
58  * \ingroup mofem_loops
59  *
60  * Struture stores context data which are set in functions run by PETSc SNES
61  * functions.
62  *
63  */
64 struct KspMethod : virtual public PetscData {
65 
66  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
67  UnknownInterface **iface) const;
68 
69  /**
70  * \brief pass information about context of KSP/DM for with finite element is
71  * computed
72  */
74 
75  KspMethod();
76 
77  virtual ~KspMethod() = default;
78 
79  /**
80  * \brief copy data form another method
81  * @param ksp ksp method
82  * @return error code
83  */
85 
86  KSPContext ksp_ctx; ///< Context
87  KSP ksp; ///< KSP solver
88 
90  Mat &ksp_A;
91  Mat &ksp_B;
92 };
93 
94 /**
95  * \brief data structure for snes (nonlinear solver) context
96  * \ingroup mofem_loops
97  *
98  * Structure stores context data which are set in functions run by PETSc SNES
99  * functions.
100  *
101  */
102 struct SnesMethod : virtual protected PetscData {
103 
104  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
105  UnknownInterface **iface) const;
106 
108 
109  SnesMethod();
110 
111  virtual ~SnesMethod() = default;
112 
113  /**
114  * \brief Copy snes data
115  */
117 
119 
120  SNES snes; ///< snes solver
121  Vec &snes_x; ///< state vector
122  Vec &snes_f; ///< residual
123  Mat &snes_A; ///< jacobian matrix
124  Mat &snes_B; ///< preconditioner of jacobian matrix
125 };
126 
127 /**
128  * \brief data structure for TS (time stepping) context
129  * \ingroup mofem_loops
130  *
131  * Structure stores context data which are set in functions run by PETSc Time
132  * Stepping functions.
133  */
134 struct TSMethod : virtual protected PetscData {
135 
136  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
137  UnknownInterface **iface) const;
138 
139  enum TSContext {
146  };
147 
148  TSMethod();
149 
150  virtual ~TSMethod() = default;
151 
152  /// \brief Copy TS solver data
154 
155  TS ts; ///< time solver
156 
158 
159  PetscInt ts_step; ///< time step number
160  PetscReal ts_a; ///< shift for U_t (see PETSc Time Solver)
161  PetscReal ts_aa; ///< shift for U_tt shift for U_tt
162  PetscReal ts_t; ///< time
163  PetscReal ts_dt; ///< time step size
164 
165  Vec &ts_u; ///< state vector
166  Vec &ts_u_t; ///< time derivative of state vector
167  Vec &ts_u_tt; ///< second time derivative of state vector
168  Vec &ts_F; ///< residual vector
169 
170  Mat &ts_A; ///< Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU +
171  ///< v*dF/dU_t + a*dF/dU_tt
172  Mat &ts_B; ///< Preconditioner for ts_A
173 };
174 
175 /**
176  * \brief Data structure to exchange data between mofem and User Loop Methods.
177  * \ingroup mofem_loops
178  *
179  * It allows to exchange data between MoFEM and user functions. It stores
180  * information about multi-indices.
181  *
182  */
184 
185  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
186  UnknownInterface **iface) const {
188  *iface = const_cast<BasicMethod *>(this);
190  }
191 
192  BasicMethod();
193  virtual ~BasicMethod() = default;
194 
195  /**
196  * @brief number currently of processed method
197  */
199 
200  /**
201  * @brief local number oe methods to process
202  */
203  int loopSize;
204 
205  /** \brief get number of evaluated element in the loop
206  */
207  inline int getNinTheLoop() const { return nInTheLoop; }
208 
209  /** \brief get loop size
210  */
211  inline int getLoopSize() const { return loopSize; }
212 
213  /**
214  * @brief Llo and hi processor rank of iterated entities
215  *
216  */
217  std::pair<int, int> loHiFERank;
218 
219  /**
220  * @brief Get lo and hi processor rank of iterated entities
221  *
222  * @return raturn std::pair<int, int> loHiFERank
223  */
224  inline auto getLoHiFERank() const { return loHiFERank; }
225 
226  /**
227  * @brief Get upper rank in loop for iterating elements
228  *
229  * @return loHiFERank.first
230  */
231  inline auto getLoFERank() const { return loHiFERank.first; }
232 
233  /**
234  * @brief Get upper rank in loop for iterating elements
235  *
236  * @return loHiFERank.first
237  */
238  inline auto getHiFERank() const { return loHiFERank.second; }
239 
240  int rAnk; ///< processor rank
241 
242  int sIze; ///< number of processors in communicator
243 
245  *refinedEntitiesPtr; ///< container of mofem dof entities
246 
248  *refinedFiniteElementsPtr; ///< container of mofem finite element entities
249 
250  const Problem *problemPtr; ///< raw pointer to problem
251 
252  const Field_multiIndex *fieldsPtr; ///< raw pointer to fields container
253 
255  *entitiesPtr; ///< raw pointer to container of field entities
256 
257  const DofEntity_multiIndex *dofsPtr; ///< raw pointer container of dofs
258 
260  *finiteElementsPtr; ///< raw pointer to container finite elements
261 
263  *finiteElementsEntitiesPtr; ///< raw pointer to container finite elements
264  ///< entities
265 
267  *adjacenciesPtr; ///< raw pointer to container to adjacencies between dofs
268  ///< and finite elements
269 
270  inline unsigned int getFieldBitNumber(std::string field_name) const {
271  if (fieldsPtr) {
272  auto field_it = fieldsPtr->get<FieldName_mi_tag>().find(field_name);
273  if (field_it != fieldsPtr->get<FieldName_mi_tag>().end())
274  return (*field_it)->getBitNumber();
275  else
276  return BITFEID_SIZE;
277  } else {
278  THROW_MESSAGE("Pointer to fields multi-index is not set");
279  return BITFEID_SIZE;
280  }
281  }
282 
283  /**
284  * @brief Copy data from other base method to this base method
285  *
286  * @param basic
287  * @return MoFEMErrorCode
288  */
290 
291  /**
292  * @brief Hook function for pre-processing
293  */
294  boost::function<MoFEMErrorCode()> preProcessHook;
295 
296  /**
297  * @brief Hook function for post-processing
298  */
299  boost::function<MoFEMErrorCode()> postProcessHook;
300 
301  /**
302  * @brief Hook function for operator
303  */
304  boost::function<MoFEMErrorCode()> operatorHook;
305 
306  /** \brief function is run at the beginning of loop
307  *
308  * It is used to zeroing matrices and vectors, calculation of shape
309  * functions on reference element, preprocessing boundary conditions, etc.
310  */
311  virtual MoFEMErrorCode preProcess();
312 
313  /** \brief function is run for every finite element
314  *
315  * It is used to calculate element local matrices and assembly. It can be
316  * used for post-processing.
317  */
318  virtual MoFEMErrorCode operator()();
319 
320  /** \brief function is run at the end of loop
321  *
322  * It is used to assembly matrices and vectors, calculating global variables,
323  * f.e. total internal energy, ect.
324  *
325  * Iterating over dofs:
326  * Example1 iterating over dofs in row by name of the field
327  * for(_IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP_(this,"DISPLACEMENT",it)) { ... }
328  *
329  *
330  */
331  virtual MoFEMErrorCode postProcess();
332 
333  /**
334  * @brief Get the cache weak ptr object
335  *
336  * \note This store problem information on entities about DOFs. Each problem
337  * store different information. If you iterate over finite elements in
338  * preprocessor of TS solve element, us TS cache in the loop. Otherwise you
339  * will create undetermined behaviour or segmentation error. This is necessary
340  * compromise over bug resilience for memory saving and performance.
341  *
342  * @return boost::weak_ptr<CacheTuple>
343  */
344  inline boost::weak_ptr<CacheTuple> getCacheWeakPtr() const {
345  return cacheWeakPtr;
346  }
347 
348  boost::movelib::unique_ptr<bool> vecAssembleSwitch;
349  boost::movelib::unique_ptr<bool> matAssembleSwitch;
350 
351  boost::weak_ptr<CacheTuple> cacheWeakPtr; // cache pointer entity data
352 };
353 
354 /**
355  * \brief structure for User Loop Methods on finite elements
356  * \ingroup mofem_loops
357  *
358  * It can be used to calculate stiffness matrices, residuals, load vectors etc.
359  * It is low level class however in some class users looking for speed and
360  * efficiency, can use it directly.
361  *
362  * This class is used with Interface::loop_finite_elements, where
363  * user overloaded operator FEMethod::operator() is executed for each element in
364  * the problem. Class have to additional methods which are overloaded by user,
365  * FEMethod::preProcess() and FEMethod::postProcess() executed at beginning and
366  * end of the loop over problem elements, respectively.
367  *
368  */
369 struct FEMethod : public BasicMethod {
370 
371  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
372  UnknownInterface **iface) const {
374  *iface = const_cast<FEMethod *>(this);
376  }
377 
378  FEMethod() = default;
379 
380  std::string feName; ///< Name of finite element
381 
382  boost::shared_ptr<const NumeredEntFiniteElement>
383  numeredEntFiniteElementPtr; ///< Pointer to finite element database
384  ///< structure
385 
386  /**
387  * @brief get finite element name
388  *
389  * @return std::string
390  */
391  inline auto getFEName() const { return feName; }
392 
393  /**
394  * @brief Tet if element to skip element
395  *
396  * If is set and return false elemnent us skiped in
397  * MoFEM::Core::loop_finite_elements
398  *
399  * \note That functionality is used to run elements on particular bit levels
400  *
401  */
402  boost::function<bool(FEMethod *fe_method_ptr)> exeTestHook;
403 
404  inline auto getDataDofsPtr() const {
405  return numeredEntFiniteElementPtr->getDataDofsPtr();
406  };
407 
408  inline auto getDataVectorDofsPtr() const {
409  return numeredEntFiniteElementPtr->getDataVectorDofsPtr();
410  };
411 
413  return numeredEntFiniteElementPtr->getDataFieldEnts();
414  }
415 
416  inline boost::shared_ptr<FieldEntity_vector_view> &
418  return const_cast<NumeredEntFiniteElement *>(
421  }
422 
423  inline auto &getRowFieldEnts() const {
424  return numeredEntFiniteElementPtr->getRowFieldEnts();
425  };
426 
427  inline auto &getRowFieldEntsPtr() const {
428  return numeredEntFiniteElementPtr->getRowFieldEntsPtr();
429  };
430 
431  inline auto &getColFieldEnts() const {
432  return numeredEntFiniteElementPtr->getColFieldEnts();
433  };
434 
435  inline auto &getColFieldEntsPtr() const {
436  return numeredEntFiniteElementPtr->getColFieldEntsPtr();
437  };
438 
439  inline auto getRowDofsPtr() const {
440  return numeredEntFiniteElementPtr->getRowDofsPtr();
441  };
442 
443  inline auto getColDofsPtr() const {
444  return numeredEntFiniteElementPtr->getColDofsPtr();
445  };
446 
447  inline auto getNumberOfNodes() const;
448 
449  inline EntityHandle getFEEntityHandle() const;
450 
451  MoFEMErrorCode getNodeData(const std::string field_name, VectorDouble &data,
452  const bool reset_dofs = true);
453 
454 };
455 
456 inline auto FEMethod::getNumberOfNodes() const {
457  return moab::CN::VerticesPerEntity(numeredEntFiniteElementPtr->getEntType());
458 };
459 
461  return numeredEntFiniteElementPtr->getEnt();
462 }
463 
464 /**
465  * \brief Data structure to exchange data between mofem and User Loop Methods on
466  * entities. \ingroup mofem_loops
467  *
468  * It allows to exchange data between MoFEM and user functions. It stores
469  * information about multi-indices.
470  */
471 struct EntityMethod : public BasicMethod {
472 
473  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
474  UnknownInterface **iface) const {
476  *iface = const_cast<EntityMethod *>(this);
478  }
479 
480  EntityMethod() = default;
481 
482  boost::shared_ptr<Field> fieldPtr;
483  boost::shared_ptr<FieldEntity> entPtr;
484 };
485 
486 /**
487  * \brief Data structure to exchange data between mofem and User Loop Methods on
488  * entities. \ingroup mofem_loops
489  *
490  * It allows to exchange data between MoFEM and user functions. It stores
491  * information about multi-indices.
492  */
493 struct DofMethod : public BasicMethod {
494 
495  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
496  UnknownInterface **iface) const {
498  *iface = const_cast<DofMethod *>(this);
500  }
501 
502  DofMethod() = default;
503 
504  boost::shared_ptr<Field> fieldPtr;
505  boost::shared_ptr<DofEntity> dofPtr;
506  boost::shared_ptr<NumeredDofEntity> dofNumeredPtr;
507 };
508 
509 /// \deprecated name changed use DofMethod insead EntMethod
511 
512 } // namespace MoFEM
513 
514 #endif // __LOOPMETHODS_HPP__
515 
516 /**
517  * \defgroup mofem_loops Loops
518  * \ingroup mofem
519  */
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::SnesMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.cpp:65
FieldEntityEntFiniteElementAdjacencyMap_multiIndex
multi_index_container< FieldEntityEntFiniteElementAdjacencyMap, indexed_by< ordered_unique< tag< Composite_Unique_mi_tag >, composite_key< FieldEntityEntFiniteElementAdjacencyMap, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getEntUniqueId >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getFeUniqueId > > >, ordered_non_unique< tag< Unique_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getEntUniqueId > >, ordered_non_unique< tag< FE_Unique_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, UId, &FieldEntityEntFiniteElementAdjacencyMap::getFeUniqueId > >, ordered_non_unique< tag< FEEnt_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, EntityHandle, &FieldEntityEntFiniteElementAdjacencyMap::getFeHandle > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FieldEntityEntFiniteElementAdjacencyMap, EntityHandle, &FieldEntityEntFiniteElementAdjacencyMap::getEntHandle > > > > FieldEntityEntFiniteElementAdjacencyMap_multiIndex
MultiIndex container keeps Adjacencies Element and dof entities adjacencies and vice versa.
Definition: AdjacencyMultiIndices.hpp:114
MoFEM::FEMethod::getRowFieldEnts
auto & getRowFieldEnts() const
Definition: LoopMethods.hpp:423
MoFEM::PetscData::Switches
std::bitset< 8 > Switches
Definition: LoopMethods.hpp:33
MoFEM::SnesMethod::SNESContext
SNESContext
Definition: LoopMethods.hpp:107
MoFEM::BasicMethod::getNinTheLoop
int getNinTheLoop() const
get number of evaluated element in the loop
Definition: LoopMethods.hpp:207
MoFEM::TSMethod::TSMethod
TSMethod()
Definition: LoopMethods.cpp:91
MoFEM::KspMethod::CTX_KSPNONE
@ CTX_KSPNONE
Definition: LoopMethods.hpp:73
MoFEM::BasicMethod::refinedEntitiesPtr
const RefEntity_multiIndex * refinedEntitiesPtr
container of mofem dof entities
Definition: LoopMethods.hpp:245
DEPRECATED
#define DEPRECATED
Definition: definitions.h:17
MoFEM::TSMethod::ts_a
PetscReal ts_a
shift for U_t (see PETSc Time Solver)
Definition: LoopMethods.hpp:160
MoFEM::PetscData::x
Vec x
Definition: LoopMethods.hpp:51
MoFEM::KspMethod::ksp_B
Mat & ksp_B
Definition: LoopMethods.hpp:91
MoFEM::FEMethod::getColFieldEntsPtr
auto & getColFieldEntsPtr() const
Definition: LoopMethods.hpp:435
MoFEM::BasicMethod::nInTheLoop
int nInTheLoop
number currently of processed method
Definition: LoopMethods.hpp:198
BITFEID_SIZE
#define BITFEID_SIZE
max number of finite elements
Definition: definitions.h:234
EntityHandle
MoFEM::FEMethod::getColDofsPtr
auto getColDofsPtr() const
Definition: LoopMethods.hpp:443
MoFEM::BasicMethod::finiteElementsEntitiesPtr
const EntFiniteElement_multiIndex * finiteElementsEntitiesPtr
Definition: LoopMethods.hpp:263
MoFEM::PetscData::CtxSetB
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38
MoFEM::EntityMethod
Data structure to exchange data between mofem and User Loop Methods on entities.
Definition: LoopMethods.hpp:471
MoFEM::TSMethod::ts_F
Vec & ts_F
residual vector
Definition: LoopMethods.hpp:168
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
MoFEM::EntityMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.hpp:473
MoFEM::FEMethod::FEMethod
FEMethod()=default
MoFEM::FEMethod::getDataDofsPtr
auto getDataDofsPtr() const
Definition: LoopMethods.hpp:404
MoFEM::PetscData::CtxSetA
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
MoFEM::FEMethod::feName
std::string feName
Name of finite element.
Definition: LoopMethods.hpp:380
MoFEM::PetscData::copyPetscData
MoFEMErrorCode copyPetscData(const PetscData &petsc_data)
Definition: LoopMethods.cpp:31
MoFEM::SnesMethod::snes_x
Vec & snes_x
state vector
Definition: LoopMethods.hpp:121
MoFEM::SnesMethod::snes_ctx
SNESContext snes_ctx
Definition: LoopMethods.hpp:118
MoFEM::RefElement_multiIndex
multi_index_container< boost::shared_ptr< RefElement >, indexed_by< ordered_unique< tag< Ent_mi_tag >, const_mem_fun< RefElement::interface_type_RefEntity, EntityHandle, &RefElement::getEnt > > > > RefElement_multiIndex
Definition: RefElementMultiIndices.hpp:180
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::FieldEntity_vector_view
std::vector< boost::weak_ptr< FieldEntity > > FieldEntity_vector_view
Definition: FieldEntsMultiIndices.hpp:478
MoFEM::SnesMethod::snes_B
Mat & snes_B
preconditioner of jacobian matrix
Definition: LoopMethods.hpp:124
MoFEM::FEMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.hpp:371
MoFEM::SnesMethod
data structure for snes (nonlinear solver) context
Definition: LoopMethods.hpp:102
MoFEM::FEMethod::getDataFieldEnts
const FieldEntity_vector_view & getDataFieldEnts() const
Definition: LoopMethods.hpp:412
MoFEM::PetscData::x_tt
Vec x_tt
Definition: LoopMethods.hpp:53
MoFEM::BasicMethod::BasicMethod
BasicMethod()
Definition: LoopMethods.cpp:109
MoFEM::PetscData::B
Mat B
Definition: LoopMethods.hpp:50
MoFEM::BasicMethod::getHiFERank
auto getHiFERank() const
Get upper rank in loop for iterating elements.
Definition: LoopMethods.hpp:238
MoFEM::FieldName_mi_tag
MultiIndex Tag for field name.
Definition: TagMultiIndices.hpp:67
MoFEM::BasicMethod::preProcessHook
boost::function< MoFEMErrorCode()> preProcessHook
Hook function for pre-processing.
Definition: LoopMethods.hpp:294
MoFEM::DofMethod::DofMethod
DofMethod()=default
MoFEM::TSMethod::CTX_TSTSMONITORSET
@ CTX_TSTSMONITORSET
Definition: LoopMethods.hpp:144
MoFEM::PetscData
Definition: LoopMethods.hpp:13
MoFEM::TSMethod::CTX_TSSETRHSJACOBIAN
@ CTX_TSSETRHSJACOBIAN
Definition: LoopMethods.hpp:141
MoFEM::DofMethod
Data structure to exchange data between mofem and User Loop Methods on entities.
Definition: LoopMethods.hpp:493
MoFEM::SnesMethod::snes_A
Mat & snes_A
jacobian matrix
Definition: LoopMethods.hpp:123
MoFEM::FEMethod::getNodeData
MoFEMErrorCode getNodeData(const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
Definition: LoopMethods.cpp:172
MoFEM::RefEntity_multiIndex
multi_index_container< boost::shared_ptr< RefEntity >, indexed_by< ordered_unique< tag< Ent_mi_tag >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getEnt > >, ordered_non_unique< tag< Ent_Ent_mi_tag >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getParentEnt > >, ordered_non_unique< tag< Composite_EntType_and_ParentEntType_mi_tag >, composite_key< RefEntity, const_mem_fun< RefEntity, EntityType, &RefEntity::getEntType >, const_mem_fun< RefEntity, EntityType, &RefEntity::getParentEntType > > >, ordered_non_unique< tag< Composite_ParentEnt_And_EntType_mi_tag >, composite_key< RefEntity, const_mem_fun< RefEntity, EntityType, &RefEntity::getEntType >, const_mem_fun< RefEntity, EntityHandle, &RefEntity::getParentEnt > > > > > RefEntity_multiIndex
Definition: RefEntsMultiIndices.hpp:760
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:574
MoFEM::KspMethod::CTX_SETFUNCTION
@ CTX_SETFUNCTION
Definition: LoopMethods.hpp:73
MoFEM::SnesMethod::CTX_SNESSETJACOBIAN
@ CTX_SNESSETJACOBIAN
Definition: LoopMethods.hpp:107
MoFEM::TSMethod::ts_A
Mat & ts_A
Definition: LoopMethods.hpp:170
MoFEM::DofMethod::dofPtr
boost::shared_ptr< DofEntity > dofPtr
Definition: LoopMethods.hpp:505
MoFEM::EntMethod
DEPRECATED typedef DofMethod EntMethod
Definition: LoopMethods.hpp:510
MoFEM::SnesMethod::CTX_SNESSETFUNCTION
@ CTX_SNESSETFUNCTION
Definition: LoopMethods.hpp:107
MoFEM::PetscData::DataContext
DataContext
Definition: LoopMethods.hpp:22
MoFEM::KspMethod::~KspMethod
virtual ~KspMethod()=default
MoFEM::PetscData::CTX_SET_X_TT
@ CTX_SET_X_TT
Definition: LoopMethods.hpp:29
MoFEM::BasicMethod::fieldsPtr
const Field_multiIndex * fieldsPtr
raw pointer to fields container
Definition: LoopMethods.hpp:252
MoFEM::TSMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.cpp:85
MoFEM::BasicMethod::adjacenciesPtr
const FieldEntityEntFiniteElementAdjacencyMap_multiIndex * adjacenciesPtr
Definition: LoopMethods.hpp:267
bool
MoFEM::DofMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.hpp:495
MoFEM::SnesMethod::snes
SNES snes
snes solver
Definition: LoopMethods.hpp:120
MoFEM::BasicMethod::sIze
int sIze
number of processors in communicator
Definition: LoopMethods.hpp:242
MoFEM::FEMethod::getRowDofsPtr
auto getRowDofsPtr() const
Definition: LoopMethods.hpp:439
MoFEM::BasicMethod::entitiesPtr
const FieldEntity_multiIndex * entitiesPtr
raw pointer to container of field entities
Definition: LoopMethods.hpp:255
MoFEM::BasicMethod::getLoopSize
int getLoopSize() const
get loop size
Definition: LoopMethods.hpp:211
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::EntityMethod::EntityMethod
EntityMethod()=default
MoFEM::BasicMethod::getLoFERank
auto getLoFERank() const
Get upper rank in loop for iterating elements.
Definition: LoopMethods.hpp:231
MoFEM::KspMethod
data structure for ksp (linear solver) context
Definition: LoopMethods.hpp:64
MoFEM::PetscData::PetscData
PetscData()
Definition: LoopMethods.cpp:27
MoFEM::PetscData::data_ctx
Switches data_ctx
Definition: LoopMethods.hpp:44
MoFEM::EntityMethod::entPtr
boost::shared_ptr< FieldEntity > entPtr
Definition: LoopMethods.hpp:483
MoFEM::SnesMethod::CTX_SNESNONE
@ CTX_SNESNONE
Definition: LoopMethods.hpp:107
MoFEM::FEMethod::getNumberOfNodes
auto getNumberOfNodes() const
Definition: LoopMethods.hpp:456
MoFEM::TSMethod::ts
TS ts
time solver
Definition: LoopMethods.hpp:155
MoFEM::KspMethod::KSPContext
KSPContext
pass information about context of KSP/DM for with finite element is computed
Definition: LoopMethods.hpp:73
MoFEM::TSMethod::CTX_TSSETIJACOBIAN
@ CTX_TSSETIJACOBIAN
Definition: LoopMethods.hpp:143
MoFEM::BasicMethod::rAnk
int rAnk
processor rank
Definition: LoopMethods.hpp:240
MoFEM::EntityMethod::fieldPtr
boost::shared_ptr< Field > fieldPtr
Definition: LoopMethods.hpp:482
MoFEM::BasicMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.hpp:185
MoFEM::SnesMethod::copySnes
MoFEMErrorCode copySnes(const SnesMethod &snes)
Copy snes data.
Definition: LoopMethods.cpp:75
MoFEM::PetscData::A
Mat A
Definition: LoopMethods.hpp:49
MoFEM::TSMethod::ts_u_t
Vec & ts_u_t
time derivative of state vector
Definition: LoopMethods.hpp:166
MoFEM::PetscData::CtxSetX_T
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:40
MoFEM::BasicMethod::preProcess
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
Definition: LoopMethods.cpp:138
MoFEM::TSMethod::CTX_TSNONE
@ CTX_TSNONE
Definition: LoopMethods.hpp:145
MoFEM::BasicMethod::postProcessHook
boost::function< MoFEMErrorCode()> postProcessHook
Hook function for post-processing.
Definition: LoopMethods.hpp:299
MoFEM::KspMethod::KspMethod
KspMethod()
Definition: LoopMethods.cpp:51
MoFEM::BasicMethod::refinedFiniteElementsPtr
const RefElement_multiIndex * refinedFiniteElementsPtr
container of mofem finite element entities
Definition: LoopMethods.hpp:248
MoFEM::TSMethod::ts_ctx
TSContext ts_ctx
Definition: LoopMethods.hpp:157
MoFEM::FEMethod::getFEName
auto getFEName() const
get finite element name
Definition: LoopMethods.hpp:391
MoFEM::PetscData::CtxSetTime
static constexpr Switches CtxSetTime
Definition: LoopMethods.hpp:42
MoFEM::PetscData::CtxSetNone
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:35
MoFEM::BasicMethod
Data structure to exchange data between mofem and User Loop Methods.
Definition: LoopMethods.hpp:183
MoFEM::BasicMethod::loopSize
int loopSize
local number oe methods to process
Definition: LoopMethods.hpp:203
MoFEM::KspMethod::ksp
KSP ksp
KSP solver.
Definition: LoopMethods.hpp:87
MoFEM::FieldEntity_multiIndex
multi_index_container< boost::shared_ptr< FieldEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, member< FieldEntity, UId, &FieldEntity::localUId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< FieldEntity::interface_type_RefEntity, EntityHandle, &FieldEntity::getEnt > > > > FieldEntity_multiIndex
Definition: FieldEntsMultiIndices.hpp:489
MoFEM::FEMethod::getDataVectorDofsPtr
auto getDataVectorDofsPtr() const
Definition: LoopMethods.hpp:408
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::PetscData::f
Vec f
Definition: LoopMethods.hpp:48
MoFEM::BasicMethod::vecAssembleSwitch
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: LoopMethods.hpp:348
MoFEM::PetscData::CtxSetX_TT
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:41
MoFEM::BasicMethod::copyBasicMethod
MoFEMErrorCode copyBasicMethod(const BasicMethod &basic)
Copy data from other base method to this base method.
Definition: LoopMethods.cpp:117
MoFEM::TSMethod::~TSMethod
virtual ~TSMethod()=default
MoFEM::PetscData::CtxSetX
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:39
MoFEM::TSMethod::ts_step
PetscInt ts_step
time step number
Definition: LoopMethods.hpp:159
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MoFEM::PetscData::CTX_SET_TIME
@ CTX_SET_TIME
Definition: LoopMethods.hpp:30
MoFEM::KspMethod::ksp_A
Mat & ksp_A
Definition: LoopMethods.hpp:90
MoFEM::PetscData::CTX_SET_B
@ CTX_SET_B
Definition: LoopMethods.hpp:26
MoFEM::SnesMethod::snes_f
Vec & snes_f
residual
Definition: LoopMethods.hpp:122
MoFEM::BasicMethod::finiteElementsPtr
const FiniteElement_multiIndex * finiteElementsPtr
raw pointer to container finite elements
Definition: LoopMethods.hpp:260
MoFEM::FEMethod::getRowFieldEntsPtr
auto & getRowFieldEntsPtr() const
Definition: LoopMethods.hpp:427
MoFEM::TSMethod::TSContext
TSContext
Definition: LoopMethods.hpp:139
MoFEM::BasicMethod::operatorHook
boost::function< MoFEMErrorCode()> operatorHook
Hook function for operator.
Definition: LoopMethods.hpp:304
MoFEM::TSMethod::ts_t
PetscReal ts_t
time
Definition: LoopMethods.hpp:162
MoFEM::PetscData::CTX_SET_A
@ CTX_SET_A
Definition: LoopMethods.hpp:25
MoFEM::BasicMethod::~BasicMethod
virtual ~BasicMethod()=default
MoFEM::BasicMethod::getFieldBitNumber
unsigned int getFieldBitNumber(std::string field_name) const
Definition: LoopMethods.hpp:270
MoFEM::TSMethod::ts_aa
PetscReal ts_aa
shift for U_tt shift for U_tt
Definition: LoopMethods.hpp:161
MoFEM::BasicMethod::dofsPtr
const DofEntity_multiIndex * dofsPtr
raw pointer container of dofs
Definition: LoopMethods.hpp:257
MoFEM::PetscData::CTX_SET_F
@ CTX_SET_F
Definition: LoopMethods.hpp:24
Field_multiIndex
multi_index_container< boost::shared_ptr< Field >, indexed_by< hashed_unique< tag< BitFieldId_mi_tag >, const_mem_fun< Field, const BitFieldId &, &Field::getId >, HashBit< BitFieldId >, EqBit< BitFieldId > >, ordered_unique< tag< Meshset_mi_tag >, member< Field, EntityHandle, &Field::meshSet > >, ordered_unique< tag< FieldName_mi_tag >, const_mem_fun< Field, boost::string_ref, &Field::getNameRef > >, ordered_non_unique< tag< BitFieldId_space_mi_tag >, const_mem_fun< Field, FieldSpace, &Field::getSpace > > > > Field_multiIndex
Field_multiIndex for Field.
Definition: FieldMultiIndices.hpp:508
MoFEM::SnesMethod::SnesMethod
SnesMethod()
Definition: LoopMethods.cpp:71
MoFEM::TSMethod::ts_u
Vec & ts_u
state vector
Definition: LoopMethods.hpp:165
MoFEM::TSMethod::copyTs
MoFEMErrorCode copyTs(const TSMethod &ts)
Copy TS solver data.
Definition: LoopMethods.cpp:96
MoFEM::TSMethod
data structure for TS (time stepping) context
Definition: LoopMethods.hpp:134
MoFEM::PetscData::CTX_SET_NONE
@ CTX_SET_NONE
Definition: LoopMethods.hpp:23
MoFEM::BasicMethod::cacheWeakPtr
boost::weak_ptr< CacheTuple > cacheWeakPtr
Definition: LoopMethods.hpp:351
MoFEM::KspMethod::copyKsp
MoFEMErrorCode copyKsp(const KspMethod &ksp)
copy data form another method
Definition: LoopMethods.cpp:55
MoFEM::DofMethod::dofNumeredPtr
boost::shared_ptr< NumeredDofEntity > dofNumeredPtr
Definition: LoopMethods.hpp:506
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::BasicMethod::problemPtr
const Problem * problemPtr
raw pointer to problem
Definition: LoopMethods.hpp:250
MoFEM::PetscData::x_t
Vec x_t
Definition: LoopMethods.hpp:52
MoFEM::FEMethod::exeTestHook
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
Tet if element to skip element.
Definition: LoopMethods.hpp:402
MoFEM::FEMethod::getColFieldEnts
auto & getColFieldEnts() const
Definition: LoopMethods.hpp:431
MoFEM::PetscData::CTX_SET_X_T
@ CTX_SET_X_T
Definition: LoopMethods.hpp:28
MoFEM::PetscData::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.cpp:21
MoFEM::BasicMethod::loHiFERank
std::pair< int, int > loHiFERank
Llo and hi processor rank of iterated entities.
Definition: LoopMethods.hpp:217
MoFEM::KspMethod::ksp_f
Vec & ksp_f
Definition: LoopMethods.hpp:89
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
MoFEM::TSMethod::CTX_TSSETIFUNCTION
@ CTX_TSSETIFUNCTION
Definition: LoopMethods.hpp:142
FiniteElement_multiIndex
multi_index_container< boost::shared_ptr< FiniteElement >, indexed_by< hashed_unique< tag< FiniteElement_Meshset_mi_tag >, member< FiniteElement, EntityHandle, &FiniteElement::meshset > >, hashed_unique< tag< BitFEId_mi_tag >, const_mem_fun< FiniteElement, BitFEId, &FiniteElement::getId >, HashBit< BitFEId >, EqBit< BitFEId > >, ordered_unique< tag< FiniteElement_name_mi_tag >, const_mem_fun< FiniteElement, boost::string_ref, &FiniteElement::getNameRef > > > > FiniteElement_multiIndex
MultiIndex for entities for FiniteElement.
Definition: FEMultiIndices.hpp:849
MoFEM::NumeredEntFiniteElement
Partitioned (Indexed) Finite Element in Problem.
Definition: FEMultiIndices.hpp:728
MoFEM::KspMethod::CTX_OPERATORS
@ CTX_OPERATORS
Definition: LoopMethods.hpp:73
MoFEM::FEMethod::getFEEntityHandle
EntityHandle getFEEntityHandle() const
Definition: LoopMethods.hpp:460
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEM::BasicMethod::postProcess
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
Definition: LoopMethods.cpp:149
MoFEM::PetscData::CtxSetF
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:36
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::field_it
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
Definition: Projection10NodeCoordsOnField.cpp:123
MoFEM::TSMethod::ts_dt
PetscReal ts_dt
time step size
Definition: LoopMethods.hpp:163
MoFEM::BasicMethod::getCacheWeakPtr
boost::weak_ptr< CacheTuple > getCacheWeakPtr() const
Get the cache weak ptr object.
Definition: LoopMethods.hpp:344
DofEntity_multiIndex
multi_index_container< boost::shared_ptr< DofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< DofEntity, UId, &DofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< DofEntity, EntityHandle, &DofEntity::getEnt > > > > DofEntity_multiIndex
MultiIndex container keeps DofEntity.
Definition: DofsMultiIndices.hpp:317
MoFEM::PetscData::CTX_SET_X
@ CTX_SET_X
Definition: LoopMethods.hpp:27
MoFEM::FEMethod::getDataFieldEntsPtr
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr() const
Definition: LoopMethods.hpp:417
MoFEM::SnesMethod::~SnesMethod
virtual ~SnesMethod()=default
MoFEM::TSMethod::ts_B
Mat & ts_B
Preconditioner for ts_A.
Definition: LoopMethods.hpp:172
MoFEM::KspMethod::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: LoopMethods.cpp:45
MoFEM::DofMethod::fieldPtr
boost::shared_ptr< Field > fieldPtr
Definition: LoopMethods.hpp:504
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
EntFiniteElement_multiIndex
multi_index_container< boost::shared_ptr< EntFiniteElement >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< EntFiniteElement, UId, &EntFiniteElement::getLocalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< EntFiniteElement::interface_type_RefEntity, EntityHandle, &EntFiniteElement::getEnt > > > > EntFiniteElement_multiIndex
MultiIndex container for EntFiniteElement.
Definition: FEMultiIndices.hpp:800
MoFEM::PetscData::~PetscData
virtual ~PetscData()=default
MoFEM::TSMethod::ts_u_tt
Vec & ts_u_tt
second time derivative of state vector
Definition: LoopMethods.hpp:167
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::BasicMethod::operator()
virtual MoFEMErrorCode operator()()
function is run for every finite element
Definition: LoopMethods.cpp:160
MoFEM::BasicMethod::getLoHiFERank
auto getLoHiFERank() const
Get lo and hi processor rank of iterated entities.
Definition: LoopMethods.hpp:224
MoFEM::TSMethod::CTX_TSSETRHSFUNCTION
@ CTX_TSSETRHSFUNCTION
Definition: LoopMethods.hpp:140
MoFEM::KspMethod::ksp_ctx
KSPContext ksp_ctx
Context.
Definition: LoopMethods.hpp:86
MoFEM::BasicMethod::matAssembleSwitch
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: LoopMethods.hpp:349