v0.14.0
Classes | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
MixedPoisson Struct Reference
Collaboration diagram for MixedPoisson:
[legend]

Classes

struct  CommonData
 
struct  OpError
 

Public Member Functions

 MixedPoisson (MoFEM::Interface &m_field)
 
MoFEMErrorCode runProblem ()
 [Run programme] More...
 

Private Member Functions

MoFEMErrorCode readMesh ()
 [Run programme] More...
 
MoFEMErrorCode setupProblem ()
 [Read mesh] More...
 
MoFEMErrorCode setIntegrationRules ()
 [Set up problem] More...
 
MoFEMErrorCode createCommonData ()
 [Set integration rule] More...
 
MoFEMErrorCode assembleSystem ()
 [Create common data] More...
 
MoFEMErrorCode solveSystem ()
 [Assemble system] More...
 
MoFEMErrorCode checkError (int iter_num=0)
 [Solve and refine loop] More...
 
MoFEMErrorCode refineOrder (int iter_num=0)
 [Solve] More...
 
MoFEMErrorCode solveRefineLoop ()
 [Refine] More...
 
MoFEMErrorCode outputResults (int iter_num=0)
 [Check error] More...
 

Static Private Member Functions

static double analyticalFunction (const double x, const double y, const double z)
 [Analytical function] More...
 
static VectorDouble analyticalFunctionGrad (const double x, const double y, const double z)
 [Analytical function] More...
 
static double sourceFunction (const double x, const double y, const double z)
 [Analytical function gradient] More...
 
static MoFEMErrorCode getTagHandle (MoFEM::Interface &m_field, const char *name, DataType type, Tag &tag_handle)
 [Source function] More...
 

Private Attributes

MoFEM::InterfacemField
 
SimplesimpleInterface
 
Range domainEntities
 
double totErrorIndicator
 
double maxErrorIndicator
 
double thetaParam
 
double indicTolerance
 
int initOrder
 
boost::shared_ptr< CommonDatacommonDataPtr
 

Detailed Description

Examples
mixed_poisson.cpp.

Definition at line 31 of file mixed_poisson.cpp.

Constructor & Destructor Documentation

◆ MixedPoisson()

MixedPoisson::MixedPoisson ( MoFEM::Interface m_field)
inline

Definition at line 33 of file mixed_poisson.cpp.

33 : mField(m_field) {}

Member Function Documentation

◆ analyticalFunction()

static double MixedPoisson::analyticalFunction ( const double  x,
const double  y,
const double  z 
)
inlinestaticprivate

[Analytical function]

Examples
mixed_poisson.cpp.

Definition at line 49 of file mixed_poisson.cpp.

50  {
51  return exp(-100. * (sqr(x) + sqr(y))) * cos(M_PI * x) * cos(M_PI * y);
52  }

◆ analyticalFunctionGrad()

static VectorDouble MixedPoisson::analyticalFunctionGrad ( const double  x,
const double  y,
const double  z 
)
inlinestaticprivate

[Analytical function]

[Analytical function gradient]

Examples
mixed_poisson.cpp.

Definition at line 56 of file mixed_poisson.cpp.

57  {
58  VectorDouble res;
59  res.resize(2);
60  res[0] = -exp(-100. * (sqr(x) + sqr(y))) *
61  (200. * x * cos(M_PI * x) + M_PI * sin(M_PI * x)) * cos(M_PI * y);
62  res[1] = -exp(-100. * (sqr(x) + sqr(y))) *
63  (200. * y * cos(M_PI * y) + M_PI * sin(M_PI * y)) * cos(M_PI * x);
64  return res;
65  }

◆ assembleSystem()

MoFEMErrorCode MixedPoisson::assembleSystem ( )
private

[Create common data]

[Assemble system]

Examples
mixed_poisson.cpp.

Definition at line 239 of file mixed_poisson.cpp.

239  {
242  pipeline_mng->getDomainLhsFE().reset();
243  pipeline_mng->getDomainRhsFE().reset();
244  pipeline_mng->getOpDomainRhsPipeline().clear();
245  pipeline_mng->getOpDomainLhsPipeline().clear();
246 
247  CHKERR AddHOOps<2, 2, 2>::add(pipeline_mng->getOpDomainLhsPipeline(), {HDIV});
248 
249  auto beta = [](const double, const double, const double) { return 1; };
250  pipeline_mng->getOpDomainLhsPipeline().push_back(
251  new OpHdivHdiv("Q", "Q", beta));
252  auto unity = []() { return 1; };
253  pipeline_mng->getOpDomainLhsPipeline().push_back(
254  new OpHdivU("Q", "U", unity, true));
255  auto source = [&](const double x, const double y, const double z) {
256  return -sourceFunction(x, y, z);
257  };
258  pipeline_mng->getOpDomainRhsPipeline().push_back(
259  new OpDomainSource("U", source));
261 }

