v0.14.0
DMMoFEM.cpp
Go to the documentation of this file.
1 
2 
3 // #undef PETSC_VERSION_RELEASE
4 // #define PETSC_VERSION_RELEASE 1
5 
6 #if PETSC_VERSION_GE(3, 6, 0)
7 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
8 // #include <petsc/private/vecimpl.h> /*I "petscdm.h" I*/
9 #else
10 #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/
11 #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/
12 #endif
13 
14 #include <DMMoFEM.hpp>
15 
16 namespace MoFEM {
17 
19  : mField_ptr(PETSC_NULL), isProblemBuild(PETSC_FALSE),
20  isPartitioned(PETSC_FALSE), isSquareMatrix(PETSC_TRUE),
21  isSubDM(PETSC_FALSE), isCompDM(PETSC_FALSE), destroyProblem(PETSC_FALSE),
22  verbosity(VERBOSE), referenceNumber(0) {
23 
24  if (!LogManager::checkIfChannelExist("DMWORLD")) {
25  auto core_log = logging::core::get();
26  core_log->add_sink(
28  core_log->add_sink(
30  core_log->add_sink(
32  LogManager::setLog("DMWORLD");
33  LogManager::setLog("DMSYNC");
34  LogManager::setLog("DMSELF");
35  MOFEM_LOG_TAG("DMWORLD", "DM");
36  MOFEM_LOG_TAG("DMSYNC", "DM");
37  MOFEM_LOG_TAG("DMSELF", "DM");
38  }
39 }
40 
41 MoFEMErrorCode DMCtx::query_interface(boost::typeindex::type_index type_index,
42  UnknownInterface **iface) const {
43  *iface = const_cast<DMCtx *>(this);
44  return 0;
45 }
46 
47 PetscErrorCode DMRegister_MoFEM(const char sname[]) {
49  CHKERR DMRegister(sname, DMCreate_MoFEM);
51 }
52 
53 PetscErrorCode DMSetOperators_MoFEM(DM dm) {
54  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
56 
57  dm->ops->createglobalvector = DMCreateGlobalVector_MoFEM;
58  dm->ops->createlocalvector = DMCreateLocalVector_MoFEM;
59  dm->ops->creatematrix = DMCreateMatrix_MoFEM;
60  dm->ops->setup = DMSetUp_MoFEM;
61  dm->ops->destroy = DMDestroy_MoFEM;
62  dm->ops->setfromoptions = DMSetFromOptions_MoFEM;
63  dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_MoFEM;
64  dm->ops->globaltolocalend = DMGlobalToLocalEnd_MoFEM;
65  dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_MoFEM;
66  dm->ops->localtoglobalend = DMLocalToGlobalEnd_MoFEM;
67  dm->ops->createfieldis = DMCreateFieldIS_MoFEM;
68 
69  // Default matrix type
70  CHKERR DMSetMatType(dm, MATMPIAIJ);
71 
73 }
74 
75 PetscErrorCode DMCreate_MoFEM(DM dm) {
76  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
78  dm->data = new DMCtx();
81 }
82 
83 PetscErrorCode DMDestroy_MoFEM(DM dm) {
84  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
85  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
87 
88  MPI_Comm comm;
89  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
90 
91  int result;
92  MPI_Comm_compare(comm, PETSC_COMM_SELF, &result);
93  if (result == MPI_IDENT)
94  MOFEM_LOG("DMSELF", Sev::noisy)
95  << "MoFEM DM destroy for problem " << dm_field->problemName
96  << " referenceNumber " << dm_field->referenceNumber;
97  else
98  MOFEM_LOG("DMWORLD", Sev::noisy)
99  << "MoFEM DM destroy for problem " << dm_field->problemName
100  << " referenceNumber " << dm_field->referenceNumber;
101 
102  if (dm_field->referenceNumber == 0) {
103  if (dm_field->destroyProblem) {
104 
105  if (dm_field->mField_ptr->check_problem(dm_field->problemName)) {
106  dm_field->mField_ptr->delete_problem(dm_field->problemName);
107  } // else problem has to be deleted by the user
108  }
109 
110  delete static_cast<DMCtx *>(dm->data);
111 
112  } else
113  --dm_field->referenceNumber;
114 
116 }
117 
118 PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr,
119  const char problem_name[],
120  const MoFEM::BitRefLevel bit_level,
121  const MoFEM::BitRefLevel bit_mask) {
123 
124  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
125  if (!dm->data) {
126  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
127  "data structure for MoFEM not yet created");
128  }
129  if (!m_field_ptr) {
130  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
131  "DM function not implemented into MoFEM");
132  }
133  dm_field->mField_ptr = m_field_ptr;
134  dm_field->problemName = problem_name;
135  if (!m_field_ptr->check_problem(dm_field->problemName)) {
136  // problem is not defined, declare problem internally set bool to
137  // destroyProblem problem with DM
138  dm_field->destroyProblem = PETSC_TRUE;
139  CHKERR dm_field->mField_ptr->add_problem(dm_field->problemName, MF_EXCL,
140  dm_field->verbosity);
141  } else {
142  dm_field->destroyProblem = PETSC_FALSE;
143  }
145  dm_field->problemName, bit_level);
147  dm_field->problemName, bit_mask);
148  dm_field->kspCtx =
149  boost::shared_ptr<KspCtx>(new KspCtx(*m_field_ptr, problem_name));
150  dm_field->snesCtx =
151  boost::shared_ptr<SnesCtx>(new SnesCtx(*m_field_ptr, problem_name));
152  dm_field->tsCtx =
153  boost::shared_ptr<TsCtx>(new TsCtx(*m_field_ptr, problem_name));
154 
155  MPI_Comm comm;
156  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
157  int result = 0;
158  MPI_Comm_compare(comm, m_field_ptr->get_comm(), &result);
159  if (result > MPI_CONGRUENT) {
160  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
161  "MoFEM and DM using different communicators");
162  }
163  MPI_Comm_size(comm, &dm_field->sIze);
164  MPI_Comm_rank(comm, &dm_field->rAnk);
165 
166  // problem structure
167  CHKERR dm_field->mField_ptr->get_problem(dm_field->problemName,
168  &dm_field->problemPtr);
169 
170  MPI_Comm_compare(comm, PETSC_COMM_SELF, &result);
171  if (result == MPI_IDENT) {
172  MOFEM_LOG("DMSELF", Sev::noisy)
173  << "MoFEM DM created for problem " << dm_field->problemName;
174  MOFEM_LOG("DMSELF", Sev::noisy) << *dm_field->problemPtr;
175  } else {
176  MOFEM_LOG("DMWORLD", Sev::noisy)
177  << "MoFEM DM created for problem " << dm_field->problemName;
178  MOFEM_LOG("DMWORLD", Sev::noisy) << *dm_field->problemPtr;
179  }
180 
182 }
183 
184 PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate) {
186 
187  if (!dm->data)
188  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
189  "data structure for MoFEM not yet created");
190 
191  if (static_cast<DMCtx *>(dm_duplicate->data)->referenceNumber == 0)
192  delete static_cast<DMCtx *>(dm_duplicate->data);
193 
194  dm_duplicate->data = dm->data;
195  ++(static_cast<DMCtx *>(dm_duplicate->data)->referenceNumber);
196 
198 }
199 
200 PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap) {
202  if (!dm->data)
203  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
204  "data structure for MoFEM not yet created on dm");
205  if (!dm_swap->data)
206  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
207  "data structure for MoFEM not yet created on swap dm");
208 
209  auto *dm_field = static_cast<DMCtx *>(dm->data);
210  auto *dm_field_swap = static_cast<DMCtx *>(dm_swap->data);
211 
212  auto tmp_field = dm_field;
213  dm_field = dm_field_swap;
214  dm_field_swap = tmp_field;
215 
217 }
218 
219 PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[]) {
221 
222  if (!dm->data) {
223  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
224  "data structure for MoFEM not yet created");
225  }
226  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
227 
228  CHKERR DMMoFEMCreateMoFEM(subdm, dm_field->mField_ptr, problem_name,
229  dm_field->problemPtr->getBitRefLevel(),
230  dm_field->problemPtr->getBitRefLevelMask());
231 
232  DMCtx *subdm_field = (DMCtx *)subdm->data;
233  subdm_field->isSubDM = PETSC_TRUE;
234  subdm_field->problemMainOfSubPtr = dm_field->problemPtr;
235  subdm_field->isPartitioned = dm_field->isPartitioned;
236  subdm_field->isSquareMatrix = PETSC_FALSE;
237  subdm->ops->setup = DMSubDMSetUp_MoFEM;
238 
240 }
241 
242 PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[]) {
243  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
245  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
246  if (!dm->data) {
247  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
248  "data structure for MoFEM not yet created");
249  }
250  if (!dm_field->isSubDM) {
251  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
252  }
253  dm_field->rowFields.push_back(field_name);
254  dm_field->mapTypeRow.erase(field_name);
256 }
257 
258 PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name) {
259  return DMMoFEMAddSubFieldRow(dm, field_name.c_str());
260 }
261 
262 PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[],
263  boost::shared_ptr<Range> r_ptr) {
264  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
266  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
267  if (!dm->data) {
268  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
269  "data structure for MoFEM not yet created");
270  }
271  if (!dm_field->isSubDM) {
272  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
273  }
274  dm_field->rowFields.push_back(field_name);
275  dm_field->mapTypeRow[field_name] = r_ptr;
277 }
278 
279 PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name,
280  boost::shared_ptr<Range> r_ptr) {
281  return DMMoFEMAddSubFieldRow(dm, field_name.c_str(), r_ptr);
282 }
283 
284 PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[]) {
285  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
287  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
288  if (!dm->data) {
289  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
290  "data structure for MoFEM not yet created");
291  }
292  if (!dm_field->isSubDM) {
293  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
294  }
295  dm_field->colFields.push_back(field_name);
296  dm_field->mapTypeCol.erase(field_name);
298 }
299 
300 PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name) {
301  return DMMoFEMAddSubFieldCol(dm, field_name.c_str());
302 }
303 
304 PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[],
305  boost::shared_ptr<Range> r_ptr) {
306  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
308  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
309  if (!dm->data) {
310  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
311  "data structure for MoFEM not yet created");
312  }
313  if (!dm_field->isSubDM) {
314  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "this is not sub-dm");
315  }
316  dm_field->colFields.push_back(field_name);
317  dm_field->mapTypeCol[field_name] = r_ptr;
319 }
320 
321 PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name,
322  boost::shared_ptr<Range> r_ptr) {
323  return DMMoFEMAddSubFieldCol(dm, field_name.c_str(), r_ptr);
324 }
325 
326 PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm) {
328  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
330  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
331  *is_sub_dm = dm_field->isSubDM;
333 }
334 
335 PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is) {
337  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
339  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
340  if (dm_field->isSubDM != PETSC_TRUE) {
341  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
342  "This DM is not created as a SubDM");
343  }
344  if (dm_field->isProblemBuild != PETSC_TRUE) {
345  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Problem is not build");
346  }
347  boost::shared_ptr<Problem::SubProblemData> sub_data =
348  dm_field->problemPtr->getSubData();
349  CHKERR sub_data->getRowIs(is);
351 }
352 
353 PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is) {
355  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
357  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
358  if (dm_field->isSubDM != PETSC_TRUE) {
359  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
360  "This DM is not created as a SubDM");
361  }
362  if (dm_field->isProblemBuild != PETSC_TRUE) {
363  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Problem is not build");
364  }
365  boost::shared_ptr<Problem::SubProblemData> sub_data =
366  dm_field->problemPtr->getSubData();
367  CHKERR sub_data->getColIs(is);
369 }
370 
371 PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[]) {
372  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
374  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
375  if (!dm->data) {
376  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
377  "data structure for MoFEM not yet created");
378  }
379  if (!dm_field->isCompDM) {
380  dm_field->isCompDM = PETSC_TRUE;
381  }
382  dm_field->rowCompPrb.push_back(prb_name);
383  if (dm_field->isSquareMatrix) {
384  dm_field->colCompPrb.push_back(prb_name);
385  }
387 }
388 
389 PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[]) {
390  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
392  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
393  if (!dm->data) {
394  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
395  "data structure for MoFEM not yet created");
396  }
397  if (!dm_field->isCompDM) {
398  dm_field->isCompDM = PETSC_TRUE;
399  }
400  if (dm_field->isSquareMatrix) {
401  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
402  "No need to add problem on column when problem block structurally "
403  "symmetric");
404  }
405  dm_field->colCompPrb.push_back(prb_name);
407 }
408 
409 PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm) {
411  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
414  *is_comp_dm = dm_field->isCompDM;
416 }
417 
418 PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr) {
419  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
421  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
422  if (!dm->data) {
423  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
424  "data structure for MoFEM not yet created");
425  }
426  *m_field_ptr = dm_field->mField_ptr;
428 }
429 
430 PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr) {
431  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
433  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
434  if (!dm->data) {
435  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
436  "data structure for MoFEM not yet created");
437  }
438  *problem_ptr = dm_field->problemPtr;
440 }
441 
442 PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem) {
444  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
446  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
447  dm_field->destroyProblem = destroy_problem;
449 }
450 
451 PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem) {
453  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
455  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
456  *destroy_problem = dm_field->destroyProblem;
458 }
459 
460 PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem) {
461  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
463  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
464  dm_field->isSquareMatrix = square_problem;
466 }
467 
468 PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name) {
469  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
471  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
473  ->resolveSharedFiniteElements(dm_field->problemPtr, fe_name);
475 }
476 
477 PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name,
478  PetscLayout *layout) {
480  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
482  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
483 
484  MPI_Comm comm;
485  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
486  CHKERR dm_field->problemPtr->getNumberOfElementsByNameAndPart(comm, fe_name,
487  layout);
489 }
490 
491 PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem) {
494  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
496  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
497  *square_problem = dm_field->isSquareMatrix;
499 }
500 
501 PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name) {
502  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
504  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
506  dm_field->problemName, fe_name);
508 }
509 
510 PetscErrorCode DMMoFEMAddElement(DM dm, std::vector<std::string> fe_name) {
512  for (auto fe : fe_name) {
513  CHKERR DMMoFEMAddElement(dm, fe);
514  }
516 }
517 
518 PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name) {
519  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
521  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
523  dm_field->problemName, fe_name);
525 }
526 
527 PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode,
528  ScatterMode scatter_mode) {
530  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
532  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
533  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setLocalGhostVector(
534  dm_field->problemPtr, COL, l, mode, scatter_mode);
535  CHKERRG(ierr);
537 }
538 
539 PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode,
540  ScatterMode scatter_mode) {
541  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
543  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
544  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setGlobalGhostVector(
545  dm_field->problemPtr, COL, g, mode, scatter_mode);
546  CHKERRG(ierr);
548 }
549 
550 PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
551  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
553  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
555  dm_field->problemPtr, *method);
556  CHKERRG(ierr);
558 }
559 
560 PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
561  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
563  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
565  dm_field->problemPtr, *method);
566  CHKERRG(ierr);
568 }
569 
570 PetscErrorCode
571 DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[],
572  MoFEM::FEMethod *method, int low_rank,
573  int up_rank, CacheTupleWeakPtr cache_ptr) {
575  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
576  ierr = dm_field->mField_ptr->loop_finite_elements(
577  dm_field->problemPtr, fe_name, *method, low_rank, up_rank, nullptr,
578  MF_EXIST, cache_ptr);
579  CHKERRG(ierr);
581 }
582 
584  DM dm, const std::string fe_name, boost::shared_ptr<MoFEM::FEMethod> method,
585  int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr) {
586  return DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name.c_str(), method.get(),
587  low_rank, up_rank, cache_ptr);
588 }
589 
590 PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[],
591  MoFEM::FEMethod *method,
592  CacheTupleWeakPtr cache_ptr) {
593  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
595  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
597  dm, fe_name, method, dm_field->rAnk, dm_field->rAnk, cache_ptr);
598  CHKERRG(ierr);
600 }
601 
602 PetscErrorCode
603 DMoFEMLoopFiniteElements(DM dm, const std::string fe_name,
604  boost::shared_ptr<MoFEM::FEMethod> method,
605  CacheTupleWeakPtr cache_ptr) {
606  return DMoFEMLoopFiniteElements(dm, fe_name.c_str(), method.get(), cache_ptr);
607 }
608 
609 PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[],
610  MoFEM::DofMethod *method) {
611  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
613  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
614  ierr =
615  dm_field->mField_ptr->loop_dofs(dm_field->problemPtr, field_name, COL,
616  *method, dm_field->rAnk, dm_field->rAnk);
617  CHKERRG(ierr);
619 }
620 
621 template <class S, class T0, class T1, class T2>
622 static PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, S fe_name, T0 method,
623  T1 pre_only, T2 post_only) {
624  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
626  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
627  if (pre_only) {
628  dm_field->kspCtx->getPreProcComputeRhs().push_back(pre_only);
629  }
630  if (method) {
631  dm_field->kspCtx->getComputeRhs().push_back(
632  PairNameFEMethodPtr(fe_name, method));
633  }
634  if (post_only) {
635  dm_field->kspCtx->getPostProcComputeRhs().push_back(post_only);
636  }
637  CHKERR DMKSPSetComputeRHS(dm, KspRhs, dm_field->kspCtx.get());
639 }
640 
641 PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[],
642  MoFEM::FEMethod *method,
643  MoFEM::BasicMethod *pre_only,
644  MoFEM::BasicMethod *post_only) {
645  return DMMoFEMKSPSetComputeRHS<const char *, MoFEM::FEMethod *,
647  dm, fe_name, method, pre_only, post_only);
648 }
649 
650 PetscErrorCode
651 DMMoFEMKSPSetComputeRHS(DM dm, const std::string fe_name,
652  boost::shared_ptr<MoFEM::FEMethod> method,
653  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
654  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
655  return DMMoFEMKSPSetComputeRHS<const std::string,
656  boost::shared_ptr<MoFEM::FEMethod>,
657  boost::shared_ptr<MoFEM::BasicMethod>,
658  boost::shared_ptr<MoFEM::BasicMethod>>(
659  dm, fe_name, method, pre_only, post_only);
660 }
661 
662 template <class S, class T0, class T1, class T2>
663 static PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, S fe_name, T0 method,
664  T1 pre_only, T2 post_only) {
665  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
667  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
668  if (pre_only) {
669  dm_field->kspCtx->getPreProcSetOperators().push_back(pre_only);
670  }
671  if (method) {
672  dm_field->kspCtx->getSetOperators().push_back(
673  PairNameFEMethodPtr(fe_name, method));
674  }
675  if (post_only) {
676  dm_field->kspCtx->getPostProcSetOperators().push_back(post_only);
677  }
678  CHKERR DMKSPSetComputeOperators(dm, KspMat, dm_field->kspCtx.get());
680 }
681 
682 PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[],
683  MoFEM::FEMethod *method,
684  MoFEM::BasicMethod *pre_only,
685  MoFEM::BasicMethod *post_only) {
686  return DMMoFEMKSPSetComputeOperators<const char *, MoFEM::FEMethod *,
689  dm, fe_name, method, pre_only, post_only);
690 }
691 
692 PetscErrorCode
693 DMMoFEMKSPSetComputeOperators(DM dm, const std::string fe_name,
694  boost::shared_ptr<MoFEM::FEMethod> method,
695  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
696  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
697  return DMMoFEMKSPSetComputeOperators<const std::string,
698  boost::shared_ptr<MoFEM::FEMethod>>(
699  dm, fe_name, method, pre_only, post_only);
700 }
701 
702 template <class S, class T0, class T1, class T2>
703 static PetscErrorCode DMMoFEMSNESSetFunction(DM dm, S fe_name, T0 method,
704  T1 pre_only, T2 post_only) {
705  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
707  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
708  if (pre_only) {
709  dm_field->snesCtx->getPreProcComputeRhs().push_back(pre_only);
710  }
711  if (method) {
712  dm_field->snesCtx->getComputeRhs().push_back(
713  PairNameFEMethodPtr(fe_name, method));
714  }
715  if (post_only) {
716  dm_field->snesCtx->getPostProcComputeRhs().push_back(post_only);
717  }
718  CHKERR DMSNESSetFunction(dm, SnesRhs, dm_field->snesCtx.get());
720 }
721 
722 PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[],
723  MoFEM::FEMethod *method,
724  MoFEM::BasicMethod *pre_only,
725  MoFEM::BasicMethod *post_only) {
726  return DMMoFEMSNESSetFunction<const char *, MoFEM::FEMethod *,
728  dm, fe_name, method, pre_only, post_only);
729 }
730 
731 PetscErrorCode
732 DMMoFEMSNESSetFunction(DM dm, const std::string fe_name,
733  boost::shared_ptr<MoFEM::FEMethod> method,
734  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
735  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
736  return DMMoFEMSNESSetFunction<const std::string,
737  boost::shared_ptr<MoFEM::FEMethod>,
738  boost::shared_ptr<MoFEM::BasicMethod>,
739  boost::shared_ptr<MoFEM::BasicMethod>>(
740  dm, fe_name, method, pre_only, post_only);
741 }
742 
743 template <class S, class T0, class T1, class T2>
744 static PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, S fe_name, T0 method,
745  T1 pre_only, T2 post_only) {
746  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
748  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
749  if (pre_only) {
750  dm_field->snesCtx->getPreProcSetOperators().push_back(pre_only);
751  }
752  if (method) {
753  dm_field->snesCtx->getSetOperators().push_back(
754  PairNameFEMethodPtr(fe_name, method));
755  }
756  if (post_only) {
757  dm_field->snesCtx->getPostProcSetOperators().push_back(post_only);
758  }
759  CHKERR DMSNESSetJacobian(dm, SnesMat, dm_field->snesCtx.get());
761 }
762 
763 PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[],
764  MoFEM::FEMethod *method,
765  MoFEM::BasicMethod *pre_only,
766  MoFEM::BasicMethod *post_only) {
767  return DMMoFEMSNESSetJacobian<const char *, MoFEM::FEMethod *,
769  dm, fe_name, method, pre_only, post_only);
770 }
771 
772 PetscErrorCode
773 DMMoFEMSNESSetJacobian(DM dm, const std::string fe_name,
774  boost::shared_ptr<MoFEM::FEMethod> method,
775  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
776  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
777  return DMMoFEMSNESSetJacobian<const std::string,
778  boost::shared_ptr<MoFEM::FEMethod>,
779  boost::shared_ptr<MoFEM::BasicMethod>,
780  boost::shared_ptr<MoFEM::BasicMethod>>(
781  dm, fe_name, method, pre_only, post_only);
782 }
783 
784 template <class S, class T0, class T1, class T2>
785 static PetscErrorCode DMMoFEMTSSetIFunction(DM dm, S fe_name, T0 method,
786  T1 pre_only, T2 post_only) {
787  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
789  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
790  if (pre_only) {
791  dm_field->tsCtx->getPreProcessIFunction().push_back(pre_only);
792  }
793  if (method) {
794  dm_field->tsCtx->getLoopsIFunction().push_back(
795  PairNameFEMethodPtr(fe_name, method));
796  }
797  if (post_only) {
798  dm_field->tsCtx->getPostProcessIFunction().push_back(post_only);
799  }
800  CHKERR DMTSSetIFunction(dm, TsSetIFunction, dm_field->tsCtx.get());
802 }
803 
804 PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[],
805  MoFEM::FEMethod *method,
806  MoFEM::BasicMethod *pre_only,
807  MoFEM::BasicMethod *post_only) {
808  return DMMoFEMTSSetIFunction<const char *, MoFEM::FEMethod *,
810  dm, fe_name, method, pre_only, post_only);
812 }
813 
814 PetscErrorCode
815 DMMoFEMTSSetIFunction(DM dm, const std::string fe_name,
816  boost::shared_ptr<MoFEM::FEMethod> method,
817  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
818  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
819  return DMMoFEMTSSetIFunction<const std::string,
820  boost::shared_ptr<MoFEM::FEMethod>,
821  boost::shared_ptr<MoFEM::BasicMethod>,
822  boost::shared_ptr<MoFEM::BasicMethod>>(
823  dm, fe_name, method, pre_only, post_only);
825 }
826 
827 template <class S, class T0, class T1, class T2>
828 static PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, S fe_name, T0 method,
829  T1 pre_only, T2 post_only) {
830  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
832  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
833  if (pre_only) {
834  dm_field->tsCtx->getPreProcessIJacobian().push_back(pre_only);
835  }
836  if (method) {
837  dm_field->tsCtx->getLoopsIJacobian().push_back(
838  PairNameFEMethodPtr(fe_name, method));
839  }
840  if (post_only) {
841  dm_field->tsCtx->getPostProcessIJacobian().push_back(post_only);
842  }
843  CHKERR DMTSSetIJacobian(dm, TsSetIJacobian, dm_field->tsCtx.get());
845 }
846 
847 PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const char fe_name[],
848  MoFEM::FEMethod *method,
849  MoFEM::BasicMethod *pre_only,
850  MoFEM::BasicMethod *post_only) {
851  return DMMoFEMTSSetIJacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
852  MoFEM::BasicMethod *>(dm, fe_name, method,
853  pre_only, post_only);
854 }
855 
856 PetscErrorCode
857 DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name,
858  boost::shared_ptr<MoFEM::FEMethod> method,
859  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
860  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
861  return DMMoFEMTSSetIJacobian<const std::string,
862  boost::shared_ptr<MoFEM::FEMethod>,
863  boost::shared_ptr<MoFEM::BasicMethod>,
864  boost::shared_ptr<MoFEM::BasicMethod>>(
865  dm, fe_name, method, pre_only, post_only);
866 }
867 
868 template <class S, class T0, class T1, class T2>
869 static PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, S fe_name, T0 method,
870  T1 pre_only, T2 post_only) {
871  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
873  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
874  if (pre_only)
875  dm_field->tsCtx->getPreProcessRHSFunction().push_back(pre_only);
876  if (method)
877  dm_field->tsCtx->getLoopsRHSFunction().push_back(
878  PairNameFEMethodPtr(fe_name, method));
879  if (post_only)
880  dm_field->tsCtx->getPostProcessRHSFunction().push_back(post_only);
881  CHKERR DMTSSetRHSFunction(dm, TsSetRHSFunction, dm_field->tsCtx.get());
883 }
884 
885 PetscErrorCode
886 DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name,
887  boost::shared_ptr<MoFEM::FEMethod> method,
888  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
889  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
890  return DMMoFEMTSSetRHSFunction<const std::string,
891  boost::shared_ptr<MoFEM::FEMethod>,
892  boost::shared_ptr<MoFEM::BasicMethod>,
893  boost::shared_ptr<MoFEM::BasicMethod>>(
894  dm, fe_name, method, pre_only, post_only);
896 }
897 
898 PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const char fe_name[],
899  MoFEM::FEMethod *method,
900  MoFEM::BasicMethod *pre_only,
901  MoFEM::BasicMethod *post_only) {
902  return DMMoFEMTSSetRHSFunction<const char *, MoFEM::FEMethod *,
904  dm, fe_name, method, pre_only, post_only);
906 }
907 
908 template <class S, class T0, class T1, class T2>
909 static PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, S fe_name, T0 method,
910  T1 pre_only, T2 post_only) {
911  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
913  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
914  if (pre_only)
915  dm_field->tsCtx->getPreProcessRHSFunction().push_back(pre_only);
916  if (method)
917  dm_field->tsCtx->getLoopsRHSFunction().push_back(
918  PairNameFEMethodPtr(fe_name, method));
919  if (post_only)
920  dm_field->tsCtx->getPostProcessRHSFunction().push_back(post_only);
921  CHKERR DMTSSetRHSJacobian(dm, TsSetRHSJacobian, dm_field->tsCtx.get());
923 }
924 
925 PetscErrorCode
926 DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name,
927  boost::shared_ptr<MoFEM::FEMethod> method,
928  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
929  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
930  return DMMoFEMTSSetRHSJacobian<const std::string,
931  boost::shared_ptr<MoFEM::FEMethod>,
932  boost::shared_ptr<MoFEM::BasicMethod>,
933  boost::shared_ptr<MoFEM::BasicMethod>>(
934  dm, fe_name, method, pre_only, post_only);
936 }
937 
938 PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const char fe_name[],
939  MoFEM::FEMethod *method,
940  MoFEM::BasicMethod *pre_only,
941  MoFEM::BasicMethod *post_only) {
942  return DMMoFEMTSSetRHSJacobian<const char *, MoFEM::FEMethod *,
944  dm, fe_name, method, pre_only, post_only);
946 }
947 
948 template <class S, class T0, class T1, class T2>
949 static PetscErrorCode DMMoFEMTSSetI2Function(DM dm, S fe_name, T0 method,
950  T1 pre_only, T2 post_only) {
951  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
953  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
954  if (pre_only) {
955  dm_field->tsCtx->getPreProcessIFunction().push_back(pre_only);
956  }
957  if (method) {
958  dm_field->tsCtx->getLoopsIFunction().push_back(
959  PairNameFEMethodPtr(fe_name, method));
960  }
961  if (post_only) {
962  dm_field->tsCtx->getPostProcessIFunction().push_back(post_only);
963  }
964  CHKERR DMTSSetI2Function(dm, TsSetI2Function, dm_field->tsCtx.get());
966 }
967 
968 PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const char fe_name[],
969  MoFEM::FEMethod *method,
970  MoFEM::BasicMethod *pre_only,
971  MoFEM::BasicMethod *post_only) {
972  return DMMoFEMTSSetI2Function<const char *, MoFEM::FEMethod *,
974  dm, fe_name, method, pre_only, post_only);
976 }
977 
978 PetscErrorCode
979 DMMoFEMTSSetI2Function(DM dm, const std::string fe_name,
980  boost::shared_ptr<MoFEM::FEMethod> method,
981  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
982  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
983  return DMMoFEMTSSetI2Function<const std::string,
984  boost::shared_ptr<MoFEM::FEMethod>,
985  boost::shared_ptr<MoFEM::BasicMethod>,
986  boost::shared_ptr<MoFEM::BasicMethod>>(
987  dm, fe_name, method, pre_only, post_only);
989 }
990 
991 template <class S, class T0, class T1, class T2>
992 static PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, S fe_name, T0 method,
993  T1 pre_only, T2 post_only) {
994  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
996  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
997  if (pre_only) {
998  dm_field->tsCtx->getPreProcessIJacobian().push_back(pre_only);
999  }
1000  if (method) {
1001  dm_field->tsCtx->getLoopsIJacobian().push_back(
1002  PairNameFEMethodPtr(fe_name, method));
1003  }
1004  if (post_only) {
1005  dm_field->tsCtx->getPostProcessIJacobian().push_back(post_only);
1006  }
1007  CHKERR DMTSSetI2Jacobian(dm, TsSetI2Jacobian, dm_field->tsCtx.get());
1009 }
1010 
1011 PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, const char fe_name[],
1012  MoFEM::FEMethod *method,
1013  MoFEM::BasicMethod *pre_only,
1014  MoFEM::BasicMethod *post_only) {
1015  return DMMoFEMTSSetI2Jacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
1016  MoFEM::BasicMethod *>(dm, fe_name, method,
1017  pre_only, post_only);
1018 }
1019 
1020 PetscErrorCode
1021 DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name,
1022  boost::shared_ptr<MoFEM::FEMethod> method,
1023  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
1024  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
1025  return DMMoFEMTSSetI2Jacobian<const std::string,
1026  boost::shared_ptr<MoFEM::FEMethod>,
1027  boost::shared_ptr<MoFEM::BasicMethod>,
1028  boost::shared_ptr<MoFEM::BasicMethod>>(
1029  dm, fe_name, method, pre_only, post_only);
1030 }
1031 
1032 template <class S, class T0, class T1, class T2>
1033 static PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, S fe_name, T0 method,
1034  T1 pre_only, T2 post_only) {
1035  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1037  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1038  if (pre_only)
1039  dm_field->tsCtx->getPreProcessMonitor().push_back(pre_only);
1040  if (method)
1041  dm_field->tsCtx->getLoopsMonitor().push_back(
1042  PairNameFEMethodPtr(fe_name, method));
1043  if (post_only)
1044  dm_field->tsCtx->getPostProcessMonitor().push_back(post_only);
1045  CHKERR TSMonitorSet(ts, TsMonitorSet, dm_field->tsCtx.get(), PETSC_NULL);
1047 }
1048 
1049 PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const char fe_name[],
1050  MoFEM::FEMethod *method,
1051  MoFEM::BasicMethod *pre_only,
1052  MoFEM::BasicMethod *post_only) {
1053  return DMMoFEMTSSetMonitor<const char *, MoFEM::FEMethod *,
1055  dm, ts, fe_name, method, pre_only, post_only);
1057 }
1058 
1059 PetscErrorCode
1060 DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name,
1061  boost::shared_ptr<MoFEM::FEMethod> method,
1062  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
1063  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
1064  return DMMoFEMTSSetMonitor<const std::string,
1065  boost::shared_ptr<MoFEM::FEMethod>,
1066  boost::shared_ptr<MoFEM::BasicMethod>,
1067  boost::shared_ptr<MoFEM::BasicMethod>>(
1068  dm, ts, fe_name, method, pre_only, post_only);
1070 }
1071 
1072 PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx) {
1073  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1075  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1076  *ksp_ctx = dm_field->kspCtx.get();
1078 }
1079 
1080 PetscErrorCode
1081 DMMoFEMGetKspCtx(DM dm, const boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx) {
1082  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1084  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1085  const_cast<boost::shared_ptr<MoFEM::KspCtx> &>(ksp_ctx) = dm_field->kspCtx;
1087 }
1088 
1089 PetscErrorCode DMMoFEMSetKspCtx(DM dm,
1090  boost::shared_ptr<MoFEM::KspCtx> ksp_ctx) {
1091  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1093  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1094  dm_field->kspCtx = ksp_ctx;
1096 }
1097 
1098 PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx) {
1099  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1101  DMCtx *dm_field = (DMCtx *)dm->data;
1102  *snes_ctx = dm_field->snesCtx.get();
1104 }
1105 
1106 PetscErrorCode
1107 DMMoFEMGetSnesCtx(DM dm, const boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx) {
1108  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1110  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1111  const_cast<boost::shared_ptr<MoFEM::SnesCtx> &>(snes_ctx) = dm_field->snesCtx;
1113 }
1114 
1115 PetscErrorCode DMMoFEMSetSnesCtx(DM dm,
1116  boost::shared_ptr<MoFEM::SnesCtx> snes_ctx) {
1117  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1119  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1120  dm_field->snesCtx = snes_ctx;
1122 }
1123 
1124 /** get if read mesh is partitioned
1125  * \ingroup dm
1126  */
1127 PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned) {
1128  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1130  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1131  dm_field->isPartitioned = is_partitioned;
1133 }
1134 
1135 /** get if read mesh is partitioned
1136  * \ingroup dm
1137  */
1138 PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned) {
1139  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1141  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1142  *is_partitioned = dm_field->isPartitioned;
1144 }
1145 
1146 PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx) {
1147  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1149  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1150  *ts_ctx = dm_field->tsCtx.get();
1152 }
1153 
1154 PetscErrorCode DMMoFEMGetTsCtx(DM dm,
1155  const boost::shared_ptr<MoFEM::TsCtx> &ts_ctx) {
1156  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1158  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1159  const_cast<boost::shared_ptr<MoFEM::TsCtx> &>(ts_ctx) = dm_field->tsCtx;
1161 }
1162 
1163 PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr<MoFEM::TsCtx> ts_ctx) {
1164  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1166  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1167  dm_field->tsCtx = ts_ctx;
1169 }
1170 
1171 PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g) {
1172  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1174  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1175  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1176  dm_field->problemName, COL, g);
1177  CHKERR VecSetDM(*g, dm);
1179 }
1180 
1181 PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, SmartPetscObj<Vec> &g_ptr) {
1182  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1184  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1185  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
1186  dm_field->problemName, COL, g_ptr);
1187  CHKERR VecSetDM(g_ptr, dm);
1189 }
1190 
1191 PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l) {
1192  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1194  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1195  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateSeq(
1196  dm_field->problemName, COL, l);
1197  CHKERR VecSetDM(*l, dm);
1199 }
1200 
1201 PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M) {
1202  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1204  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1205  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1207  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1208  M);
1209  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1211  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1212  M);
1213  } else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1215  ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1216  dm_field->problemName, M);
1217  } else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1219  ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1220  dm_field->problemName, M);
1221  } else {
1222  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1223  "Matrix type not implemented");
1224  }
1225  CHKERR MatSetDM(*M, dm);
1227 }
1228 
1229 PetscErrorCode DMCreateMatrix_MoFEM(DM dm, SmartPetscObj<Mat> &M) {
1230  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1232  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1233  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
1235  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
1236  M);
1237  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
1239  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
1240  M);
1241  } else if (strcmp(dm->mattype, MATAIJCUSPARSE) == 0) {
1243  ->createMPIAIJCUSPARSEWithArrays<PetscGlobalIdx_mi_tag>(
1244  dm_field->problemName, M);
1245  } else if (strcmp(dm->mattype, MATSEQAIJCUSPARSE) == 0) {
1247  ->createSeqAIJCUSPARSEWithArrays<PetscLocalIdx_mi_tag>(
1248  dm_field->problemName, M);
1249  } else {
1250  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1251  "Matrix type not implemented");
1252  }
1253  CHKERR MatSetDM(M, dm);
1255 }
1256 
1257 #if PETSC_VERSION_GE(3, 7, 0)
1258 PetscErrorCode DMSetFromOptions_MoFEM(PetscOptionItems *PetscOptionsObject,
1259  DM dm) {
1260 #elif PETSC_VERSION_GE(3, 5, 3)
1261 PetscErrorCode DMSetFromOptions_MoFEM(PetscOptions *PetscOptionsObject, DM dm) {
1262 #else
1263 PetscErrorCode DMSetFromOptions_MoFEM(DM dm) {
1264 #endif
1265 
1266  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1268  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1269 #if PETSC_VERSION_GE(3, 5, 3)
1270  ierr = PetscOptionsHead(PetscOptionsObject, "DMMoFEM Options");
1271  CHKERRG(ierr);
1272 #else
1273  ierr = PetscOptionsHead("DMMoFEM Options");
1274  CHKERRG(ierr);
1275 #endif
1276  ierr = PetscOptionsBool("-dm_is_partitioned",
1277  "set if mesh is partitioned (works which native MOAB "
1278  "file format, i.e. h5m",
1279  "DMSetUp", dm_field->isPartitioned,
1280  &dm_field->isPartitioned, NULL);
1281  CHKERRG(ierr);
1283 }
1284 
1285 PetscErrorCode DMSetUp_MoFEM(DM dm) {
1286  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1287  ProblemsManager *prb_mng_ptr;
1289  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1290  CHKERR dm_field->mField_ptr->getInterface(prb_mng_ptr);
1291 
1292  if (dm_field->isCompDM) {
1293  // It is composite probelm
1294  CHKERR prb_mng_ptr->buildComposedProblem(
1295  dm_field->problemName, dm_field->rowCompPrb, dm_field->colCompPrb,
1296  dm_field->isSquareMatrix == PETSC_TRUE, dm_field->verbosity);
1297  } else {
1298  if (dm_field->isPartitioned) {
1300  dm_field->problemName, dm_field->isSquareMatrix == PETSC_TRUE,
1301  dm_field->verbosity);
1302  } else {
1303  CHKERR prb_mng_ptr->buildProblem(dm_field->problemName,
1304  dm_field->isSquareMatrix == PETSC_TRUE,
1305  dm_field->verbosity);
1306  CHKERR prb_mng_ptr->partitionProblem(dm_field->problemName,
1307  dm_field->verbosity);
1308  }
1309  }
1310 
1311  // Partition finite elements
1312  if (dm_field->isPartitioned) {
1313  CHKERR prb_mng_ptr->partitionFiniteElements(
1314  dm_field->problemName, true, 0, dm_field->sIze, dm_field->verbosity);
1316  dm_field->problemName, dm_field->verbosity);
1317  } else {
1318  // partition finite elemnets
1319  CHKERR prb_mng_ptr->partitionFiniteElements(dm_field->problemName, false,
1320  -1, -1, dm_field->verbosity);
1321  // Get ghost DOFs
1322  CHKERR prb_mng_ptr->partitionGhostDofs(dm_field->problemName,
1323  dm_field->verbosity);
1324  }
1325 
1326  // Set flag that problem is build and partitioned
1327  dm_field->isProblemBuild = PETSC_TRUE;
1328 
1330 }
1331 
1332 PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm) {
1333  PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
1334  ProblemsManager *prb_mng_ptr;
1336 
1337  DMCtx *subdm_field = static_cast<DMCtx *>(subdm->data);
1338 
1339  // build sub dm problem
1340  CHKERR subdm_field->mField_ptr->getInterface(prb_mng_ptr);
1341 
1342  map<std::string, boost::shared_ptr<Range>> *entity_map_row = nullptr;
1343  map<std::string, boost::shared_ptr<Range>> *entity_map_col = nullptr;
1344 
1345  if (subdm_field->mapTypeRow.size())
1346  entity_map_row = &subdm_field->mapTypeRow;
1347  if (subdm_field->mapTypeCol.size())
1348  entity_map_col = &subdm_field->mapTypeCol;
1349 
1350  CHKERR prb_mng_ptr->buildSubProblem(
1351  subdm_field->problemName, subdm_field->rowFields, subdm_field->colFields,
1352  subdm_field->problemMainOfSubPtr->getName(),
1353  subdm_field->isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1354  subdm_field->verbosity);
1355 
1356  // partition problem
1357  subdm_field->isPartitioned = subdm_field->isPartitioned;
1358  if (subdm_field->isPartitioned) {
1359  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, true,
1360  0, subdm_field->sIze,
1361  subdm_field->verbosity);
1362  // set ghost nodes
1364  subdm_field->problemName, subdm_field->verbosity);
1365  } else {
1366  // partition finite elements
1367  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, false,
1368  -1, -1, subdm_field->verbosity);
1369  // set ghost nodes
1370  CHKERR prb_mng_ptr->partitionGhostDofs(subdm_field->problemName,
1371  subdm_field->verbosity);
1372  }
1373 
1374  subdm_field->isProblemBuild = PETSC_TRUE;
1375 
1377 }
1378 
1379 PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec g, InsertMode mode,
1380  Vec l) {
1381  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1383  CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
1385 }
1386 
1387 PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec g, InsertMode mode, Vec l) {
1389  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1391 
1392  CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
1393 
1394  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1395  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1396  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1397 
1398  double *array_loc, *array_glob;
1399  CHKERR VecGetArray(l, &array_loc);
1400  CHKERR VecGetArray(g, &array_glob);
1401  switch (mode) {
1402  case INSERT_VALUES:
1403  cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1404  break;
1405  case ADD_VALUES:
1406  cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1407  break;
1408  default:
1409  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1410  }
1411  CHKERR VecRestoreArray(l, &array_loc);
1412  CHKERR VecRestoreArray(g, &array_glob);
1414 }
1415 
1416 PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM dm, Vec l, InsertMode mode,
1417  Vec g) {
1418 
1419  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1421 
1422  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1423  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1424  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1425 
1426  double *array_loc, *array_glob;
1427  CHKERR VecGetArray(l, &array_loc);
1428  CHKERR VecGetArray(g, &array_glob);
1429  switch (mode) {
1430  case INSERT_VALUES:
1431  cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1432  break;
1433  case ADD_VALUES:
1434  cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1435  break;
1436  default:
1437  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1438  }
1439  CHKERR VecRestoreArray(l, &array_loc);
1440  CHKERR VecRestoreArray(g, &array_glob);
1441 
1443 }
1444 
1445 PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec l, InsertMode mode, Vec g) {
1446  //
1449 }
1450 
1451 PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields,
1452  char ***fieldNames, IS **fields) {
1453  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1455 
1456  if (numFields) {
1457  *numFields = 0;
1458  }
1459  if (fieldNames) {
1460  *fieldNames = NULL;
1461  }
1462  if (fields) {
1463  *fields = NULL;
1464  }
1465 
1466  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1467  auto fields_ptr = dm_field->mField_ptr->get_fields();
1468  Field_multiIndex::iterator fit, hi_fit;
1469  fit = fields_ptr->begin();
1470  hi_fit = fields_ptr->end();
1471  *numFields = std::distance(fit, hi_fit);
1472 
1473  if (fieldNames) {
1474  CHKERR PetscMalloc1(*numFields, fieldNames);
1475  }
1476  if (fields) {
1477  CHKERR PetscMalloc1(*numFields, fields);
1478  }
1479 
1480  for (int f = 0; fit != hi_fit; fit++, f++) {
1481  if (fieldNames) {
1482  CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1483  (char **)&(*fieldNames)[f]);
1484  }
1485  if (fields) {
1486  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1487  ->isCreateProblemFieldAndRank(
1488  dm_field->problemPtr->getName(), ROW, fit->get()->getName(), 0,
1489  fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1490  }
1491  }
1492 
1494 }
1495 
1496 PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[],
1497  IS *is) {
1498  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1500  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1501  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1502  ->isCreateProblemFieldAndRank(dm_field->problemPtr->getName(), ROW,
1503  field_name, 0, 1000, is);
1505 }
1506 
1507 PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb) {
1508  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1510  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1511  dm_field->verbosity = verb;
1513 }
1514 
1515 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::DMCtx::destroyProblem
PetscBool destroyProblem
If true destroy problem with DM.
Definition: DMMoFEM.hpp:974
MoFEM::KspRhs
PetscErrorCode KspRhs(KSP ksp, Vec f, void *ctx)
Run over elements in the lists.
Definition: KspCtx.cpp:21
MoFEM::DMSetFromOptions_MoFEM
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition: DMMoFEM.cpp:1263
g
constexpr double g
Definition: shallow_wave.cpp:63
MoFEM::DMMoFEMSetSnesCtx
PetscErrorCode DMMoFEMSetSnesCtx(DM dm, boost::shared_ptr< MoFEM::SnesCtx > snes_ctx)
Set MoFEM::SnesCtx data structure.
Definition: DMMoFEM.cpp:1115
MoFEM::LogManager::checkIfChannelExist
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
Definition: LogManager.cpp:404
MoFEM::CoreInterface::problem_basic_method_postProcess
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::CoreInterface::loop_dofs
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.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::DMMoFEMSetKspCtx
PetscErrorCode DMMoFEMSetKspCtx(DM dm, boost::shared_ptr< MoFEM::KspCtx > ksp_ctx)
set MoFEM::KspCtx data structure
Definition: DMMoFEM.cpp:1089
MoFEM::DMMoFEMKSPSetComputeRHS
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set KSP right hand side evaluation function
Definition: DMMoFEM.cpp:641
MoFEM::DMCtx::referenceNumber
int referenceNumber
Definition: DMMoFEM.hpp:980
MoFEM::DMMoFEMAddSubFieldCol
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:284
MoFEM::DMMoFEMGetDestroyProblem
PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem)
Definition: DMMoFEM.cpp:451
MoFEM::Problem::getNumberOfElementsByNameAndPart
MoFEMErrorCode getNumberOfElementsByNameAndPart(MPI_Comm comm, const std::string name, PetscLayout *layout) const
Get number of finite elements by name on processors.
Definition: ProblemsMultiIndices.cpp:98
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::DMCtx::isCompDM
PetscBool isCompDM
Definition: DMMoFEM.hpp:968
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
MoFEM::TsSetIFunction
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:56
MoFEM::DMMoFEMSetSquareProblem
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition: DMMoFEM.cpp:460
MoFEM::DMCtx::problemPtr
const Problem * problemPtr
pinter to problem data structure
Definition: DMMoFEM.hpp:960
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::DMGlobalToLocalBegin_MoFEM
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMoFEM.cpp:1379
MoFEM::DMMoFEMGetSubColIS
PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is)
get sub problem is
Definition: DMMoFEM.cpp:353
MoFEM::DMMoFEMAddColCompositeProblem
PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on col.
Definition: DMMoFEM.cpp:389
MoFEM::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
DMMoFEM.hpp
Discrete manager interface for MoFEM.
MoFEM::DMSetUp_MoFEM
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMoFEM.cpp:1285
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:527
MoFEM::TsSetI2Jacobian
PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, PetscReal a, PetscReal aa, Mat A, Mat B, void *ctx)
Calculation Jacobian for second order PDE in time.
Definition: TsCtx.cpp:511
MoFEM::Problem::getNbGhostDofsRow
DofIdx getNbGhostDofsRow() const
Definition: ProblemsMultiIndices.hpp:380
MoFEM::DMMoFEMTSSetRHSJacobian
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: DMMoFEM.cpp:926
MoFEM::DMMoFEMGetIsSubDM
PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm)
Definition: DMMoFEM.cpp:326
MoFEM::DofMethod
Data structure to exchange data between mofem and User Loop Methods on entities.
Definition: LoopMethods.hpp:493
MoFEM::DMMoFEMKSPSetComputeOperators
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: DMMoFEM.cpp:682
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
MoFEM::DMLocalToGlobalBegin_MoFEM
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
Definition: DMMoFEM.cpp:1416
MoFEM::DMCtx::sIze
int sIze
Definition: DMMoFEM.hpp:957
MoFEM::DMMoFEMSwapDMCtx
PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap)
Swap internal data struture.
Definition: DMMoFEM.cpp:200
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:501
MoFEM::CoreInterface::check_problem
virtual bool check_problem(const std::string name)=0
check if problem exist
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
ROW
@ ROW
Definition: definitions.h:123
MoFEM::DMMoFEMGetIsCompDM
PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm)
Get if this DM is composite DM.
Definition: DMMoFEM.cpp:409
MoFEM::TsMonitorSet
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:259
MoFEM::DMCreateGlobalVector_MoFEM
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1171
MoFEM::DMCtx::mapTypeRow
std::map< std::string, boost::shared_ptr< Range > > mapTypeRow
Definition: DMMoFEM.hpp:971
MoFEM::DMMoFEMSetDestroyProblem
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
Definition: DMMoFEM.cpp:442
MoFEM::DMCtx::mapTypeCol
std::map< std::string, boost::shared_ptr< Range > > mapTypeCol
Definition: DMMoFEM.hpp:972
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
VERBOSE
@ VERBOSE
Definition: definitions.h:209
MoFEM::CoreInterface::problem_basic_method_preProcess
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
MoFEM::DMCreateMatrix_MoFEM
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1201
MoFEM::DMMoFEMAddRowCompositeProblem
PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on row.
Definition: DMMoFEM.cpp:371
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::DMCtx::rAnk
int rAnk
Definition: DMMoFEM.hpp:957
MoFEM::DMMoFEMTSSetIJacobian
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: DMMoFEM.cpp:857
MoFEM::DMoFEMLoopFiniteElementsUpAndLowRank
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:571
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEM::DMMoFEMTSSetI2Jacobian
PetscErrorCode DMMoFEMTSSetI2Jacobian(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: DMMoFEM.cpp:1021
MoFEM::DMDestroy_MoFEM
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
Definition: DMMoFEM.cpp:83
MoFEM::ProblemsManager::partitionGhostDofsOnDistributedMesh
MoFEMErrorCode partitionGhostDofsOnDistributedMesh(const std::string name, int verb=VERBOSE)
determine ghost nodes on distributed meshes
Definition: ProblemsManager.cpp:2473
MoFEM::DMSubDMSetUp_MoFEM
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
Definition: DMMoFEM.cpp:1332
MoFEM::DMMoFEMTSSetI2Function
PetscErrorCode DMMoFEMTSSetI2Function(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 implicit function evaluation function
Definition: DMMoFEM.cpp:979
MoFEM::DMCtx::verbosity
int verbosity
verbosity
Definition: DMMoFEM.hpp:979
MoFEM::DMMoFEMCreateSubDM
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition: DMMoFEM.cpp:219
RowColData
RowColData
RowColData.
Definition: definitions.h:123
MoFEM::DMMoFEMGetSubRowIS
PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is)
get sub problem is
Definition: DMMoFEM.cpp:335
MoFEM::PairNameFEMethodPtr
Definition: AuxPETSc.hpp:12
MoFEM::Problem::getSubData
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
Definition: ProblemsMultiIndices.hpp:119
MoFEM::LogManager::getStrmSync
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
MoFEM::DMCtx::kspCtx
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
Definition: DMMoFEM.hpp:982
MoFEM::DMCtx::colCompPrb
std::vector< std::string > colCompPrb
Definition: DMMoFEM.hpp:970
MoFEM::DMoFEMMeshToGlobalVector
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition: DMMoFEM.cpp:539
COL
@ COL
Definition: definitions.h:123
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::CoreInterface::modify_problem_unset_finite_element
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...
MoFEM::DMMoFEMGetSnesCtx
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition: DMMoFEM.cpp:1098
MoFEM::TsCtx
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:17
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
MoFEM::DMCtx::isProblemBuild
PetscBool isProblemBuild
True if problem is build.
Definition: DMMoFEM.hpp:950
MoFEM::DMCtx::tsCtx
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
Definition: DMMoFEM.hpp:984
MoFEM::DMMoFEMSNESSetJacobian
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: DMMoFEM.cpp:763
MoFEM::TsSetIJacobian
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobian in TS solver.
Definition: TsCtx.cpp:165
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
MoFEM::DMoFEMPreProcessFiniteElements
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:550
MoFEM::DMCtx
PETSc Discrete Manager data structure.
Definition: DMMoFEM.hpp:944
MoFEM::DMCtx::isSquareMatrix
PetscBool isSquareMatrix
true if rows equals to cols
Definition: DMMoFEM.hpp:955
MoFEM::DMMoFEMGetSquareProblem
PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem)
get squared problem
Definition: DMMoFEM.cpp:491
MoFEM::DMGlobalToLocalEnd_MoFEM
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMoFEM.cpp:1387
MoFEM::BasicMethod
Data structure to exchange data between mofem and User Loop Methods.
Definition: LoopMethods.hpp:183
MoFEM::DMoFEMPostProcessFiniteElements
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:560
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::Problem::getNbLocalDofsRow
DofIdx getNbLocalDofsRow() const
Definition: ProblemsMultiIndices.hpp:378
MoFEM::DMMoFEMCreateMoFEM
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: DMMoFEM.cpp:118
MoFEM::ProblemsManager::partitionFiniteElements
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
Definition: ProblemsManager.cpp:2167
MoFEM::DMCtx::isPartitioned
PetscBool isPartitioned
true if read mesh is on parts
Definition: DMMoFEM.hpp:954
MoFEM::DMCtx::DMCtx
DMCtx()
Definition: DMMoFEM.cpp:18
MoFEM::DMSetOperators_MoFEM
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition: DMMoFEM.cpp:53
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
MoFEM::DMCreateFieldIS_MoFEM
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition: DMMoFEM.cpp:1451
MoFEM::ProblemsManager::buildComposedProblem
MoFEMErrorCode buildComposedProblem(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
Definition: ProblemsManager.cpp:1250
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
MoFEM::DMCtx::snesCtx
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
Definition: DMMoFEM.hpp:983
MoFEM::DMMoFEMResolveSharedFiniteElements
PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name)
Resolve shared entities.
Definition: DMMoFEM.cpp:468
MoFEM::DMMoFEMGetKspCtx
PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx)
get MoFEM::KspCtx data structure
Definition: DMMoFEM.cpp:1072
MoFEM::CoreInterface::get_fields
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MoFEM::DMMoFEMSetTsCtx
PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr< MoFEM::TsCtx > ts_ctx)
Set MoFEM::TsCtx data structure.
Definition: DMMoFEM.cpp:1163
MoFEM::SnesRhs
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:27
MoFEM::DMMoFEMTSSetRHSFunction
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: DMMoFEM.cpp:886
MoFEM::DMCreate_MoFEM
PetscErrorCode DMCreate_MoFEM(DM dm)
Create dm data structure with MoFEM data structure.
Definition: DMMoFEM.cpp:75
MoFEM::DMCtx::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: DMMoFEM.cpp:41
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
MoFEM::DMMoFEMGetFieldIS
PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[], IS *is)
get field is in the problem
Definition: DMMoFEM.cpp:1496
MoFEM::DMCtx::mField_ptr
Interface * mField_ptr
MoFEM interface.
Definition: DMMoFEM.hpp:949
MoFEM::DMoFEMGetInterfacePtr
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMoFEM.cpp:418
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::KspMat
PetscErrorCode KspMat(KSP ksp, Mat A, Mat B, void *ctx)
Run over elements in the list.
Definition: KspCtx.cpp:86
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::DMMoFEMTSSetMonitor
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: DMMoFEM.cpp:1060
MoFEM::DMMoFEMTSSetIFunction
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: DMMoFEM.cpp:804
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ProblemsManager::buildProblemOnDistributedMesh
MoFEMErrorCode buildProblemOnDistributedMesh(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures, assuming that mesh is distributed (collective)
Definition: ProblemsManager.cpp:290
MoFEM::DMMoFEMGetIsPartitioned
PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned)
Definition: DMMoFEM.cpp:1138
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:430
MoFEM::DMMoFEMUnSetElement
PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name)
unset element from dm
Definition: DMMoFEM.cpp:518
MoFEM::CacheTupleWeakPtr
boost::weak_ptr< CacheTuple > CacheTupleWeakPtr
Definition: FEMultiIndices.hpp:494
MoFEM::Problem::getBitRefLevelMask
BitRefLevel getBitRefLevelMask() const
Definition: ProblemsMultiIndices.hpp:384
MoFEM::DMCtx::colFields
std::vector< std::string > colFields
Definition: DMMoFEM.hpp:965
MoFEM::DMMoFEMGetTsCtx
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition: DMMoFEM.cpp:1146
MoFEM::Problem::getBitRefLevel
BitRefLevel getBitRefLevel() const
Definition: ProblemsMultiIndices.hpp:383
MoFEM::KspCtx
Interface for linear (KSP) solver.
Definition: KspCtx.hpp:13
MoFEM::TsSetRHSJacobian
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition: TsCtx.cpp:422
MoFEM::CoreInterface::modify_problem_add_finite_element
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
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::DMCtx::rowFields
std::vector< std::string > rowFields
Definition: DMMoFEM.hpp:964
MoFEM::CoreInterface::modify_problem_mask_ref_level_add_bit
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
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::SnesCtx
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:13
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::TsSetRHSFunction
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition: TsCtx.cpp:323
MoFEM::DMMoFEMAddSubFieldRow
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:242
MoFEM::DMCtx::isSubDM
PetscBool isSubDM
Definition: DMMoFEM.hpp:963
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::SmartPetscObj< Vec >
MoFEM::TsSetI2Function
PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, Vec F, void *ctx)
Calculation the right hand side for second order PDE in time.
Definition: TsCtx.cpp:612
MoFEM::SnesMat
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:139
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::DMMoFEMGetProblemFiniteElementLayout
PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name, PetscLayout *layout)
Get finite elements layout in the problem.
Definition: DMMoFEM.cpp:477
MoFEM::DMMoFEMDuplicateDMCtx
PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate)
Duplicate internal data struture.
Definition: DMMoFEM.cpp:184
MF_EXCL
@ MF_EXCL
Definition: definitions.h:99
MoFEM::Problem::getName
auto getName() const
Definition: ProblemsMultiIndices.hpp:372
MoFEM::DMLocalToGlobalEnd_MoFEM
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition: DMMoFEM.cpp:1445
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:590
MF_EXIST
@ MF_EXIST
Definition: definitions.h:100
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::DMoFEMLoopDofs
PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[], MoFEM::DofMethod *method)
execute method for dofs on field in problem
Definition: DMMoFEM.cpp:609
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::CoreInterface::delete_problem
virtual MoFEMErrorCode delete_problem(const std::string name)=0
Delete problem.
MoFEM::DMCreateLocalVector_MoFEM
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVector.
Definition: DMMoFEM.cpp:1191
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
MoFEM::DMCtx::problemName
std::string problemName
Problem name.
Definition: DMMoFEM.hpp:951
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
MoFEM::DMCtx::rowCompPrb
std::vector< std::string > rowCompPrb
Definition: DMMoFEM.hpp:969
MoFEM::ProblemsManager::partitionProblem
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
Definition: ProblemsManager.cpp:1683
MoFEM::DMCtx::problemMainOfSubPtr
const Problem * problemMainOfSubPtr
pinter to main problem to sub-problem
Definition: DMMoFEM.hpp:966
MoFEM::DMMoFEMSetVerbosity
PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb)
Set verbosity level.
Definition: DMMoFEM.cpp:1507
MoFEM::ProblemsManager::buildSubProblem
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, boost::shared_ptr< Range >> *entityMapRow=nullptr, const map< std::string, boost::shared_ptr< Range >> *entityMapCol=nullptr, int verb=VERBOSE)
build sub problem
Definition: ProblemsManager.cpp:981
MoFEM::DMMoFEMSetIsPartitioned
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1127
MoFEM::DMMoFEMSNESSetFunction
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: DMMoFEM.cpp:722