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