◆ checkError()

MoFEMErrorCode MixedPoisson::checkError ( int  iter_num = 0)
private

[Solve and refine loop]

[Check error]

Examples
mixed_poisson.cpp.

Definition at line 371 of file mixed_poisson.cpp.

371  {
374  pipeline_mng->getDomainLhsFE().reset();
375  pipeline_mng->getDomainRhsFE().reset();
376  pipeline_mng->getOpDomainRhsPipeline().clear();
377 
379  {HDIV, L2});
380 
381  pipeline_mng->getOpDomainRhsPipeline().push_back(
382  new OpCalculateScalarFieldValues("U", commonDataPtr->approxVals));
383  pipeline_mng->getOpDomainRhsPipeline().push_back(
385  commonDataPtr->approxValsGrad));
386  pipeline_mng->getOpDomainRhsPipeline().push_back(
387  new OpCalculateHVecVectorField<3>("Q", commonDataPtr->approxFlux));
388  pipeline_mng->getOpDomainRhsPipeline().push_back(
390 
391  commonDataPtr->maxErrorIndicator = 0;
392  CHKERR VecZeroEntries(commonDataPtr->petscVec);
393  CHKERR pipeline_mng->loopFiniteElements();
394  CHKERR VecAssemblyBegin(commonDataPtr->petscVec);
395  CHKERR VecAssemblyEnd(commonDataPtr->petscVec);
396  CHKERR VecGhostUpdateBegin(commonDataPtr->petscVec, INSERT_VALUES,
397  SCATTER_FORWARD);
398  CHKERR VecGhostUpdateEnd(commonDataPtr->petscVec, INSERT_VALUES,
399  SCATTER_FORWARD);
400 
401  MPI_Allreduce(&commonDataPtr->maxErrorIndicator, &maxErrorIndicator, 1,
402  MPI_DOUBLE, MPI_MAX, mField.get_comm());
403 
404  const double *array;
405  CHKERR VecGetArrayRead(commonDataPtr->petscVec, &array);
406  MOFEM_LOG("EXAMPLE", Sev::inform)
407  << "Global error indicator (max): " << commonDataPtr->maxErrorIndicator;
408  MOFEM_LOG("EXAMPLE", Sev::inform)
409  << "Global error indicator (total): "
410  << std::sqrt(array[CommonData::ERROR_INDICATOR_TOTAL]);
411  MOFEM_LOG("EXAMPLE", Sev::inform)
412  << "Global error L2 norm: "
413  << std::sqrt(array[CommonData::ERROR_L2_NORM]);
414  MOFEM_LOG("EXAMPLE", Sev::inform)
415  << "Global error H1 seminorm: "
416  << std::sqrt(array[CommonData::ERROR_H1_SEMINORM]);
417 
419  CHKERR VecRestoreArrayRead(commonDataPtr->petscVec, &array);
420 
421  std::vector<Tag> tag_handles;
422  tag_handles.resize(4);
423  CHKERR getTagHandle(mField, "ERROR_L2_NORM", MB_TYPE_DOUBLE, tag_handles[0]);
424  CHKERR getTagHandle(mField, "ERROR_H1_SEMINORM", MB_TYPE_DOUBLE,
425  tag_handles[1]);
426  CHKERR getTagHandle(mField, "ERROR_INDICATOR", MB_TYPE_DOUBLE,
427  tag_handles[2]);
428  CHKERR getTagHandle(mField, "ORDER", MB_TYPE_INTEGER, tag_handles[3]);
429 
430  ParallelComm *pcomm =
431  ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
432  if (pcomm == NULL)
433  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "Communicator not set");
434 
435  tag_handles.push_back(pcomm->part_tag());
436  std::ostringstream strm;
437  strm << "error_" << iter_num << ".h5m";
438  CHKERR mField.get_moab().write_file(strm.str().c_str(), "MOAB",
439  "PARALLEL=WRITE_PART", 0, 0,
440  tag_handles.data(), tag_handles.size());
442 }

