v0.8.23
DMMMoFEM.cpp
Go to the documentation of this file.
1 /* This file is part of MoFEM.
2  * MoFEM is free software: you can redistribute it and/or modify it under
3  * the terms of the GNU Lesser General Public License as published by the
4  * Free Software Foundation, either version 3 of the License, or (at your
5  * option) any later version.
6  *
7  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
8  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
9  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
10  * License for more details.
11  *
12  * You should have received a copy of the GNU Lesser General Public
13  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
14 
15 // #undef PETSC_VERSION_RELEASE
16 // #define PETSC_VERSION_RELEASE 1
17 
18 #if PETSC_VERSION_GE(3, 6, 0)
19 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
20 // #include <petsc/private/vecimpl.h> /*I "petscdm.h" I*/
21 #else
22 #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/
23 #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/
24 #endif
25 
26 #include <DMMoFEM.hpp>
27 
28 namespace MoFEM {
29 
31  : mField_ptr(PETSC_NULL), isProblemBuild(PETSC_FALSE),
32  isPartitioned(PETSC_FALSE), isSquareMatrix(PETSC_TRUE),
33  isSubDM(PETSC_FALSE), isCompDM(PETSC_FALSE), destroyProblem(PETSC_FALSE),
34  verbosity(VERBOSE), referenceNumber(0) {}
35 
37  // cerr << "Snes " << snesCtx.use_count() << endl;
38  // cerr << "Destroy DMCtx" << endl;
39 }
40 
42  UnknownInterface **iface) const {
44  *iface = NULL;
45  if (uuid == IDD_DMCTX) {
46  *iface = const_cast<DMCtx *>(this);
48  }
49  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown interface");
51 }
52 
53 PetscErrorCode DMRegister_MoFEM(const char sname[]) {
55  CHKERR DMRegister(sname, DMCreate_MoFEM);
57 }
58 
59 PetscErrorCode DMSetOperators_MoFEM(DM dm) {
60  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
62 
63  dm->ops->createglobalvector = DMCreateGlobalVector_MoFEM;
64  dm->ops->createlocalvector = DMCreateLocalVector_MoFEM;
65  dm->ops->creatematrix = DMCreateMatrix_MoFEM;
66  dm->ops->setup = DMSetUp_MoFEM;
67  dm->ops->destroy = DMDestroy_MoFEM;
68  dm->ops->setfromoptions = DMSetFromOptions_MoFEM;
69  dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_MoFEM;
70  dm->ops->globaltolocalend = DMGlobalToLocalEnd_MoFEM;
71  dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_MoFEM;
72  dm->ops->localtoglobalend = DMLocalToGlobalEnd_MoFEM;
73  dm->ops->createfieldis = DMCreateFieldIS_MoFEM;
74 
75  // Default matrix type
76  CHKERR DMSetMatType(dm, MATMPIAIJ);
77 
79 }
80 
81 PetscErrorCode DMCreate_MoFEM(DM dm) {
82  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
84  dm->data = new DMCtx();
87 }
88 
89 PetscErrorCode DMDestroy_MoFEM(DM dm) {
90  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
91  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
93  if (((DMCtx *)dm->data)->referenceNumber == 0) {
94  if (dm_field->destroyProblem) {
95  if (dm_field->mField_ptr->check_problem(dm_field->problemName)) {
96  dm_field->mField_ptr->delete_problem(dm_field->problemName);
97  } // else problem has to be deleted by the user
98  }
99  // cerr << "Destroy " << dm_field->problemName << endl;
100  delete ((DMCtx *)dm->data);
101  } else {
102  // cerr << "Dereference " << dm_field->problemName << " " <<
103  // ((DMCtx*)dm->data)->referenceNumber << endl;
104  (((DMCtx *)dm->data)->referenceNumber)--;
105  }
107 }
108 
109 PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr,
110  const char problem_name[],
111  const MoFEM::BitRefLevel bit_level,
112  const MoFEM::BitRefLevel bit_mask) {
114 
115  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
116  if (!dm->data) {
117  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
118  "data structure for MoFEM not yet created");
119  }
120  if (!m_field_ptr) {
121  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
122  "DM function not implemented into MoFEM");
123  }
124  dm_field->mField_ptr = m_field_ptr;
125  dm_field->problemName = problem_name;
126  if (!m_field_ptr->check_problem(dm_field->problemName)) {
127  // problem is not defined, declare problem internally set bool to
128  // destroyProblem problem with DM
129  dm_field->destroyProblem = PETSC_TRUE;
130  CHKERR dm_field->mField_ptr->add_problem(dm_field->problemName, MF_EXCL,
131  dm_field->verbosity);
132  } else {
133  dm_field->destroyProblem = PETSC_FALSE;
134  }
136  dm_field->problemName, bit_level);
138  dm_field->problemName, bit_mask);
139  dm_field->kspCtx =
140  boost::shared_ptr<KspCtx>(new KspCtx(*m_field_ptr, problem_name));
141  dm_field->snesCtx =
142  boost::shared_ptr<SnesCtx>(new SnesCtx(*m_field_ptr, problem_name));
143  dm_field->tsCtx =
144  boost::shared_ptr<TsCtx>(new TsCtx(*m_field_ptr, problem_name));
145 
146  MPI_Comm comm;
147  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
148  int result = 0;
149  MPI_Comm_compare(comm, m_field_ptr->get_comm(), &result);
150  // std::cerr << result << " " << MPI_IDENT << " " << MPI_CONGRUENT << " " <<
151  // MPI_SIMILAR << " " << MPI_UNEQUAL << std::endl;
152  if (result > MPI_CONGRUENT) {
153  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
154  "MoFEM and DM using different communicators");
155  }
156  MPI_Comm_size(comm, &dm_field->sIze);
157  MPI_Comm_rank(comm, &dm_field->rAnk);
158 
159  // problem structure
160  CHKERR dm_field->mField_ptr->get_problem(dm_field->problemName,
161  &dm_field->problemPtr);
162 
164 }
165 
166 PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[]) {
168 
169  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
170  if (!dm->data) {
171  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
172  "data structure for MoFEM not yet created");
173  }
174  CHKERR DMMoFEMCreateMoFEM(subdm, dm_field->mField_ptr, problem_name,
175  dm_field->problemPtr->getBitRefLevel());
176 
177  DMCtx *subdm_field = (DMCtx *)subdm->data;
178  subdm_field->isSubDM = PETSC_TRUE;
179  subdm_field->problemMainOfSubPtr = dm_field->problemPtr;
180  subdm_field->isPartitioned = dm_field->isPartitioned;
181  subdm_field->isSquareMatrix = PETSC_FALSE;
182  subdm->ops->setup = DMSubDMSetUp_MoFEM;
183 
185 }
186 
187 PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[],
188  EntityType lo_type, EntityType hi_type) {
189  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
191  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
192  if (!dm->data) {
193  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
194  "data structure for MoFEM not yet created");
195  }
196  if (!dm_field->isSubDM) {
197  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
198  }
199  dm_field->rowFields.push_back(field_name);
200  if (lo_type != MBVERTEX || hi_type != MBMAXTYPE) {
201  if (!dm_field->mapTypeRow)
202  dm_field->mapTypeRow = boost::make_shared<
203  std::map<std::string, std::pair<EntityType, EntityType>>>();
204  (*dm_field->mapTypeRow)[field_name] =
205  std::pair<EntityType, EntityType>(lo_type, hi_type);
206  }
208 }
209 
210 PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[],
211  EntityType lo_type, EntityType hi_type) {
212  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
214  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
215  if (!dm->data) {
216  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
217  "data structure for MoFEM not yet created");
218  }
219  if (!dm_field->isSubDM) {
220  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
221  }
222  dm_field->colFields.push_back(field_name);
223  if (lo_type != MBVERTEX || hi_type != MBMAXTYPE) {
224  if (!dm_field->mapTypeCol)
225  dm_field->mapTypeCol = boost::make_shared<
226  std::map<std::string, std::pair<EntityType, EntityType>>>();
227  (*dm_field->mapTypeCol)[field_name] =
228  std::pair<EntityType, EntityType>(lo_type, hi_type);
229  }
231 }
232 
233 PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm) {
235  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
237  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
238  *is_sub_dm = dm_field->isSubDM;
240 }
241 
242 PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is) {
244  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
246  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
247  if (dm_field->isSubDM != PETSC_TRUE) {
248  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
249  "This DM is not created as a SubDM");
250  }
251  if (dm_field->isProblemBuild != PETSC_TRUE) {
252  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Problem is not build");
253  }
254  boost::shared_ptr<Problem::SubProblemData> sub_data =
255  dm_field->problemPtr->getSubData();
256  CHKERR sub_data->getRowIs(is);
258 }
259 
260 PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is) {
262  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
264  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
265  if (dm_field->isSubDM != PETSC_TRUE) {
266  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
267  "This DM is not created as a SubDM");
268  }
269  if (dm_field->isProblemBuild != PETSC_TRUE) {
270  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Problem is not build");
271  }
272  boost::shared_ptr<Problem::SubProblemData> sub_data =
273  dm_field->problemPtr->getSubData();
274  CHKERR sub_data->getColIs(is);
276 }
277 
278 PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[]) {
279  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
281  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
282  if (!dm->data) {
283  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
284  "data structure for MoFEM not yet created");
285  }
286  if (!dm_field->isCompDM) {
287  dm_field->isCompDM = PETSC_TRUE;
288  }
289  dm_field->rowCompPrb.push_back(prb_name);
290  if (dm_field->isSquareMatrix) {
291  dm_field->colCompPrb.push_back(prb_name);
292  }
294 }
295 
296 PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[]) {
297  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
299  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
300  if (!dm->data) {
301  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
302  "data structure for MoFEM not yet created");
303  }
304  if (!dm_field->isCompDM) {
305  dm_field->isCompDM = PETSC_TRUE;
306  }
307  if (dm_field->isSquareMatrix) {
308  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
309  "No need to add problem on column when problem block structurally "
310  "symmetric");
311  }
312  dm_field->colCompPrb.push_back(prb_name);
314 }
315 
316 PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm) {
318  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
320  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
321  *is_comp_dm = dm_field->isCompDM;
323 }
324 
325 PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr) {
326  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
328  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
329  if (!dm->data) {
330  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
331  "data structure for MoFEM not yet created");
332  }
333  *m_field_ptr = dm_field->mField_ptr;
335 }
336 
337 PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr) {
338  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
340  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
341  if (!dm->data) {
342  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
343  "data structure for MoFEM not yet created");
344  }
345  *problem_ptr = dm_field->problemPtr;
347 }
348 
349 PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem) {
351  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
353  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
354  dm_field->destroyProblem = destroy_problem;
356 }
357 
358 PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem) {
360  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
362  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
363  *destroy_problem = dm_field->destroyProblem;
365 }
366 
367 PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem) {
369  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
371  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
372  dm_field->isSquareMatrix = square_problem;
374 }
375 
376 PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, const char fe_name[]) {
378  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
380  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
382  ->resolveSharedFiniteElements(dm_field->problemPtr, fe_name);
384 }
385 
386 PetscErrorCode DMMoFEMResolveSharedEntities(DM dm, const char fe_name[]) {
387  return DMMoFEMResolveSharedFiniteElements(dm, fe_name);
388 }
389 
390 PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, const char fe_name[],
391  PetscLayout *layout) {
393  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
395  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
396 
397  MPI_Comm comm;
398  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
399  CHKERR dm_field->problemPtr->getNumberOfElementsByNameAndPart(comm, fe_name,
400  layout);
402 }
403 
404 PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem) {
407  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
409  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
410  *square_problem = dm_field->isSquareMatrix;
412 }
413 
414 PetscErrorCode DMMoFEMAddElement(DM dm, const char fe_name[]) {
415  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
417  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
419  dm_field->problemName, fe_name);
420  CHKERRG(ierr);
422 }
423 
424 PetscErrorCode DMMoFEMUnSetElement(DM dm, const char fe_name[]) {
425  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
427  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
429  dm_field->problemName, fe_name);
430  CHKERRG(ierr);
432 }
433 
434 PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode,
435  ScatterMode scatter_mode) {
437  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
439  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
440  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setLocalGhostVector(
441  dm_field->problemPtr, COL, l, mode, scatter_mode);
442  CHKERRG(ierr);
444 }
445 
446 PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode,
447  ScatterMode scatter_mode) {
448  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
450  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
451  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setGlobalGhostVector(
452  dm_field->problemPtr, COL, g, mode, scatter_mode);
453  CHKERRG(ierr);
455 }
456 
457 PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
458  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
460  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
462  dm_field->problemPtr, *method);
463  CHKERRG(ierr);
465 }
466 
467 PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
468  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
470  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
472  dm_field->problemPtr, *method);
473  CHKERRG(ierr);
475 }
476 
477 PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[],
478  MoFEM::FEMethod *method,
479  int low_rank, int up_rank) {
481  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
482  ierr = dm_field->mField_ptr->loop_finite_elements(
483  dm_field->problemPtr, fe_name, *method, low_rank, up_rank);
484  CHKERRG(ierr);
486 }
487 
488 PetscErrorCode
489 DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const std::string fe_name,
490  boost::shared_ptr<MoFEM::FEMethod> method,
491  int low_rank, int up_rank) {
492  return DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name.c_str(), method.get(),
493  low_rank, up_rank);
494 }
495 
496 PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[],
497  MoFEM::FEMethod *method) {
498  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
500  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
501  ierr = DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name, method,
502  dm_field->rAnk, dm_field->rAnk);
503  CHKERRG(ierr);
505 }
506 
507 PetscErrorCode
508 DMoFEMLoopFiniteElements(DM dm, const std::string fe_name,
509  boost::shared_ptr<MoFEM::FEMethod> method) {
510  return DMoFEMLoopFiniteElements(dm, fe_name.c_str(), method.get());
511 }
512 
513 PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[],
514  MoFEM::DofMethod *method) {
515  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
517  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
518  ierr =
519  dm_field->mField_ptr->loop_dofs(dm_field->problemPtr, field_name, COL,
520  *method, dm_field->rAnk, dm_field->rAnk);
521  CHKERRG(ierr);
523 }
524 
525 template <class S, class T0, class T1, class T2>
526 static PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, S fe_name, T0 method,
527  T1 pre_only, T2 post_only) {
528  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
530  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
531  if (pre_only) {
532  dm_field->kspCtx->get_preProcess_to_do_Rhs().push_back(pre_only);
533  }
534  if (method) {
535  dm_field->kspCtx->get_loops_to_do_Rhs().push_back(
536  PairNameFEMethodPtr(fe_name, method));
537  }
538  if (post_only) {
539  dm_field->kspCtx->get_postProcess_to_do_Rhs().push_back(post_only);
540  }
541  CHKERR DMKSPSetComputeRHS(dm, KspRhs, dm_field->kspCtx.get());
543 }
544 
545 PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[],
546  MoFEM::FEMethod *method,
547  MoFEM::BasicMethod *pre_only,
548  MoFEM::BasicMethod *post_only) {
549  return DMMoFEMKSPSetComputeRHS<const char *, MoFEM::FEMethod *,
551  dm, fe_name, method, pre_only, post_only);
552 }
553 
554 PetscErrorCode
555 DMMoFEMKSPSetComputeRHS(DM dm, const std::string fe_name,
556  boost::shared_ptr<MoFEM::FEMethod> method,
557  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
558  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
559  return DMMoFEMKSPSetComputeRHS<const std::string,
560  boost::shared_ptr<MoFEM::FEMethod>,
561  boost::shared_ptr<MoFEM::BasicMethod>,
562  boost::shared_ptr<MoFEM::BasicMethod>>(
563  dm, fe_name, method, pre_only, post_only);
564 }
565 
566 template <class S, class T0, class T1, class T2>
567 static PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, S fe_name, T0 method,
568  T1 pre_only, T2 post_only) {
569  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
571  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
572  if (pre_only) {
573  dm_field->kspCtx->get_preProcess_to_do_Mat().push_back(pre_only);
574  }
575  if (method) {
576  dm_field->kspCtx->get_loops_to_do_Mat().push_back(
577  PairNameFEMethodPtr(fe_name, method));
578  }
579  if (post_only) {
580  dm_field->kspCtx->get_postProcess_to_do_Mat().push_back(post_only);
581  }
582  CHKERR DMKSPSetComputeOperators(dm, KspMat, dm_field->kspCtx.get());
584 }
585 
586 PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[],
587  MoFEM::FEMethod *method,
588  MoFEM::BasicMethod *pre_only,
589  MoFEM::BasicMethod *post_only) {
590  return DMMoFEMKSPSetComputeOperators<const char *, MoFEM::FEMethod *,
593  dm, fe_name, method, pre_only, post_only);
594 }
595 
596 PetscErrorCode
597 DMMoFEMKSPSetComputeOperators(DM dm, const std::string fe_name,
598  boost::shared_ptr<MoFEM::FEMethod> method,
599  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
600  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
601  return DMMoFEMKSPSetComputeOperators<const std::string,
602  boost::shared_ptr<MoFEM::FEMethod>>(
603  dm, fe_name, method, pre_only, post_only);
604 }
605 
606 template <class S, class T0, class T1, class T2>
607 static PetscErrorCode DMMoFEMSNESSetFunction(DM dm, S fe_name, T0 method,
608  T1 pre_only, T2 post_only) {
609  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
611  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
612  if (pre_only) {
613  dm_field->snesCtx->get_preProcess_to_do_Rhs().push_back(pre_only);
614  }
615  if (method) {
616  dm_field->snesCtx->get_loops_to_do_Rhs().push_back(
617  PairNameFEMethodPtr(fe_name, method));
618  }
619  if (post_only) {
620  dm_field->snesCtx->get_postProcess_to_do_Rhs().push_back(post_only);
621  }
622  CHKERR DMSNESSetFunction(dm, SnesRhs, dm_field->snesCtx.get());
624 }
625 
626 PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[],
627  MoFEM::FEMethod *method,
628  MoFEM::BasicMethod *pre_only,
629  MoFEM::BasicMethod *post_only) {
630  return DMMoFEMSNESSetFunction<const char *, MoFEM::FEMethod *,
632  dm, fe_name, method, pre_only, post_only);
633 }
634 
635 PetscErrorCode
636 DMMoFEMSNESSetFunction(DM dm, const std::string fe_name,
637  boost::shared_ptr<MoFEM::FEMethod> method,
638  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
639  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
640  return DMMoFEMSNESSetFunction<const std::string,
641  boost::shared_ptr<MoFEM::FEMethod>,
642  boost::shared_ptr<MoFEM::BasicMethod>,
643  boost::shared_ptr<MoFEM::BasicMethod>>(
644  dm, fe_name, method, pre_only, post_only);
645 }
646 
647 template <class S, class T0, class T1, class T2>
648 static PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, S fe_name, T0 method,
649  T1 pre_only, T2 post_only) {
650  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
652  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
653  if (pre_only) {
654  dm_field->snesCtx->get_preProcess_to_do_Mat().push_back(pre_only);
655  }
656  if (method) {
657  dm_field->snesCtx->get_loops_to_do_Mat().push_back(
658  PairNameFEMethodPtr(fe_name, method));
659  }
660  if (post_only) {
661  dm_field->snesCtx->get_postProcess_to_do_Mat().push_back(post_only);
662  }
663  CHKERR DMSNESSetJacobian(dm, SnesMat, dm_field->snesCtx.get());
665 }
666 
667 PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[],
668  MoFEM::FEMethod *method,
669  MoFEM::BasicMethod *pre_only,
670  MoFEM::BasicMethod *post_only) {
671  return DMMoFEMSNESSetJacobian<const char *, MoFEM::FEMethod *,
673  dm, fe_name, method, pre_only, post_only);
674 }
675 
676 PetscErrorCode
677 DMMoFEMSNESSetJacobian(DM dm, const std::string fe_name,
678  boost::shared_ptr<MoFEM::FEMethod> method,
679  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
680  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
681  return DMMoFEMSNESSetJacobian<const std::string,
682  boost::shared_ptr<MoFEM::FEMethod>,
683  boost::shared_ptr<MoFEM::BasicMethod>,
684  boost::shared_ptr<MoFEM::BasicMethod>>(
685  dm, fe_name, method, pre_only, post_only);
686 }
687 
688 template <class S, class T0, class T1, class T2>
689 static PetscErrorCode DMMoFEMTSSetIFunction(DM dm, S fe_name, T0 method,
690  T1 pre_only, T2 post_only) {
691  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
693  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
694  if (pre_only) {
695  dm_field->tsCtx->get_preProcess_to_do_IFunction().push_back(pre_only);
696  }
697  if (method) {
698  dm_field->tsCtx->get_loops_to_do_IFunction().push_back(
699  PairNameFEMethodPtr(fe_name, method));
700  }
701  if (post_only) {
702  dm_field->tsCtx->get_postProcess_to_do_IFunction().push_back(post_only);
703  }
704  CHKERR DMTSSetIFunction(dm, TsSetIFunction, dm_field->tsCtx.get());
706 }
707 
708 PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[],
709  MoFEM::FEMethod *method,
710  MoFEM::BasicMethod *pre_only,
711  MoFEM::BasicMethod *post_only) {
712  return DMMoFEMTSSetIFunction<const char *, MoFEM::FEMethod *,
714  dm, fe_name, method, pre_only, post_only);
716 }
717 
718 PetscErrorCode
719 DMMoFEMTSSetIFunction(DM dm, const std::string fe_name,
720  boost::shared_ptr<MoFEM::FEMethod> method,
721  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
722  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
723  return DMMoFEMTSSetIFunction<const std::string,
724  boost::shared_ptr<MoFEM::FEMethod>,
725  boost::shared_ptr<MoFEM::BasicMethod>,
726  boost::shared_ptr<MoFEM::BasicMethod>>(
727  dm, fe_name, method, pre_only, post_only);
729 }
730 
731 template <class S, class T0, class T1, class T2>
732 static PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, S fe_name, T0 method,
733  T1 pre_only, T2 post_only) {
734  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
736  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
737  if (pre_only) {
738  dm_field->tsCtx->get_preProcess_to_do_IJacobian().push_back(pre_only);
739  }
740  if (method) {
741  dm_field->tsCtx->get_loops_to_do_IJacobian().push_back(
742  PairNameFEMethodPtr(fe_name, method));
743  }
744  if (post_only) {
745  dm_field->tsCtx->get_postProcess_to_do_IJacobian().push_back(post_only);
746  }
747  CHKERR DMTSSetIJacobian(dm, TsSetIJacobian, dm_field->tsCtx.get());
749 }
750 
751 PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const char fe_name[],
752  MoFEM::FEMethod *method,
753  MoFEM::BasicMethod *pre_only,
754  MoFEM::BasicMethod *post_only) {
755  return DMMoFEMTSSetIJacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
756  MoFEM::BasicMethod *>(dm, fe_name, method,
757  pre_only, post_only);
758 }
759 
760 PetscErrorCode
761 DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name,
762  boost::shared_ptr<MoFEM::FEMethod> method,
763  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
764  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
765  return DMMoFEMTSSetIJacobian<const std::string,
766  boost::shared_ptr<MoFEM::FEMethod>,
767  boost::shared_ptr<MoFEM::BasicMethod>,
768  boost::shared_ptr<MoFEM::BasicMethod>>(
769  dm, fe_name, method, pre_only, post_only);
770 }
771 
772 template <class S, class T0, class T1, class T2>
773 static PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, S fe_name, T0 method,
774  T1 pre_only, T2 post_only) {
775  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
777  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
778  if (pre_only)
779  dm_field->tsCtx->get_preProcess_to_do_RHSFunction().push_back(pre_only);
780  if (method)
781  dm_field->tsCtx->get_loops_to_do_RHSFunction().push_back(
782  PairNameFEMethodPtr(fe_name, method));
783  if (post_only)
784  dm_field->tsCtx->get_postProcess_to_do_RHSFunction().push_back(post_only);
785  CHKERR DMTSSetRHSFunction(dm, TsSetRHSFunction, dm_field->tsCtx.get());
787 }
788 
789 PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const char fe_name[],
790  MoFEM::FEMethod *method,
791  MoFEM::BasicMethod *pre_only,
792  MoFEM::BasicMethod *post_only) {
793  return DMMoFEMTSSetRHSFunction<const char *, MoFEM::FEMethod *,
795  dm, fe_name, method, pre_only, post_only);
797 }
798 
799 PetscErrorCode
800 DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name,
801  boost::shared_ptr<MoFEM::FEMethod> method,
802  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
803  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
804  return DMMoFEMTSSetRHSFunction<const std::string,
805  boost::shared_ptr<MoFEM::FEMethod>,
806  boost::shared_ptr<MoFEM::BasicMethod>,
807  boost::shared_ptr<MoFEM::BasicMethod>>(
808  dm, fe_name, method, pre_only, post_only);
810 }
811 
812 template <class S, class T0, class T1, class T2>
813 static PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, S fe_name, T0 method,
814  T1 pre_only, T2 post_only) {
815  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
817  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
818  if (pre_only)
819  dm_field->tsCtx->get_preProcess_to_do_RHSFunction().push_back(pre_only);
820  if (method)
821  dm_field->tsCtx->get_loops_to_do_RHSFunction().push_back(
822  PairNameFEMethodPtr(fe_name, method));
823  if (post_only)
824  dm_field->tsCtx->get_postProcess_to_do_RHSFunction().push_back(post_only);
825  CHKERR DMTSSetRHSJacobian(dm, TsSetRHSJacobian, dm_field->tsCtx.get());
827 }
828 
829 PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const char fe_name[],
830  MoFEM::FEMethod *method,
831  MoFEM::BasicMethod *pre_only,
832  MoFEM::BasicMethod *post_only) {
833  return DMMoFEMTSSetRHSJacobian<const char *, MoFEM::FEMethod *,
835  dm, fe_name, method, pre_only, post_only);
837 }
838 
839 PetscErrorCode
840 DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name,
841  boost::shared_ptr<MoFEM::FEMethod> method,
842  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
843  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
844  return DMMoFEMTSSetRHSJacobian<const std::string,
845  boost::shared_ptr<MoFEM::FEMethod>,
846  boost::shared_ptr<MoFEM::BasicMethod>,
847  boost::shared_ptr<MoFEM::BasicMethod>>(
848  dm, fe_name, method, pre_only, post_only);
850 }
851 
852 template <class S, class T0, class T1, class T2>
853 static PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, S fe_name, T0 method,
854  T1 pre_only, T2 post_only) {
855  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
857  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
858  if (pre_only)
859  dm_field->tsCtx->get_preProcess_to_do_Monitor().push_back(pre_only);
860  if (method)
861  dm_field->tsCtx->get_loops_to_do_Monitor().push_back(
862  PairNameFEMethodPtr(fe_name, method));
863  if (post_only)
864  dm_field->tsCtx->get_postProcess_to_do_Monitor().push_back(post_only);
865  CHKERR TSMonitorSet(ts, TsMonitorSet, dm_field->tsCtx.get(), PETSC_NULL);
867 }
868 
869 PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const char fe_name[],
870  MoFEM::FEMethod *method,
871  MoFEM::BasicMethod *pre_only,
872  MoFEM::BasicMethod *post_only) {
873  return DMMoFEMTSSetMonitor<const char *, MoFEM::FEMethod *,
875  dm, ts, fe_name, method, pre_only, post_only);
877 }
878 
879 PetscErrorCode
880 DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name,
881  boost::shared_ptr<MoFEM::FEMethod> method,
882  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
883  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
884  return DMMoFEMTSSetMonitor<const std::string,
885  boost::shared_ptr<MoFEM::FEMethod>,
886  boost::shared_ptr<MoFEM::BasicMethod>,
887  boost::shared_ptr<MoFEM::BasicMethod>>(
888  dm, ts, fe_name, method, pre_only, post_only);
890 }
891 
892 PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx) {
893  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
895  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
896  *ksp_ctx = dm_field->kspCtx.get();
898 }
899 
900 PetscErrorCode
901 DMMoFEMGetKspCtx(DM dm, const boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx) {
902  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
904  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
905  const_cast<boost::shared_ptr<MoFEM::KspCtx> &>(ksp_ctx) = dm_field->kspCtx;
907 }
908 
909 PetscErrorCode DMMoFEMSetKspCtx(DM dm,
910  boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx) {
911  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
913  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
914  dm_field->kspCtx = ksp_ctx;
916 }
917 
918 PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx) {
919  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
921  DMCtx *dm_field = (DMCtx *)dm->data;
922  *snes_ctx = dm_field->snesCtx.get();
924 }
925 
926 PetscErrorCode
927 DMMoFEMGetSnesCtx(DM dm, const boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx) {
928  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
930  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
931  const_cast<boost::shared_ptr<MoFEM::SnesCtx> &>(snes_ctx) = dm_field->snesCtx;
933 }
934 
935 PetscErrorCode DMMoFEMSetSnesCtx(DM dm,
936  boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx) {
937  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
939  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
940  dm_field->snesCtx = snes_ctx;
942 }
943 
944 /** get if read mesh is partitioned
945  * \ingroup dm
946  */
947 PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned) {
948  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
950  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
951  dm_field->isPartitioned = is_partitioned;
953 }
954 
955 /** get if read mesh is partitioned
956  * \ingroup dm
957  */
958 PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned) {
959  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
961  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
962  *is_partitioned = dm_field->isPartitioned;
964 }
965 
966 PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx) {
967  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
969  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
970  *ts_ctx = dm_field->tsCtx.get();
972 }
973 
974 PetscErrorCode DMMoFEMGetTsCtx(DM dm,
975  const boost::shared_ptr<MoFEM::TsCtx> &ts_ctx) {
976  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
978  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
979  const_cast<boost::shared_ptr<MoFEM::TsCtx> &>(ts_ctx) = dm_field->tsCtx;
981 }
982 
983 PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr<MoFEM::TsCtx> &ts_ctx) {
984  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
986  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
987  dm_field->tsCtx = ts_ctx;
989 }
990 
991 PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g) {
992  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
994  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
995  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
996  dm_field->problemName, COL, g);
998 }
999 
1000 PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, SmartPetscObj<Vec> &g_ptr) {
1001  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1003  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1004  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1005  dm_field->problemName, COL, g_ptr);
1007 }
1008 
1009 PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l) {
1010  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1012  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1013  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateSeq(
1014  dm_field->problemName, COL, l);
1016 }
1017 
1018 PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M) {
1019  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1021  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1022  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1024  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1025  M);
1026  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1028  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1029  M);
1030  } else {
1031  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1032  "Matrix type not implemented");
1033  }
1035 }
1036 
1037 PetscErrorCode DMCreateMatrix_MoFEM(DM dm, SmartPetscObj<Mat> &M) {
1038  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1040  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1041  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1043  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1044  M);
1045  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1047  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1048  M);
1049  } else {
1050  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1051  "Matrix type not implemented");
1052  }
1054 }
1055 
1056 #if PETSC_VERSION_GE(3, 7, 0)
1057 PetscErrorCode DMSetFromOptions_MoFEM(PetscOptionItems *PetscOptionsObject,
1058  DM dm) {
1059 #elif PETSC_VERSION_GE(3, 5, 3)
1060 PetscErrorCode DMSetFromOptions_MoFEM(PetscOptions *PetscOptionsObject, DM dm) {
1061 #else
1062 PetscErrorCode DMSetFromOptions_MoFEM(DM dm) {
1063 #endif
1064 
1065  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1067  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1068 #if PETSC_VERSION_GE(3, 5, 3)
1069  ierr = PetscOptionsHead(PetscOptionsObject, "DMMoFEM Options");
1070  CHKERRG(ierr);
1071 #else
1072  ierr = PetscOptionsHead("DMMoFEM Options");
1073  CHKERRG(ierr);
1074 #endif
1075  ierr = PetscOptionsBool("-dm_is_partitioned",
1076  "set if mesh is partitioned (works which native MOAB "
1077  "file format, i.e. h5m",
1078  "DMSetUp", dm_field->isPartitioned,
1079  &dm_field->isPartitioned, NULL);
1080  CHKERRG(ierr);
1082 }
1083 
1084 PetscErrorCode DMSetUp_MoFEM(DM dm) {
1085  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1086  ProblemsManager *prb_mng_ptr;
1088  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1089  CHKERR dm_field->mField_ptr->getInterface(prb_mng_ptr);
1090 
1091  if (dm_field->isCompDM) {
1092  // It is composite probelm
1093  CHKERR prb_mng_ptr->buildCompsedProblem(
1094  dm_field->problemName, dm_field->rowCompPrb, dm_field->colCompPrb,
1095  dm_field->isSquareMatrix == PETSC_TRUE, dm_field->verbosity);
1096  } else {
1097  if (dm_field->isPartitioned) {
1099  dm_field->problemName, dm_field->isSquareMatrix == PETSC_TRUE,
1100  dm_field->verbosity);
1101  } else {
1102  CHKERR prb_mng_ptr->buildProblem(dm_field->problemName,
1103  dm_field->isSquareMatrix == PETSC_TRUE,
1104  dm_field->verbosity);
1105  CHKERR prb_mng_ptr->partitionProblem(dm_field->problemName,
1106  dm_field->verbosity);
1107  }
1108  }
1109 
1110  // Partition finite elements
1111  if (dm_field->isPartitioned) {
1112  CHKERR prb_mng_ptr->partitionFiniteElements(
1113  dm_field->problemName, true, 0, dm_field->sIze, dm_field->verbosity);
1115  dm_field->problemName, dm_field->verbosity);
1116  } else {
1117  // partition finite elemnets
1118  CHKERR prb_mng_ptr->partitionFiniteElements(dm_field->problemName, false,
1119  -1, -1, dm_field->verbosity);
1120  // Get ghost DOFs
1121  CHKERR prb_mng_ptr->partitionGhostDofs(dm_field->problemName,
1122  dm_field->verbosity);
1123  }
1124 
1125  // Set flag that problem is build and partitioned
1126  dm_field->isProblemBuild = PETSC_TRUE;
1127 
1129 }
1130 
1131 PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm) {
1132  PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
1133  ProblemsManager *prb_mng_ptr;
1135 
1136  DMCtx *subdm_field = static_cast<DMCtx *>(subdm->data);
1137 
1138  // build sub dm problem
1139  CHKERR subdm_field->mField_ptr->getInterface(prb_mng_ptr);
1140 
1141  map<std::string, std::pair<EntityType, EntityType>> *entity_map_row = nullptr;
1142  map<std::string, std::pair<EntityType, EntityType>> *entity_map_col = nullptr;
1143 
1144  if (subdm_field->mapTypeRow)
1145  entity_map_row = subdm_field->mapTypeRow.get();
1146  if (subdm_field->mapTypeCol)
1147  entity_map_row = subdm_field->mapTypeCol.get();
1148 
1149  CHKERR prb_mng_ptr->buildSubProblem(
1150  subdm_field->problemName, subdm_field->rowFields, subdm_field->colFields,
1151  subdm_field->problemMainOfSubPtr->getName(),
1152  subdm_field->isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1153  subdm_field->verbosity);
1154 
1155  // partition problem
1156  subdm_field->isPartitioned = subdm_field->isPartitioned;
1157  if (subdm_field->isPartitioned) {
1158  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, true,
1159  0, subdm_field->sIze,
1160  subdm_field->verbosity);
1161  // set ghost nodes
1163  subdm_field->problemName, subdm_field->verbosity);
1164  } else {
1165  // partition finite elements
1166  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, false,
1167  -1, -1, subdm_field->verbosity);
1168  // set ghost nodes
1169  CHKERR prb_mng_ptr->partitionGhostDofs(subdm_field->problemName,
1170  subdm_field->verbosity);
1171  }
1172 
1173  subdm_field->isProblemBuild = PETSC_TRUE;
1174 
1176 }
1177 
1178 PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec g, InsertMode mode,
1179  Vec l) {
1180  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1182  CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
1184 }
1185 
1186 PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec g, InsertMode mode, Vec l) {
1188  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1190 
1191  CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
1192 
1193  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1194  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1195  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1196 
1197  double *array_loc, *array_glob;
1198  CHKERR VecGetArray(l, &array_loc);
1199  CHKERR VecGetArray(g, &array_glob);
1200  switch (mode) {
1201  case INSERT_VALUES:
1202  cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1203  break;
1204  case ADD_VALUES:
1205  cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1206  break;
1207  default:
1208  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1209  }
1210  CHKERR VecRestoreArray(l, &array_loc);
1211  CHKERR VecRestoreArray(g, &array_glob);
1213 }
1214 
1215 PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM dm, Vec l, InsertMode mode,
1216  Vec g) {
1217 
1218  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1220 
1221  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1222  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1223  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1224 
1225  double *array_loc, *array_glob;
1226  CHKERR VecGetArray(l, &array_loc);
1227  CHKERR VecGetArray(g, &array_glob);
1228  switch (mode) {
1229  case INSERT_VALUES:
1230  cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1231  break;
1232  case ADD_VALUES:
1233  cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1234  break;
1235  default:
1236  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1237  }
1238  CHKERR VecRestoreArray(l, &array_loc);
1239  CHKERR VecRestoreArray(g, &array_glob);
1240 
1242 }
1243 
1244 PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec l, InsertMode mode, Vec g) {
1245  //
1248 }
1249 
1250 PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields,
1251  char ***fieldNames, IS **fields) {
1252  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1254 
1255  if (numFields) {
1256  *numFields = 0;
1257  }
1258  if (fieldNames) {
1259  *fieldNames = NULL;
1260  }
1261  if (fields) {
1262  *fields = NULL;
1263  }
1264 
1265  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1266  const Field_multiIndex *fields_ptr;
1267  CHKERR dm_field->mField_ptr->get_fields(&fields_ptr);
1268  Field_multiIndex::iterator fit, hi_fit;
1269  fit = fields_ptr->begin();
1270  hi_fit = fields_ptr->end();
1271  *numFields = std::distance(fit, hi_fit);
1272 
1273  if (fieldNames) {
1274  CHKERR PetscMalloc1(*numFields, fieldNames);
1275  }
1276  if (fields) {
1277  CHKERR PetscMalloc1(*numFields, fields);
1278  }
1279 
1280  for (int f = 0; fit != hi_fit; fit++, f++) {
1281  if (fieldNames) {
1282  CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1283  (char **)&(*fieldNames)[f]);
1284  }
1285  if (fields) {
1286  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1287  ->isCreateProblemFieldAndRank(
1288  dm_field->problemPtr->getName(), ROW, fit->get()->getName(), 0,
1289  fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1290  }
1291  }
1292 
1294 }
1295 
1296 PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[],
1297  IS *is) {
1298  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1300  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1301  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1302  ->isCreateProblemFieldAndRank(dm_field->problemPtr->getName(), ROW,
1303  field_name, 0, 1000, is);
1305 }
1306 
1307 PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb) {
1308  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1310  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1311  dm_field->verbosity = verb;
1313 }
1314 
1315 } // namespace MoFEM
std::string problemName
Problem name.
Definition: DMMoFEM.hpp:893
PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is)
get sub problem is
Definition: DMMMoFEM.cpp:242
structure for User Loop Methods on finite elementsIt can be used to calculate stiffness matrices,...
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition: DMMMoFEM.cpp:626
MoFEM interface unique ID.
Problem manager is used to build and partition problems.
PetscBool isSquareMatrix
true if rows equals to cols
Definition: DMMoFEM.hpp:897
virtual bool check_problem(const std::string name)=0
check if problem exist
PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb)
Set verbosity level.
Definition: DMMMoFEM.cpp:1307
int verbosity
verbosity
Definition: DMMoFEM.hpp:923
Managing BitRefLevels.
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set compute operator for KSP solver via sub-matrix and IS.
Definition: DMMMoFEM.cpp:545
DEPRECATED PetscErrorCode DMMoFEMResolveSharedEntities(DM dm, const char fe_name[])
Definition: DMMMoFEM.cpp:386
virtual MoFEMErrorCode delete_problem(const std::string name)=0
Delete problem.
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1178
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
boost::shared_ptr< std::map< std::string, std::pair< EntityType, EntityType > > > mapTypeRow
Definition: DMMoFEM.hpp:914
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
void cblas_daxpy(const int N, const double alpha, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_daxpy.c:11
PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is)
get sub problem is
Definition: DMMMoFEM.cpp:260
PetscBool isSubDM
Definition: DMMoFEM.hpp:905
Section manager is used to create indexes and sectionsFIXME: ISManager is not properly testsed by ato...
Definition: ISManager.hpp:36
PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, const char fe_name[])
Resolve shared entities.
Definition: DMMMoFEM.cpp:376
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMMoFEM.cpp:434
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problemIt if true is assumed that matrix has the same indexing on rows and columns....
Definition: DMMMoFEM.cpp:367
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:27
PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm)
Definition: DMMMoFEM.cpp:233
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[], EntityType lo_type=MBVERTEX, EntityType hi_type=MBMAXTYPE)
Definition: DMMMoFEM.cpp:210
MoFEMErrorCode buildCompsedProblem(const std::string out_name, const std::vector< std::string > add_row_problems, const std::vector< std::string > add_col_problems, const bool square_matrix=true, int verb=1)
build composite problem
PetscErrorCode DMMoFEMAddElement(DM dm, const char fe_name[])
add element to dm
Definition: DMMMoFEM.cpp:414
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition: SnesCtx.cpp:121
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string &name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Interface * mField_ptr
MoFEM interface.
Definition: DMMoFEM.hpp:891
base class for all interface classes
PetscBool isPartitioned
true if read mesh is on parts
Definition: DMMoFEM.hpp:896
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition: TsCtx.cpp:241
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMMoFEM.cpp:1084
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:544
PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on colThis create block on col with DOFs from problem of given name.
Definition: DMMMoFEM.cpp:296
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.
MoFEMErrorCode query_interface(const MOFEMuuid &uuid, UnknownInterface **iface) const
Definition: DMMMoFEM.cpp:41
DofIdx getNbLocalDofsRow() const
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:190
const Problem * problemPtr
pinter to problem data structure
Definition: DMMoFEM.hpp:902
std::vector< std::string > colCompPrb
Definition: DMMoFEM.hpp:912
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
Data structure to exchange data between mofem and User Loop Methods on entities.It allows to exchange...
virtual ~DMCtx()
Definition: DMMMoFEM.cpp:36
keeps basic data about problemThis is low level structure with information about problem,...
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1215
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMMoFEM.cpp:1018
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Definition: DMMMoFEM.cpp:586
DofIdx getNbGhostDofsRow() const
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition: DMMMoFEM.cpp:761
RowColData
RowColData.
Definition: definitions.h:186
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:496
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
Discrete manager interface for MoFEM.
PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm)
Get if this DM is composite DM.
Definition: DMMMoFEM.cpp:316
PetscBool destroyProblem
If true destroy problem with DM.
Definition: DMMoFEM.hpp:918
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
Definition: DMMMoFEM.cpp:89
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVectorsets the routine to create a global vector associated with the shell DM.
Definition: DMMMoFEM.cpp:991
PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, const char fe_name[], PetscLayout *layout)
Get finite elements layout in the problem.
Definition: DMMMoFEM.cpp:390
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
Definition: DMMoFEM.hpp:928
PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side function
Definition: DMMMoFEM.cpp:800
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition: DMMMoFEM.cpp:166
PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on rowThis create block on row with DOFs from problem of given name.
Definition: DMMMoFEM.cpp:278
PetscErrorCode KspRhs(KSP ksp, Vec f, void *ctx)
Run over elements in the lists.
Definition: KspCtx.cpp:23
std::vector< std::string > rowFields
Definition: DMMoFEM.hpp:906
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[], EntityType lo_type=MBVERTEX, EntityType hi_type=MBMAXTYPE)
Definition: DMMMoFEM.cpp:187
PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr< MoFEM::TsCtx > &ts_ctx)
Set MoFEM::TsCtx data structureIt take over pointer, do not delete it, DM will destroy pointer when i...
Definition: DMMMoFEM.cpp:983
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:17
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
Definition: DMMMoFEM.cpp:1131
MoFEMErrorCode partitionGhostDofsOnDistributedMesh(const std::string name, int verb=VERBOSE)
determine ghost nodes on distributed meshes
virtual MoFEMErrorCode get_problem(const std::string &problem_name, const Problem **problem_ptr) const =0
Get problem database (data structure)
PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem)
get squared problemIt if true is assumed that matrix has the same indexing on rows and columns....
Definition: DMMMoFEM.cpp:404
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition: DMMMoFEM.cpp:918
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:66
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition: DMMMoFEM.cpp:446
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
Definition: DMMoFEM.hpp:926
PetscBool isProblemBuild
True if problem is build.
Definition: DMMoFEM.hpp:892
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:27
const Problem * problemMainOfSubPtr
pinter to main problem to sub-problem
Definition: DMMoFEM.hpp:908
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:53
PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side jacobian
Definition: DMMMoFEM.cpp:840
PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem)
Definition: DMMMoFEM.cpp:358
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
Definition: DMMoFEM.hpp:927
PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx)
get MoFEM::KspCtx data structure
Definition: DMMMoFEM.cpp:892
virtual MoFEMErrorCode modify_problem_unset_finite_element(const std::string &name_problem, const std::string &fe_name)=0
unset finite element from problem, this remove entities assigned to finite element to a particular pr...
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition: DMMMoFEM.cpp:667
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
PetscErrorCode DMMoFEMUnSetElement(DM dm, const char fe_name[])
unset element from dm
Definition: DMMMoFEM.cpp:424
PetscErrorCode KspMat(KSP ksp, Mat A, Mat B, void *ctx)
Run over elenents in the list.
Definition: KspCtx.cpp:71
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMMoFEM.cpp:457
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:36
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1186
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMMoFEM.cpp:467
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMMoFEM.cpp:109
PetscErrorCode DMCreate_MoFEM(DM dm)
Create dm data structure with MoFEM data structure.
Definition: DMMMoFEM.cpp:81
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition: DMMMoFEM.cpp:59
#define CHKERR
Inline error check.
Definition: definitions.h:596
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition: DMMMoFEM.cpp:1250
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition: DMMMoFEM.cpp:880
MoFEMErrorCode buildSubProblem(const std::string out_name, const std::vector< std::string > &fields_row, const std::vector< std::string > &fields_col, const std::string main_problem, const bool square_matrix=true, const map< std::string, std::pair< EntityType, EntityType >> *entityMapRow=nullptr, const map< std::string, std::pair< EntityType, EntityType >> *entityMapCol=nullptr, int verb=VERBOSE)
build sub problem
std::vector< std::string > colFields
Definition: DMMoFEM.hpp:907
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode buildProblemOnDistributedMesh(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures, assuming that mesh is distributed (collective)Mesh is distributed,...
BitRefLevel getBitRefLevel() const
PetscBool isCompDM
Definition: DMMoFEM.hpp:910
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
std::vector< std::string > rowCompPrb
Definition: DMMoFEM.hpp:911
Interface for linear (KSP) solver.
Definition: KspCtx.hpp:27
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition: DMMMoFEM.cpp:966
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[], IS *is)
get field is in the problem
Definition: DMMMoFEM.cpp:1296
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition: TsCtx.cpp:325
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMMoFEM.cpp:325
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elementsFunction which partition finite elements based on dofs partitioning....
PETSc Discrete Manager data structure.
Definition: DMMoFEM.hpp:886
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
Definition: DMMMoFEM.cpp:349
PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set TS implicit function evaluation function
Definition: DMMMoFEM.cpp:708
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMMoFEM.cpp:337
PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned)
Definition: DMMMoFEM.cpp:958
static const MOFEMuuid IDD_DMCTX
Definition: DMMoFEM.hpp:871
virtual MPI_Comm & get_comm() const =0
std::string getName() const
virtual MoFEMErrorCode get_fields(const Field_multiIndex **fields_ptr) const =0
Get fields multi-index from database.
PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[], MoFEM::DofMethod *method)
execute method for dofs on field in problem
Definition: DMMMoFEM.cpp:513
MoFEMErrorCode getNumberOfElementsByNameAndPart(MPI_Comm comm, const std::string name, PetscLayout *layout) const
Get number of finite elements by name on processors.
intrusive_ptr for managing petsc objects
Definition: AuxPETSc.hpp:124
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:477
boost::shared_ptr< std::map< std::string, std::pair< EntityType, EntityType > > > mapTypeCol
Definition: DMMoFEM.hpp:916
PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition: SnesCtx.cpp:23
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethodThis function set data about problem, adjacencies and other multi-indices in ...
PetscErrorCode DMMoFEMSetKspCtx(DM dm, boost::shared_ptr< MoFEM::KspCtx > &ksp_ctx)
set MoFEM::KspCtx data structure
Definition: DMMMoFEM.cpp:909
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMMoFEM.cpp:947
Data structure to exchange data between mofem and User Loop Methods.It allows to exchange data betwee...
Matrix manager is used to build and partition problems.
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1244
virtual MoFEMErrorCode modify_problem_mask_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
set dof mask ref level for problem
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobina in TS solver.
Definition: TsCtx.cpp:108
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition: DMMMoFEM.cpp:1062
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVectorsets the routine to create a local vector associated with the shell DM.
Definition: DMMMoFEM.cpp:1009
PetscErrorCode DMMoFEMSetSnesCtx(DM dm, boost::shared_ptr< MoFEM::SnesCtx > &snes_ctx)
Set MoFEM::SnesCtx data structure.
Definition: DMMMoFEM.cpp:935