v0.15.0
Loading...
Searching...
No Matches
LoopMethods.cpp
Go to the documentation of this file.
1/** \file LoopMethods.cpp
2\brief methods for managing loops
3*/
4
5
6#include <MoFEM.hpp>
7
8namespace MoFEM {
9
18
19// PetscData
21PetscData::query_interface(boost::typeindex::type_index type_index,
22 UnknownInterface **iface) const {
23 *iface = const_cast<PetscData *>(this);
24 return 0;
25}
26
28 : data_ctx(PetscData::CtxSetNone), f(PETSC_NULLPTR), A(PETSC_NULLPTR),
29 B(PETSC_NULLPTR), x(PETSC_NULLPTR), dx(PETSC_NULLPTR), x_t(PETSC_NULLPTR),
30 x_tt(PETSC_NULLPTR) {}
31
33 this->data_ctx = petsc_data.data_ctx;
34 this->f = petsc_data.f;
35 this->A = petsc_data.A;
36 this->B = petsc_data.B;
37 this->x = petsc_data.x;
38 this->dx = petsc_data.dx;
39 this->x_t = petsc_data.x_t;
40 this->x_tt = petsc_data.x_tt;
41 return 0;
42}
43
44// KSP
46KspMethod::query_interface(boost::typeindex::type_index type_index,
47 UnknownInterface **iface) const {
48 *iface = const_cast<KspMethod *>(this);
49 return 0;
50};
51
53 : ksp_ctx(CTX_KSPNONE), ksp(PETSC_NULLPTR), ksp_f(PetscData::f),
54 ksp_A(PetscData::A), ksp_B(PetscData::B) {}
55
63
64// SNES
66SnesMethod::query_interface(boost::typeindex::type_index type_index,
67 UnknownInterface **iface) const {
68 *iface = const_cast<SnesMethod *>(this);
69 return 0;
70}
71
73 : snes_ctx(CTX_SNESNONE), snes_x(PetscData::x), snes_dx(PetscData::dx),
74 snes_f(PetscData::f), snes_A(PetscData::A), snes_B(PetscData::B) {}
75
83
84// TAO
86TaoMethod::query_interface(boost::typeindex::type_index type_index,
87 UnknownInterface **iface) const {
88 *iface = const_cast<TaoMethod *>(this);
89 return 0;
90}
91
93 : tao_ctx(CTX_TAO_NONE), tao_x(PetscData::x), tao_f(PetscData::f),
94 tao_A(PetscData::A), tao_B(PetscData::B) {}
95
99 this->tao_ctx = tao.tao_ctx;
100 this->tao = tao.tao;
102}
103
104// TS
106TSMethod::query_interface(boost::typeindex::type_index type_index,
107 UnknownInterface **iface) const {
108 *iface = const_cast<TSMethod *>(this);
109 return 0;
110}
111
113 : ts_ctx(CTX_TSNONE), ts_step(-1), ts_a(0), ts_t(0), ts_u(PetscData::x),
114 ts_u_t(PetscData::x_t), ts_u_tt(PetscData::x_tt), ts_F(PetscData::f),
115 ts_A(PetscData::A), ts_B(PetscData::B) {}
116
120 this->ts_ctx = ts.ts_ctx;
121 this->ts = ts.ts;
122 this->ts_step = ts.ts_step;
123 this->ts_a = ts.ts_a;
124 this->ts_t = ts.ts_t;
125 this->ts_dt = ts.ts_dt;
127}
128
129// BasicMethod
132 nInTheLoop(0), loopSize(0), rAnk(-1), sIze(-1),
133 refinedEntitiesPtr(nullptr), refinedFiniteElementsPtr(nullptr),
134 problemPtr(nullptr), fieldsPtr(nullptr), entitiesPtr(nullptr),
135 dofsPtr(nullptr), finiteElementsPtr(nullptr),
136 finiteElementsEntitiesPtr(nullptr), adjacenciesPtr(nullptr) {}
137
140
141 this->nInTheLoop = basic.nInTheLoop;
142 this->loopSize = basic.loopSize;
143 this->rAnk = basic.rAnk;
144 this->sIze = basic.sIze;
147 this->problemPtr = basic.problemPtr;
148 this->fieldsPtr = basic.fieldsPtr;
149 this->entitiesPtr = basic.entitiesPtr;
150 this->dofsPtr = basic.dofsPtr;
153 this->adjacenciesPtr = basic.adjacenciesPtr;
154 this->cacheWeakPtr = basic.cacheWeakPtr;
155
157}
158
161 if (preProcessHook) {
163 CHKERRG(ierr);
164 } else {
165 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
166 "should be implemented by user in derived class (preProcess)");
167 }
169}
172 if (postProcessHook) {
174 CHKERRG(ierr);
175 } else {
176 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
177 "should be implemented by user in derived class (postProcess)");
178 }
180}
183 if (operatorHook) {
184 ierr = operatorHook();
185 CHKERRG(ierr);
186 } else {
187 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
188 "should be implemented by user in derived class (operator)");
189 }
191}
192
194 VectorDouble &data,
195 const bool reset_dofs) {
197
198 // TODO: [CORE-60] Fix implementation not to use DOFs
199
200 auto get_nodes_field_data =
201 [&](boost::shared_ptr<FEDofEntity_multiIndex> &&dofs,
202 VectorDouble &nodes_data) {
204
205 auto bit_number = getFieldBitNumber(field_name);
206 auto &dofs_by_uid = dofs->get<Unique_mi_tag>();
207 auto dit =
208 dofs_by_uid.lower_bound(FieldEntity::getLocalUniqueIdCalculate(
209 bit_number, get_id_for_min_type<MBVERTEX>()));
210 auto hi_dit =
211 dofs_by_uid.upper_bound(FieldEntity::getLocalUniqueIdCalculate(
212 bit_number, get_id_for_max_type<MBVERTEX>()));
213
214 if (dit != hi_dit) {
215 auto &first_dof = **dit;
216 const int num_nodes = getNumberOfNodes();
217 const int nb_dof_idx = first_dof.getNbOfCoeffs();
218 const int max_nb_dofs = nb_dof_idx * num_nodes;
219 nodes_data.resize(max_nb_dofs, false);
220 nodes_data.clear();
221
222 for (; dit != hi_dit;) {
223 const auto &dof_ptr = *dit;
224 const auto &dof = *dof_ptr;
225 const auto &sn = *dof.getSideNumberPtr();
226 const int side_number = sn.side_number;
227
228 int pos = side_number * nb_dof_idx;
229 auto ent_filed_data_vec = dof.getEntFieldData();
230 for (int ii = 0; ii != nb_dof_idx; ++ii) {
231 nodes_data[pos] = ent_filed_data_vec[ii];
232 ++pos;
233 ++dit;
234 }
235 }
236
237 } else if (reset_dofs) {
238 nodes_data.resize(0, false);
239 }
240
242 };
243
244 return get_nodes_field_data(getDataDofsPtr(), data);
245
247}
248
249} // namespace MoFEM
#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 CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
MoFEM::TsCtx * ts_ctx
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
EntityHandle get_id_for_max_type()
EntityHandle get_id_for_min_type()
constexpr auto field_name
Data structure to exchange data between mofem and User Loop Methods.
int loopSize
local number oe methods to process
const FieldEntity_multiIndex * entitiesPtr
raw pointer to container of field entities
unsigned int getFieldBitNumber(std::string field_name) const
MoFEMErrorCode copyBasicMethod(const BasicMethod &basic)
Copy data from other base method to this base method.
const FiniteElement_multiIndex * finiteElementsPtr
raw pointer to container finite elements
const RefElement_multiIndex * refinedFiniteElementsPtr
container of mofem finite element entities
boost::function< MoFEMErrorCode()> preProcessHook
Hook function for pre-processing.
const DofEntity_multiIndex * dofsPtr
raw pointer container of dofs
int rAnk
processor rank
const Field_multiIndex * fieldsPtr
raw pointer to fields container
int nInTheLoop
number currently of processed method
const RefEntity_multiIndex * refinedEntitiesPtr
container of mofem dof entities
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
boost::function< MoFEMErrorCode()> operatorHook
Hook function for operator.
int sIze
number of processors in communicator
boost::weak_ptr< CacheTuple > cacheWeakPtr
boost::function< MoFEMErrorCode()> postProcessHook
Hook function for post-processing.
MoFEMErrorCode getNodeData(const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
auto getDataDofsPtr() const
auto getNumberOfNodes() const
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
data structure for ksp (linear solver) context
MoFEMErrorCode copyKsp(const KspMethod &ksp)
copy data form another method
KSPContext ksp_ctx
Context.
KSP ksp
KSP solver.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
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
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)
data structure for snes (nonlinear solver) context
MoFEMErrorCode copySnes(const SnesMethod &snes)
Copy snes data.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
SNESContext snes_ctx
SNES snes
snes solver
data structure for TS (time stepping) context
TS ts
time solver
PetscReal ts_t
time
PetscReal ts_dt
time step size
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
PetscReal ts_a
shift for U_t (see PETSc Time Solver)
PetscInt ts_step
time step number
MoFEMErrorCode copyTs(const TSMethod &ts)
Copy TS solver data.
Tao tao
tao solver
MoFEMErrorCode copyTao(const TaoMethod &tao)
Copy TAO data.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
base class for all interface classes