◆ createCommonData()

MoFEMErrorCode MixedPoisson::createCommonData ( )
private

[Set integration rule]

[Create common data]

Examples
mixed_poisson.cpp.

Definition at line 221 of file mixed_poisson.cpp.

221  {
223  commonDataPtr = boost::make_shared<CommonData>();
224  PetscInt ghosts[3] = {0, 1, 2};
225  if (!mField.get_comm_rank())
226  commonDataPtr->petscVec =
227  createGhostVector(mField.get_comm(), 3, 3, 0, ghosts);
228  else
229  commonDataPtr->petscVec =
230  createGhostVector(mField.get_comm(), 0, 3, 3, ghosts);
231  commonDataPtr->approxVals = boost::make_shared<VectorDouble>();
232  commonDataPtr->approxValsGrad = boost::make_shared<MatrixDouble>();
233  commonDataPtr->approxFlux = boost::make_shared<MatrixDouble>();
235 }

◆ getTagHandle()

static MoFEMErrorCode MixedPoisson::getTagHandle ( MoFEM::Interface m_field,
const char *  name,
DataType  type,
Tag &  tag_handle 
)
inlinestaticprivate

[Source function]

Examples
mixed_poisson.cpp.

Definition at line 79 of file mixed_poisson.cpp.

81  {
83  int int_val = 0;
84  double double_val = 0;
85  switch (type) {
86  case MB_TYPE_INTEGER:
87  CHKERR m_field.get_moab().tag_get_handle(
88  name, 1, type, tag_handle, MB_TAG_CREAT | MB_TAG_SPARSE, &int_val);
89  break;
90  case MB_TYPE_DOUBLE:
91  CHKERR m_field.get_moab().tag_get_handle(
92  name, 1, type, tag_handle, MB_TAG_CREAT | MB_TAG_SPARSE, &double_val);
93  break;
94  default:
95  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
96  "Wrong data type %d for tag", type);
97  }
99  }

◆ outputResults()

MoFEMErrorCode MixedPoisson::outputResults ( int  iter_num = 0)
private

[Check error]

[Output results]

Examples
mixed_poisson.cpp.

Definition at line 446 of file mixed_poisson.cpp.

446  {
449  pipeline_mng->getDomainLhsFE().reset();
450 
452 
453  auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
454 
455  CHKERR AddHOOps<2, 2, 2>::add(post_proc_fe->getOpPtrVector(), {HDIV});
456 
457  auto u_ptr = boost::make_shared<VectorDouble>();
458  auto flux_ptr = boost::make_shared<MatrixDouble>();
459  post_proc_fe->getOpPtrVector().push_back(
460  new OpCalculateScalarFieldValues("U", u_ptr));
461  post_proc_fe->getOpPtrVector().push_back(
462  new OpCalculateHVecVectorField<3>("Q", flux_ptr));
463 
465 
466  post_proc_fe->getOpPtrVector().push_back(
467 
468  new OpPPMap(post_proc_fe->getPostProcMesh(),
469  post_proc_fe->getMapGaussPts(),
470 
471  OpPPMap::DataMapVec{{"U", u_ptr}},
472 
473  OpPPMap::DataMapMat{{"Q", flux_ptr}},
474 
476 
478 
479  )
480 
481  );
482 
483  pipeline_mng->getDomainRhsFE() = post_proc_fe;
484  CHKERR pipeline_mng->loopFiniteElements();
485 
486  std::ostringstream strm;
487  strm << "out_" << iter_num << ".h5m";
488  CHKERR post_proc_fe->writeFile(strm.str().c_str());
490 }

◆ readMesh()

MoFEMErrorCode MixedPoisson::readMesh ( )
private

[Run programme]

[Read mesh]

Examples
mixed_poisson.cpp.

Definition at line 155 of file mixed_poisson.cpp.

◆ refineOrder()

MoFEMErrorCode MixedPoisson::refineOrder ( int  iter_num = 0)
private

[Solve]

[Refine]

Examples
mixed_poisson.cpp.

Definition at line 286 of file mixed_poisson.cpp.

