v0.8.21
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 DMMoFEMResolveSharedEntities(DM dm, const char fe_name[]) {
378  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
380  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
381  CHKERR dm_field->mField_ptr->resolve_shared_ents(dm_field->problemPtr,
382  fe_name);
384 }
385 
386 PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, const char fe_name[],
387  PetscLayout *layout) {
389  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
391  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
392 
393  MPI_Comm comm;
394  CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
395  CHKERR dm_field->problemPtr->getNumberOfElementsByNameAndPart(comm, fe_name,
396  layout);
398 }
399 
400 PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem) {
403  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
405  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
406  *square_problem = dm_field->isSquareMatrix;
408 }
409 
410 PetscErrorCode DMMoFEMAddElement(DM dm, const char fe_name[]) {
411  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
415  dm_field->problemName, fe_name);
416  CHKERRG(ierr);
418 }
419 
420 PetscErrorCode DMMoFEMUnSetElement(DM dm, const char fe_name[]) {
421  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
423  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
425  dm_field->problemName, fe_name);
426  CHKERRG(ierr);
428 }
429 
430 PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode,
431  ScatterMode scatter_mode) {
433  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
435  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
436  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setLocalGhostVector(
437  dm_field->problemPtr, COL, l, mode, scatter_mode);
438  CHKERRG(ierr);
440 }
441 
442 PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode,
443  ScatterMode scatter_mode) {
444  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
446  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
447  ierr = dm_field->mField_ptr->getInterface<VecManager>()->setGlobalGhostVector(
448  dm_field->problemPtr, COL, g, mode, scatter_mode);
449  CHKERRG(ierr);
451 }
452 
453 PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
454  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
456  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
458  dm_field->problemPtr, *method);
459  CHKERRG(ierr);
461 }
462 
463 PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method) {
464  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
466  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
468  dm_field->problemPtr, *method);
469  CHKERRG(ierr);
471 }
472 
473 PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[],
474  MoFEM::FEMethod *method,
475  int low_rank, int up_rank) {
477  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
478  ierr = dm_field->mField_ptr->loop_finite_elements(
479  dm_field->problemPtr, fe_name, *method, low_rank, up_rank);
480  CHKERRG(ierr);
482 }
483 
484 PetscErrorCode
485 DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const std::string &fe_name,
486  boost::shared_ptr<MoFEM::FEMethod> method,
487  int low_rank, int up_rank) {
488  return DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name.c_str(), method.get(),
489  low_rank, up_rank);
490 }
491 
492 PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[],
493  MoFEM::FEMethod *method) {
494  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
496  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
497  ierr = DMoFEMLoopFiniteElementsUpAndLowRank(dm, fe_name, method,
498  dm_field->rAnk, dm_field->rAnk);
499  CHKERRG(ierr);
501 }
502 
503 PetscErrorCode
504 DMoFEMLoopFiniteElements(DM dm, const std::string &fe_name,
505  boost::shared_ptr<MoFEM::FEMethod> method) {
506  return DMoFEMLoopFiniteElements(dm, fe_name.c_str(), method.get());
507 }
508 
509 PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[],
510  MoFEM::DofMethod *method) {
511  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
513  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
514  ierr =
515  dm_field->mField_ptr->loop_dofs(dm_field->problemPtr, field_name, COL,
516  *method, dm_field->rAnk, dm_field->rAnk);
517  CHKERRG(ierr);
519 }
520 
521 template <class S, class T0, class T1, class T2>
522 static PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, S fe_name, T0 method,
523  T1 pre_only, T2 post_only) {
524  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
526  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
527  if (pre_only) {
528  dm_field->kspCtx->get_preProcess_to_do_Rhs().push_back(pre_only);
529  }
530  if (method) {
531  dm_field->kspCtx->get_loops_to_do_Rhs().push_back(
532  PairNameFEMethodPtr(fe_name, method));
533  }
534  if (post_only) {
535  dm_field->kspCtx->get_postProcess_to_do_Rhs().push_back(post_only);
536  }
537  CHKERR DMKSPSetComputeRHS(dm, KspRhs, dm_field->kspCtx.get());
539 }
540 
541 PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[],
542  MoFEM::FEMethod *method,
543  MoFEM::BasicMethod *pre_only,
544  MoFEM::BasicMethod *post_only) {
545  return DMMoFEMKSPSetComputeRHS<const char *, MoFEM::FEMethod *,
547  dm, fe_name, method, pre_only, post_only);
548 }
549 
550 PetscErrorCode
551 DMMoFEMKSPSetComputeRHS(DM dm, const std::string &fe_name,
552  boost::shared_ptr<MoFEM::FEMethod> method,
553  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
554  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
555  return DMMoFEMKSPSetComputeRHS<const std::string &,
556  boost::shared_ptr<MoFEM::FEMethod>,
557  boost::shared_ptr<MoFEM::BasicMethod>,
558  boost::shared_ptr<MoFEM::BasicMethod>>(
559  dm, fe_name, method, pre_only, post_only);
560 }
561 
562 template <class S, class T0, class T1, class T2>
563 static PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, S fe_name, T0 method,
564  T1 pre_only, T2 post_only) {
565  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
567  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
568  if (pre_only) {
569  dm_field->kspCtx->get_preProcess_to_do_Mat().push_back(pre_only);
570  }
571  if (method) {
572  dm_field->kspCtx->get_loops_to_do_Mat().push_back(
573  PairNameFEMethodPtr(fe_name, method));
574  }
575  if (post_only) {
576  dm_field->kspCtx->get_postProcess_to_do_Mat().push_back(post_only);
577  }
578  CHKERR DMKSPSetComputeOperators(dm, KspMat, dm_field->kspCtx.get());
580 }
581 
582 PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[],
583  MoFEM::FEMethod *method,
584  MoFEM::BasicMethod *pre_only,
585  MoFEM::BasicMethod *post_only) {
586  return DMMoFEMKSPSetComputeOperators<const char *, MoFEM::FEMethod *,
589  dm, fe_name, method, pre_only, post_only);
590 }
591 
592 PetscErrorCode
593 DMMoFEMKSPSetComputeOperators(DM dm, const std::string &fe_name,
594  boost::shared_ptr<MoFEM::FEMethod> method,
595  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
596  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
597  return DMMoFEMKSPSetComputeOperators<const std::string &,
598  boost::shared_ptr<MoFEM::FEMethod>>(
599  dm, fe_name, method, pre_only, post_only);
600 }
601 
602 template <class S, class T0, class T1, class T2>
603 static PetscErrorCode DMMoFEMSNESSetFunction(DM dm, S fe_name, T0 method,
604  T1 pre_only, T2 post_only) {
605  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
607  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
608  if (pre_only) {
609  dm_field->snesCtx->get_preProcess_to_do_Rhs().push_back(pre_only);
610  }
611  if (method) {
612  dm_field->snesCtx->get_loops_to_do_Rhs().push_back(
613  PairNameFEMethodPtr(fe_name, method));
614  }
615  if (post_only) {
616  dm_field->snesCtx->get_postProcess_to_do_Rhs().push_back(post_only);
617  }
618  CHKERR DMSNESSetFunction(dm, SnesRhs, dm_field->snesCtx.get());
620 }
621 
622 PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[],
623  MoFEM::FEMethod *method,
624  MoFEM::BasicMethod *pre_only,
625  MoFEM::BasicMethod *post_only) {
626  return DMMoFEMSNESSetFunction<const char *, MoFEM::FEMethod *,
628  dm, fe_name, method, pre_only, post_only);
629 }
630 
631 PetscErrorCode
632 DMMoFEMSNESSetFunction(DM dm, const std::string &fe_name,
633  boost::shared_ptr<MoFEM::FEMethod> method,
634  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
635  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
636  return DMMoFEMSNESSetFunction<const std::string &,
637  boost::shared_ptr<MoFEM::FEMethod>,
638  boost::shared_ptr<MoFEM::BasicMethod>,
639  boost::shared_ptr<MoFEM::BasicMethod>>(
640  dm, fe_name, method, pre_only, post_only);
641 }
642 
643 template <class S, class T0, class T1, class T2>
644 static PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, S fe_name, T0 method,
645  T1 pre_only, T2 post_only) {
646  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
648  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
649  if (pre_only) {
650  dm_field->snesCtx->get_preProcess_to_do_Mat().push_back(pre_only);
651  }
652  if (method) {
653  dm_field->snesCtx->get_loops_to_do_Mat().push_back(
654  PairNameFEMethodPtr(fe_name, method));
655  }
656  if (post_only) {
657  dm_field->snesCtx->get_postProcess_to_do_Mat().push_back(post_only);
658  }
659  CHKERR DMSNESSetJacobian(dm, SnesMat, dm_field->snesCtx.get());
661 }
662 
663 PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[],
664  MoFEM::FEMethod *method,
665  MoFEM::BasicMethod *pre_only,
666  MoFEM::BasicMethod *post_only) {
667  return DMMoFEMSNESSetJacobian<const char *, MoFEM::FEMethod *,
669  dm, fe_name, method, pre_only, post_only);
670 }
671 
672 PetscErrorCode
673 DMMoFEMSNESSetJacobian(DM dm, const std::string &fe_name,
674  boost::shared_ptr<MoFEM::FEMethod> method,
675  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
676  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
677  return DMMoFEMSNESSetJacobian<const std::string &,
678  boost::shared_ptr<MoFEM::FEMethod>,
679  boost::shared_ptr<MoFEM::BasicMethod>,
680  boost::shared_ptr<MoFEM::BasicMethod>>(
681  dm, fe_name, method, pre_only, post_only);
682 }
683 
684 template <class S, class T0, class T1, class T2>
685 static PetscErrorCode DMMoFEMTSSetIFunction(DM dm, S fe_name, T0 method,
686  T1 pre_only, T2 post_only) {
687  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
689  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
690  if (pre_only) {
691  dm_field->tsCtx->get_preProcess_to_do_IFunction().push_back(pre_only);
692  }
693  if (method) {
694  dm_field->tsCtx->get_loops_to_do_IFunction().push_back(
695  PairNameFEMethodPtr(fe_name, method));
696  }
697  if (post_only) {
698  dm_field->tsCtx->get_postProcess_to_do_IFunction().push_back(post_only);
699  }
700  CHKERR DMTSSetIFunction(dm, TsSetIFunction, dm_field->tsCtx.get());
702 }
703 
704 PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[],
705  MoFEM::FEMethod *method,
706  MoFEM::BasicMethod *pre_only,
707  MoFEM::BasicMethod *post_only) {
708  return DMMoFEMTSSetIFunction<const char *, MoFEM::FEMethod *,
710  dm, fe_name, method, pre_only, post_only);
712 }
713 
714 PetscErrorCode
715 DMMoFEMTSSetIFunction(DM dm, const std::string &fe_name,
716  boost::shared_ptr<MoFEM::FEMethod> method,
717  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
718  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
719  return DMMoFEMTSSetIFunction<const std::string,
720  boost::shared_ptr<MoFEM::FEMethod>,
721  boost::shared_ptr<MoFEM::BasicMethod>,
722  boost::shared_ptr<MoFEM::BasicMethod>>(
723  dm, fe_name, method, pre_only, post_only);
725 }
726 
727 template <class S, class T0, class T1, class T2>
728 static PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, S fe_name, T0 method,
729  T1 pre_only, T2 post_only) {
730  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
732  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
733  if (pre_only) {
734  dm_field->tsCtx->get_preProcess_to_do_IJacobian().push_back(pre_only);
735  }
736  if (method) {
737  dm_field->tsCtx->get_loops_to_do_IJacobian().push_back(
738  PairNameFEMethodPtr(fe_name, method));
739  }
740  if (post_only) {
741  dm_field->tsCtx->get_postProcess_to_do_IJacobian().push_back(post_only);
742  }
743  CHKERR DMTSSetIJacobian(dm, TsSetIJacobian, dm_field->tsCtx.get());
745 }
746 
747 PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const char fe_name[],
748  MoFEM::FEMethod *method,
749  MoFEM::BasicMethod *pre_only,
750  MoFEM::BasicMethod *post_only) {
751  return DMMoFEMTSSetIJacobian<const char *, FEMethod *, MoFEM::BasicMethod *,
752  MoFEM::BasicMethod *>(dm, fe_name, method,
753  pre_only, post_only);
754 }
755 
756 PetscErrorCode
757 DMMoFEMTSSetIJacobian(DM dm, const std::string &fe_name,
758  boost::shared_ptr<MoFEM::FEMethod> method,
759  boost::shared_ptr<MoFEM::BasicMethod> pre_only,
760  boost::shared_ptr<MoFEM::BasicMethod> post_only) {
761  return DMMoFEMTSSetIJacobian<const std::string &,
762  boost::shared_ptr<MoFEM::FEMethod>,
763  boost::shared_ptr<MoFEM::BasicMethod>,
764  boost::shared_ptr<MoFEM::BasicMethod>>(
765  dm, fe_name, method, pre_only, post_only);
766 }
767 
768 PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx) {
769  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
771  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
772  *ksp_ctx = dm_field->kspCtx.get();
774 }
775 
776 PetscErrorCode
777 DMMoFEMGetKspCtx(DM dm, const boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx) {
778  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
780  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
781  const_cast<boost::shared_ptr<MoFEM::KspCtx> &>(ksp_ctx) = dm_field->kspCtx;
783 }
784 
785 PetscErrorCode DMMoFEMSetKspCtx(DM dm,
786  boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx) {
787  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
789  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
790  dm_field->kspCtx = ksp_ctx;
792 }
793 
794 PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx) {
795  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
797  DMCtx *dm_field = (DMCtx *)dm->data;
798  *snes_ctx = dm_field->snesCtx.get();
800 }
801 
802 PetscErrorCode
803 DMMoFEMGetSnesCtx(DM dm, const boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx) {
804  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
806  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
807  const_cast<boost::shared_ptr<MoFEM::SnesCtx> &>(snes_ctx) = dm_field->snesCtx;
809 }
810 
811 PetscErrorCode DMMoFEMSetSnesCtx(DM dm,
812  boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx) {
813  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
815  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
816  dm_field->snesCtx = snes_ctx;
818 }
819 
820 /** get if read mesh is partitioned
821  * \ingroup dm
822  */
823 PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned) {
824  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
826  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
827  dm_field->isPartitioned = is_partitioned;
829 }
830 
831 /** get if read mesh is partitioned
832  * \ingroup dm
833  */
834 PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned) {
835  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
837  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
838  *is_partitioned = dm_field->isPartitioned;
840 }
841 
842 PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx) {
843  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
845  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
846  *ts_ctx = dm_field->tsCtx.get();
848 }
849 
850 PetscErrorCode DMMoFEMGetTsCtx(DM dm,
851  const boost::shared_ptr<MoFEM::TsCtx> &ts_ctx) {
852  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
854  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
855  const_cast<boost::shared_ptr<MoFEM::TsCtx> &>(ts_ctx) = dm_field->tsCtx;
857 }
858 
859 PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr<MoFEM::TsCtx> &ts_ctx) {
860  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
862  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
863  dm_field->tsCtx = ts_ctx;
865 }
866 
867 PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g) {
868  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
870  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
871  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateGhost(
872  dm_field->problemName, COL, g);
874 }
875 
876 PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l) {
877  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
879  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
880  CHKERR dm_field->mField_ptr->getInterface<VecManager>()->vecCreateSeq(
881  dm_field->problemName, COL, l);
883 }
884 
885 PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M) {
886  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
888  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
889  if (strcmp(dm->mattype, MATMPIAIJ) == 0) {
891  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(dm_field->problemName,
892  M);
893  } else if (strcmp(dm->mattype, MATAIJ) == 0) {
895  ->createSeqAIJWithArrays<PetscLocalIdx_mi_tag>(dm_field->problemName,
896  M);
897  } else {
898  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
899  "Matrix type not implemented");
900  }
902 }
903 
904 #if PETSC_VERSION_GE(3, 7, 0)
905 PetscErrorCode DMSetFromOptions_MoFEM(PetscOptionItems *PetscOptionsObject,
906  DM dm) {
907 #elif PETSC_VERSION_GE(3, 5, 3)
908 PetscErrorCode DMSetFromOptions_MoFEM(PetscOptions *PetscOptionsObject, DM dm) {
909 #else
910 PetscErrorCode DMSetFromOptions_MoFEM(DM dm) {
911 #endif
912 
913  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
915  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
916 #if PETSC_VERSION_GE(3, 5, 3)
917  ierr = PetscOptionsHead(PetscOptionsObject, "DMMoFEM Options");
918  CHKERRG(ierr);
919 #else
920  ierr = PetscOptionsHead("DMMoFEM Options");
921  CHKERRG(ierr);
922 #endif
923  ierr = PetscOptionsBool("-dm_is_partitioned",
924  "set if mesh is partitioned (works which native MOAB "
925  "file format, i.e. h5m",
926  "DMSetUp", dm_field->isPartitioned,
927  &dm_field->isPartitioned, NULL);
928  CHKERRG(ierr);
930 }
931 
932 PetscErrorCode DMSetUp_MoFEM(DM dm) {
933  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
934  ProblemsManager *prb_mng_ptr;
936  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
937  CHKERR dm_field->mField_ptr->getInterface(prb_mng_ptr);
938 
939  if (dm_field->isCompDM) {
940  // It is composite probelm
941  CHKERR prb_mng_ptr->buildCompsedProblem(
942  dm_field->problemName, dm_field->rowCompPrb, dm_field->colCompPrb,
943  dm_field->isSquareMatrix == PETSC_TRUE, dm_field->verbosity);
944  } else {
945  if (dm_field->isPartitioned) {
947  dm_field->problemName, dm_field->isSquareMatrix == PETSC_TRUE,
948  dm_field->verbosity);
949  } else {
950  CHKERR prb_mng_ptr->buildProblem(dm_field->problemName,
951  dm_field->isSquareMatrix == PETSC_TRUE,
952  dm_field->verbosity);
953  CHKERR prb_mng_ptr->partitionProblem(dm_field->problemName,
954  dm_field->verbosity);
955  }
956  }
957 
958  // Partition finite elements
959  if (dm_field->isPartitioned) {
960  CHKERR prb_mng_ptr->partitionFiniteElements(
961  dm_field->problemName, true, 0, dm_field->sIze, dm_field->verbosity);
963  dm_field->problemName, dm_field->verbosity);
964  } else {
965  // partition finite elemnets
966  CHKERR prb_mng_ptr->partitionFiniteElements(dm_field->problemName, false,
967  -1, -1, dm_field->verbosity);
968  // Get ghost DOFs
969  CHKERR prb_mng_ptr->partitionGhostDofs(dm_field->problemName,
970  dm_field->verbosity);
971  }
972 
973  // Set flag that problem is build and partitioned
974  dm_field->isProblemBuild = PETSC_TRUE;
975 
977 }
978 
979 PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm) {
980  PetscValidHeaderSpecific(subdm, DM_CLASSID, 1);
981  ProblemsManager *prb_mng_ptr;
983 
984  DMCtx *subdm_field = static_cast<DMCtx *>(subdm->data);
985 
986  // build sub dm problem
987  CHKERR subdm_field->mField_ptr->getInterface(prb_mng_ptr);
988 
989  map<std::string, std::pair<EntityType, EntityType>> *entity_map_row = nullptr;
990  map<std::string, std::pair<EntityType, EntityType>> *entity_map_col = nullptr;
991 
992  if(subdm_field->mapTypeRow)
993  entity_map_row = subdm_field->mapTypeRow.get();
994  if(subdm_field->mapTypeCol)
995  entity_map_row = subdm_field->mapTypeCol.get();
996 
997  CHKERR prb_mng_ptr->buildSubProblem(
998  subdm_field->problemName, subdm_field->rowFields, subdm_field->colFields,
999  subdm_field->problemMainOfSubPtr->getName(),
1000  subdm_field->isSquareMatrix == PETSC_TRUE, entity_map_row, entity_map_col,
1001  subdm_field->verbosity);
1002 
1003  // partition problem
1004  subdm_field->isPartitioned = subdm_field->isPartitioned;
1005  if (subdm_field->isPartitioned) {
1006  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, true,
1007  0, subdm_field->sIze,
1008  subdm_field->verbosity);
1009  // set ghost nodes
1011  subdm_field->problemName, subdm_field->verbosity);
1012  } else {
1013  // partition finite elements
1014  CHKERR prb_mng_ptr->partitionFiniteElements(subdm_field->problemName, false,
1015  -1, -1, subdm_field->verbosity);
1016  // set ghost nodes
1017  CHKERR prb_mng_ptr->partitionGhostDofs(subdm_field->problemName,
1018  subdm_field->verbosity);
1019  }
1020 
1021  subdm_field->isProblemBuild = PETSC_TRUE;
1022 
1024 }
1025 
1026 PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec g, InsertMode mode,
1027  Vec l) {
1028  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1030  CHKERR VecGhostUpdateBegin(g, INSERT_VALUES, SCATTER_FORWARD);
1032 }
1033 
1034 PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec g, InsertMode mode, Vec l) {
1036  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1038 
1039  CHKERR VecGhostUpdateEnd(g, INSERT_VALUES, SCATTER_FORWARD);
1040 
1041  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1042  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1043  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1044 
1045  double *array_loc, *array_glob;
1046  CHKERR VecGetArray(l, &array_loc);
1047  CHKERR VecGetArray(g, &array_glob);
1048  switch (mode) {
1049  case INSERT_VALUES:
1050  cblas_dcopy(nb_dofs + nb_ghost, array_glob, 1, array_loc, 1);
1051  break;
1052  case ADD_VALUES:
1053  cblas_daxpy(nb_dofs + nb_ghost, 1, array_glob, 1, array_loc, 1);
1054  break;
1055  default:
1056  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1057  }
1058  CHKERR VecRestoreArray(l, &array_loc);
1059  CHKERR VecRestoreArray(g, &array_glob);
1061 }
1062 
1063 PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM dm, Vec l, InsertMode mode,
1064  Vec g) {
1065 
1066  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1068 
1069  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1070  int nb_dofs = dm_field->problemPtr->getNbLocalDofsRow();
1071  int nb_ghost = dm_field->problemPtr->getNbGhostDofsRow();
1072 
1073  double *array_loc, *array_glob;
1074  CHKERR VecGetArray(l, &array_loc);
1075  CHKERR VecGetArray(g, &array_glob);
1076  switch (mode) {
1077  case INSERT_VALUES:
1078  cblas_dcopy(nb_dofs + nb_ghost, array_loc, 1, array_glob, 1);
1079  break;
1080  case ADD_VALUES:
1081  cblas_daxpy(nb_dofs + nb_ghost, 1, array_loc, 1, array_glob, 1);
1082  break;
1083  default:
1084  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1085  }
1086  CHKERR VecRestoreArray(l, &array_loc);
1087  CHKERR VecRestoreArray(g, &array_glob);
1088 
1090 }
1091 
1092 PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec l, InsertMode mode, Vec g) {
1093  //
1096 }
1097 
1098 PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields,
1099  char ***fieldNames, IS **fields) {
1100  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1102 
1103  if (numFields) {
1104  *numFields = 0;
1105  }
1106  if (fieldNames) {
1107  *fieldNames = NULL;
1108  }
1109  if (fields) {
1110  *fields = NULL;
1111  }
1112 
1113  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1114  const Field_multiIndex *fields_ptr;
1115  CHKERR dm_field->mField_ptr->get_fields(&fields_ptr);
1116  Field_multiIndex::iterator fit, hi_fit;
1117  fit = fields_ptr->begin();
1118  hi_fit = fields_ptr->end();
1119  *numFields = std::distance(fit, hi_fit);
1120 
1121  if (fieldNames) {
1122  CHKERR PetscMalloc1(*numFields, fieldNames);
1123  }
1124  if (fields) {
1125  CHKERR PetscMalloc1(*numFields, fields);
1126  }
1127 
1128  for (int f = 0; fit != hi_fit; fit++, f++) {
1129  if (fieldNames) {
1130  CHKERR PetscStrallocpy(fit->get()->getName().c_str(),
1131  (char **)&(*fieldNames)[f]);
1132  }
1133  if (fields) {
1134  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1135  ->isCreateProblemFieldAndRank(
1136  dm_field->problemPtr->getName(), ROW, fit->get()->getName(), 0,
1137  fit->get()->getNbOfCoeffs(), &(*fields)[f]);
1138  }
1139  }
1140 
1142 }
1143 
1144 PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[],
1145  IS *is) {
1146  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1148  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1149  CHKERR dm_field->mField_ptr->getInterface<ISManager>()
1150  ->isCreateProblemFieldAndRank(dm_field->problemPtr->getName(), ROW,
1151  field_name, 0, 1000, is);
1153 }
1154 
1155 PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb) {
1156  PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1158  DMCtx *dm_field = static_cast<DMCtx *>(dm->data);
1159  dm_field->verbosity = verb;
1161 }
1162 
1163 } // namespace MoFEM
std::string problemName
Problem name.
Definition: DMMoFEM.hpp:748
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:622
MoFEM interface unique ID.
Problem manager is used to build and partition problems.
PetscBool isSquareMatrix
true if rows equals to cols
Definition: DMMoFEM.hpp:752
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:1155
int verbosity
verbosity
Definition: DMMoFEM.hpp:778
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:541
virtual MoFEMErrorCode delete_problem(const std::string name)=0
Delete problem.
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1026
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:769
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:760
Section manager is used to create indexes and sectionsFIXME: ISManager is not properly testsed by ato...
Definition: ISManager.hpp:36
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
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:430
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
PetscErrorCode DMMoFEMResolveSharedEntities(DM dm, const char fe_name[])
Resolve shared entities.
Definition: DMMMoFEM.cpp:376
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:410
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:110
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:746
base class for all interface classes
PetscBool isPartitioned
true if read mesh is on parts
Definition: DMMoFEM.hpp:751
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:475
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMMoFEM.cpp:932
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:542
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
const Problem * problemPtr
pinter to problem data structure
Definition: DMMoFEM.hpp:757
std::vector< std::string > colCompPrb
Definition: DMMoFEM.hpp:767
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:506
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:1063
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMMoFEM.cpp:885
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:582
DofIdx getNbGhostDofsRow() const
RowColData
RowColData.
Definition: definitions.h:184
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:492
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:773
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:867
PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, const char fe_name[], PetscLayout *layout)
Get finite elements layout in the problem.
Definition: DMMMoFEM.cpp:386
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
Definition: DMMoFEM.hpp:783
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:761
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:859
virtual MoFEMErrorCode loop_finite_elements(const std::string &problem_name, const std::string &fe_name, FEMethod &method, MoFEMTypes bh=MF_EXIST, int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
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:23
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
Definition: DMMMoFEM.cpp:979
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:400
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition: DMMMoFEM.cpp:794
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:442
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
Definition: DMMoFEM.hpp:781
PetscBool isProblemBuild
True if problem is build.
Definition: DMMoFEM.hpp:747
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:27
const Problem * problemMainOfSubPtr
pinter to main problem to sub-problem
Definition: DMMoFEM.hpp:763
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:53
PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem)
Definition: DMMMoFEM.cpp:358
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
Definition: DMMoFEM.hpp:782
PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx)
get MoFEM::KspCtx data structure
Definition: DMMMoFEM.cpp:768
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:663
virtual MoFEMErrorCode resolve_shared_ents(const Problem *problem_ptr, const std::string &fe_name, int verb=DEFAULT_VERBOSITY)=0
resolve shared entities for finite elements in the problem
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:420
PetscErrorCode KspMat(KSP ksp, Mat A, Mat B, void *ctx)
Run over elenents in the list.
Definition: KspCtx.cpp:59
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:453
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:35
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition: DMMMoFEM.cpp:1034
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMMoFEM.cpp:463
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:594
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition: DMMMoFEM.cpp:1098
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:762
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:765
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:766
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:842
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:1144
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:741
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:704
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405
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:834
static const MOFEMuuid IDD_DMCTX
Definition: DMMoFEM.hpp:726
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:509
MoFEMErrorCode getNumberOfElementsByNameAndPart(MPI_Comm comm, const std::string name, PetscLayout *layout) const
Get number of finite elements by name on processors.
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:473
boost::shared_ptr< std::map< std::string, std::pair< EntityType, EntityType > > > mapTypeCol
Definition: DMMoFEM.hpp:771
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:785
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMMoFEM.cpp:823
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:757
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:1092
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:103
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition: DMMMoFEM.cpp:910
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVectorsets the routine to create a local vector associated with the shell DM.
Definition: DMMMoFEM.cpp:876
PetscErrorCode DMMoFEMSetSnesCtx(DM dm, boost::shared_ptr< MoFEM::SnesCtx > &snes_ctx)
Set MoFEM::SnesCtx data structure.
Definition: DMMMoFEM.cpp:811