286  {
288  Tag th_error_ind, th_order;
289  CHKERR getTagHandle(mField, "ERROR_INDICATOR", MB_TYPE_DOUBLE, th_error_ind);
290  CHKERR getTagHandle(mField, "ORDER", MB_TYPE_INTEGER, th_order);
291 
292  std::vector<Range> refinement_levels;
293  refinement_levels.resize(iter_num + 1);
294  for (auto ent : domainEntities) {
295  double err_indic = 0;
296  CHKERR mField.get_moab().tag_get_data(th_error_ind, &ent, 1, &err_indic);
297 
298  int order, new_order;
299  CHKERR mField.get_moab().tag_get_data(th_order, &ent, 1, &order);
300  new_order = order + 1;
301  Range refined_ents;
302  if (err_indic > thetaParam * maxErrorIndicator) {
303  refined_ents.insert(ent);
304  Range adj;
305  CHKERR mField.get_moab().get_adjacencies(&ent, 1, 1, false, adj,
306  moab::Interface::UNION);
307  refined_ents.merge(adj);
308  refinement_levels[new_order - initOrder].merge(refined_ents);
309  CHKERR mField.get_moab().tag_set_data(th_order, &ent, 1, &new_order);
310  }
311  }
312 
313  for (int ll = 1; ll < refinement_levels.size(); ll++) {
314  CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
315  refinement_levels[ll]);
316 
317  if (initOrder + ll > 8) {
318  MOFEM_LOG("EXAMPLE", Sev::warning)
319  << "setting approximation order higher than 8 is not currently "
320  "supported"
321  << endl;
322  } else {
323  CHKERR mField.set_field_order(refinement_levels[ll], "U", initOrder + ll);
324  CHKERR mField.set_field_order(refinement_levels[ll], "Q",
325  initOrder + ll + 1);
326  }
327  }
328 
329  CHKERR mField.getInterface<CommInterface>()->synchroniseFieldEntities("Q");
330  CHKERR mField.getInterface<CommInterface>()->synchroniseFieldEntities("U");
334  mField.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
337 }

◆ runProblem()

MoFEMErrorCode MixedPoisson::runProblem ( )

[Run programme]

Examples
mixed_poisson.cpp.

Definition at line 144 of file mixed_poisson.cpp.

144  {
146  CHKERR readMesh();
151 }

◆ setIntegrationRules()

MoFEMErrorCode MixedPoisson::setIntegrationRules ( )
private

[Set up problem]

[Set integration rule]

Examples
mixed_poisson.cpp.

Definition at line 207 of file mixed_poisson.cpp.

207  {
209 
210  auto rule = [](int, int, int p) -> int { return 2 * p; };
211 
213  CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
214  CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
215 
217 }

◆ setupProblem()

MoFEMErrorCode MixedPoisson::setupProblem ( )
private

[Read mesh]

[Set up problem]

Examples
mixed_poisson.cpp.

Definition at line 165 of file mixed_poisson.cpp.

165  {
167 
169 
170  int nb_quads = 0;
171  CHKERR mField.get_moab().get_number_entities_by_type(0, MBQUAD, nb_quads);
172  auto base = AINSWORTH_LEGENDRE_BASE;
173  if (nb_quads) {
174  // AINSWORTH_LEGENDRE_BASE is not implemented for HDIV/HCURL space on quads
175  base = DEMKOWICZ_JACOBI_BASE;
176  }
177 
178  // Note that in 2D case HDIV and HCURL spaces are isomorphic, and therefore
179  // only base for HCURL has been implemented in 2D. Base vectors for HDIV space
180  // are be obtained after rotation of HCURL base vectors by a right angle
181  CHKERR simpleInterface->addDomainField("Q", HCURL, base, 1);
182 
183  thetaParam = 0.5;
184  CHKERR PetscOptionsGetReal(PETSC_NULL, "", "-theta", &thetaParam, PETSC_NULL);
185 
186  indicTolerance = 1e-3;
187  CHKERR PetscOptionsGetReal(PETSC_NULL, "", "-indic_tol", &indicTolerance,
188  PETSC_NULL);
189 
190  initOrder = 2;
191  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &initOrder, PETSC_NULL);
195 
196  CHKERR mField.get_moab().get_entities_by_dimension(0, 2, domainEntities);
197  Tag th_order;
198  CHKERR getTagHandle(mField, "ORDER", MB_TYPE_INTEGER, th_order);
199  for (auto ent : domainEntities) {
200  CHKERR mField.get_moab().tag_set_data(th_order, &ent, 1, &initOrder);
201  }
203 }

◆ solveRefineLoop()

MoFEMErrorCode MixedPoisson::solveRefineLoop ( )
private

[Refine]

[Solve and refine loop]

Examples
mixed_poisson.cpp.

Definition at line 341 of file mixed_poisson.cpp.

341  {
346  CHKERR checkError();
348 
349  int iter_num = 1;
350  while (fabs(indicTolerance) > DBL_EPSILON &&
352  MOFEM_LOG("EXAMPLE", Sev::inform) << "Refinement iteration " << iter_num;
353 
354  CHKERR refineOrder(iter_num);
358  CHKERR checkError(iter_num);
359  CHKERR outputResults(iter_num);
360 
361  iter_num++;
362  if (iter_num > 100)
363  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
364  "Too many refinement iterations");
365  }
367 }

◆ solveSystem()

MoFEMErrorCode MixedPoisson::solveSystem ( )
private

[Assemble system]

[Solve]

Examples
mixed_poisson.cpp.

Definition at line 265 of file mixed_poisson.cpp.

265  {
268  auto solver = pipeline_mng->createKSP();
269  CHKERR KSPSetFromOptions(solver);
270  CHKERR KSPSetUp(solver);
271 
272  auto dm = simpleInterface->getDM();
273  auto D = createDMVector(dm);
274  auto F = vectorDuplicate(D);
275  CHKERR VecZeroEntries(D);
276 
277  CHKERR KSPSolve(solver, F, D);
278  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
279  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
280  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
282 }

◆ sourceFunction()

static double MixedPoisson::sourceFunction ( const double  x,
const double  y,
const double  z 
)
inlinestaticprivate

[Analytical function gradient]

[Source function]

Definition at line 69 of file mixed_poisson.cpp.

69  {
70  return -exp(-100. * (sqr(x) + sqr(y))) *
71  (400. * M_PI *
72  (x * cos(M_PI * y) * sin(M_PI * x) +
73  y * cos(M_PI * x) * sin(M_PI * y)) +
74  2. * (20000. * (sqr(x) + sqr(y)) - 200. - sqr(M_PI)) *
75  cos(M_PI * x) * cos(M_PI * y));
76  }

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<CommonData> MixedPoisson::commonDataPtr
private

Definition at line 128 of file mixed_poisson.cpp.

◆ domainEntities

Range MixedPoisson::domainEntities
private

Definition at line 40 of file mixed_poisson.cpp.

◆ indicTolerance

double MixedPoisson::indicTolerance
private

Definition at line 45 of file mixed_poisson.cpp.

◆ initOrder

int MixedPoisson::initOrder
private

Definition at line 46 of file mixed_poisson.cpp.

◆ maxErrorIndicator

double MixedPoisson::maxErrorIndicator
private

Definition at line 42 of file mixed_poisson.cpp.

◆ mField

MoFEM::Interface& MixedPoisson::mField
private

Definition at line 37 of file mixed_poisson.cpp.

◆ simpleInterface

Simple* MixedPoisson::simpleInterface
private

Definition at line 38 of file mixed_poisson.cpp.

◆ thetaParam

double MixedPoisson::thetaParam
private

Definition at line 44 of file mixed_poisson.cpp.

◆ totErrorIndicator

double MixedPoisson::totErrorIndicator
private

Definition at line 41 of file mixed_poisson.cpp.


The documentation for this struct was generated from the following file:
MixedPoisson::domainEntities
Range domainEntities
Definition: mixed_poisson.cpp:40
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::PipelineManager::getDomainRhsFE
boost::shared_ptr< FEMethod > & getDomainRhsFE()
Definition: PipelineManager.hpp:405
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
OpError
Definition: initial_diffusion.cpp:61
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::OpPostProcMapInMoab::DataMapVec
std::map< std::string, boost::shared_ptr< VectorDouble > > DataMapVec
Definition: PostProcBrokenMeshInMoabBase.hpp:700
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::PipelineManager::loopFiniteElements
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
Definition: PipelineManager.cpp:19
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM::DMSetUp_MoFEM
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMoFEM.cpp:1303
MixedPoisson::readMesh
MoFEMErrorCode readMesh()
[Run programme]
Definition: mixed_poisson.cpp:155
MixedPoisson::CommonData::ERROR_INDICATOR_TOTAL
@ ERROR_INDICATOR_TOTAL
Definition: mixed_poisson.cpp:123
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:523
MixedPoisson::initOrder
int initOrder
Definition: mixed_poisson.cpp:46
order
constexpr int order
Definition: dg_projection.cpp:18
MixedPoisson::outputResults
MoFEMErrorCode outputResults(int iter_num=0)
[Check error]
Definition: mixed_poisson.cpp:446
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MixedPoisson::CommonData::ERROR_L2_NORM
@ ERROR_L2_NORM
Definition: mixed_poisson.cpp:121
MixedPoisson::mField
MoFEM::Interface & mField
Definition: mixed_poisson.cpp:37
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
MixedPoisson::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: mixed_poisson.cpp:128
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:747
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::createGhostVector
auto createGhostVector(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[])
Create smart ghost vector.
Definition: PetscSmartObj.hpp:179
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MixedPoisson::sourceFunction
static double sourceFunction(const double x, const double y, const double z)
[Analytical function gradient]
Definition: mixed_poisson.cpp:69
MixedPoisson::indicTolerance
double indicTolerance
Definition: mixed_poisson.cpp:45
MixedPoisson::refineOrder
MoFEMErrorCode refineOrder(int iter_num=0)
[Solve]
Definition: mixed_poisson.cpp:286
MixedPoisson::setupProblem
MoFEMErrorCode setupProblem()
[Read mesh]
Definition: mixed_poisson.cpp:165
MixedPoisson::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
[Set up problem]
Definition: mixed_poisson.cpp:207
MoFEM::PetscOptionsGetReal
PetscErrorCode PetscOptionsGetReal(PetscOptions *, const char pre[], const char name[], PetscReal *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:152
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MixedPoisson::solveSystem
MoFEMErrorCode solveSystem()
[Assemble system]
Definition: mixed_poisson.cpp:265
MoFEM::OpCalculateHVecVectorField
Get vector field for H-div approximation.
Definition: UserDataOperators.hpp:2115
double
convert.type
type
Definition: convert.py:64
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::Simple::addDomainField
MoFEMErrorCode addDomainField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEM::PipelineManager::createKSP
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
Definition: PipelineManager.cpp:49
MoFEM::PipelineManager::setDomainRhsIntegrationRule
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:530
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
MoFEM::OpPostProcMapInMoab::DataMapMat
std::map< std::string, boost::shared_ptr< MatrixDouble > > DataMapMat
Definition: PostProcBrokenMeshInMoabBase.hpp:701
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:545
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:503
MoFEM::PipelineManager::setDomainLhsIntegrationRule
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:503
MixedPoisson::maxErrorIndicator
double maxErrorIndicator
Definition: mixed_poisson.cpp:42
Range
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
OpHdivHdiv
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 3, 2 > OpHdivHdiv
Definition: mixed_poisson.cpp:23
MixedPoisson::simpleInterface
Simple * simpleInterface
Definition: mixed_poisson.cpp:38
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
OpHdivU
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< 2 > OpHdivU
Definition: mixed_poisson.cpp:25
MixedPoisson::thetaParam
double thetaParam
Definition: mixed_poisson.cpp:44
MixedPoisson::totErrorIndicator
double totErrorIndicator
Definition: mixed_poisson.cpp:41
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::Simple::getBitRefLevel
BitRefLevel & getBitRefLevel()
Get the BitRefLevel.
Definition: Simple.hpp:354
MixedPoisson::solveRefineLoop
MoFEMErrorCode solveRefineLoop()
[Refine]
Definition: mixed_poisson.cpp:341
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MixedPoisson::checkError
MoFEMErrorCode checkError(int iter_num=0)
[Solve and refine loop]
Definition: mixed_poisson.cpp:371
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
sqr
double sqr(double x)
Definition: mixed_poisson.cpp:29
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MixedPoisson::assembleSystem
MoFEMErrorCode assembleSystem()
[Create common data]
Definition: mixed_poisson.cpp:239
source
auto source
Definition: poisson_2d_dis_galerkin.cpp:61
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:683
MixedPoisson::CommonData::ERROR_H1_SEMINORM
@ ERROR_H1_SEMINORM
Definition: mixed_poisson.cpp:122
MixedPoisson::createCommonData
MoFEMErrorCode createCommonData()
[Set integration rule]
Definition: mixed_poisson.cpp:221
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpDomainSource
Definition: mixed_poisson.cpp:27
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MixedPoisson::getTagHandle
static MoFEMErrorCode getTagHandle(MoFEM::Interface &m_field, const char *name, DataType type, Tag &tag_handle)
[Source function]
Definition: mixed_poisson.cpp:79
F
@ F
Definition: free_surface.cpp:394
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698