v0.14.0
poisson_2d_homogeneous.cpp
Go to the documentation of this file.
1 /**
2  * \file poisson_2d_homogeneous.cpp
3  * \example poisson_2d_homogeneous.cpp
4  *
5  * Solution of poisson equation. Direct implementation of User Data Operators
6  * for teaching proposes.
7  *
8  * \note In practical application we suggest use form integrators to generalise
9  * and simplify code. However, here we like to expose user to ways how to
10  * implement data operator from scratch.
11  */
12 
13 constexpr auto field_name = "U";
14 
15 constexpr int SPACE_DIM =
16  EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
17 
19 
20 using namespace MoFEM;
21 using namespace Poisson2DHomogeneousOperators;
22 
24 
25 static char help[] = "...\n\n";
26 
28 public:
30 
31  // Declaration of the main function to run analysis
32  MoFEMErrorCode runProgram();
33 
34 private:
35  // Declaration of other main functions called in runProgram()
36  MoFEMErrorCode readMesh();
37  MoFEMErrorCode setupProblem();
38  MoFEMErrorCode boundaryCondition();
39  MoFEMErrorCode assembleSystem();
40  MoFEMErrorCode setIntegrationRules();
41  MoFEMErrorCode solveSystem();
42  MoFEMErrorCode outputResults();
43  MoFEMErrorCode checkResults();
44 
45  // MoFEM interfaces
47  // Field name
49  // Approximation order
50  int oRder;
51  // Function to calculate the Source term
52  static double sourceTermFunction(const double x, const double y,
53  const double z) {
54  return 2. * M_PI * M_PI * sin(M_PI * x) * sin(M_PI * y);
55  }
56  // PETSc vector for storing norms
58  int atom_test = 0;
59  enum NORMS { NORM = 0, LAST_NORM };
60 };
61 
63  : mField(m_field) {}
64 
65 //! [Read mesh]
68 
72 
74 }
75 //! [Read mesh]
76 
77 //! [Setup problem]
80  Range domain_ents;
81  CHKERR mField.get_moab().get_entities_by_dimension(0, SPACE_DIM, domain_ents,
82  true);
83  auto get_ents_by_dim = [&](const auto dim) {
84  if (dim == SPACE_DIM) {
85  return domain_ents;
86  } else {
87  Range ents;
88  if (dim == 0)
89  CHKERR mField.get_moab().get_connectivity(domain_ents, ents, true);
90  else
91  CHKERR mField.get_moab().get_entities_by_dimension(0, dim, ents, true);
92  return ents;
93  }
94  };
95 
96  // Select base
97  auto get_base = [&]() {
98  auto domain_ents = get_ents_by_dim(SPACE_DIM);
99  if (domain_ents.empty())
100  CHK_THROW_MESSAGE(MOFEM_NOT_FOUND, "Empty mesh");
101  const auto type = type_from_handle(domain_ents[0]);
102  switch (type) {
103  case MBQUAD:
104  return DEMKOWICZ_JACOBI_BASE;
105  case MBHEX:
106  return DEMKOWICZ_JACOBI_BASE;
107  case MBTRI:
109  case MBTET:
111  default:
112  CHK_THROW_MESSAGE(MOFEM_NOT_FOUND, "Element type not handled");
113  }
114  return NOBASE;
115  };
116  auto base = get_base();
118 
119  int oRder = 3;
120  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &oRder, PETSC_NULL);
121  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-atom_test", &atom_test,
122  PETSC_NULL);
124 
126 
128 }
129 //! [Setup problem]
130 
131 //! [Boundary condition]
134 
135  auto bc_mng = mField.getInterface<BcManager>();
136 
137  // Remove BCs from blockset name "BOUNDARY_CONDITION" or SETU, note that you
138  // can use regular expression to put list of blocksets;
139  CHKERR bc_mng->removeBlockDOFsOnEntities<BcScalarMeshsetType<BLOCKSET>>(
140  simpleInterface->getProblemName(), "(BOUNDARY_CONDITION|SETU)",
141  std::string(field_name), true);
142 
144 }
145 //! [Boundary condition]
146 
147 //! [Assemble system]
150 
151  auto pipeline_mng = mField.getInterface<PipelineManager>();
152 
153  { // Push operators to the Pipeline that is responsible for calculating LHS
155  pipeline_mng->getOpDomainLhsPipeline(), {H1});
156  pipeline_mng->getOpDomainLhsPipeline().push_back(
158  }
159 
160  { // Push operators to the Pipeline that is responsible for calculating RHS
161 
162  auto set_values_to_bc_dofs = [&](auto &fe) {
163  auto get_bc_hook = [&]() {
165  return hook;
166  };
167  fe->preProcessHook = get_bc_hook();
168  };
169 
170  // you can skip that if boundary condition is prescribing zero
171  auto calculate_residual_from_set_values_on_bc = [&](auto &pipeline) {
172  using DomainEle =
175  using OpInternal = FormsIntegrators<DomainEleOp>::Assembly<
177 
178  auto grad_u_vals_ptr = boost::make_shared<MatrixDouble>();
179  pipeline_mng->getOpDomainRhsPipeline().push_back(
181  grad_u_vals_ptr));
182  pipeline_mng->getOpDomainRhsPipeline().push_back(
183  new OpInternal(field_name, grad_u_vals_ptr,
184  [](double, double, double) constexpr { return -1; }));
185  };
186 
188  pipeline_mng->getOpDomainRhsPipeline(), {H1});
189  set_values_to_bc_dofs(pipeline_mng->getDomainRhsFE());
190  calculate_residual_from_set_values_on_bc(
191  pipeline_mng->getOpDomainRhsPipeline());
192  pipeline_mng->getOpDomainRhsPipeline().push_back(
194  }
195 
197 }
198 //! [Assemble system]
199 
200 //! [Set integration rules]
203 
204  auto rule_lhs = [](int, int, int p) -> int { return 2 * (p - 1); };
205  auto rule_rhs = [](int, int, int p) -> int { return p; };
206 
207  auto pipeline_mng = mField.getInterface<PipelineManager>();
208  CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule_lhs);
209  CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule_rhs);
210 
212 }
213 //! [Set integration rules]
214 
215 //! [Solve system]
218 
219  auto pipeline_mng = mField.getInterface<PipelineManager>();
220 
221  auto ksp_solver = pipeline_mng->createKSP();
222  CHKERR KSPSetFromOptions(ksp_solver);
223  CHKERR KSPSetUp(ksp_solver);
224 
225  // Create RHS and solution vectors
226  auto dm = simpleInterface->getDM();
227  auto F = createDMVector(dm);
228  auto D = vectorDuplicate(F);
229 
230  // Solve the system
231  CHKERR KSPSolve(ksp_solver, F, D);
232 
233  // Scatter result data on the mesh
234  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
235  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
236  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
237 
239 }
240 //! [Solve system]
241 
242 //! [Output results]
245 
246  auto pipeline_mng = mField.getInterface<PipelineManager>();
247  pipeline_mng->getDomainLhsFE().reset();
248 
249  auto post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
251  post_proc_fe->getOpPtrVector(), {H1});
252 
253  auto u_ptr = boost::make_shared<VectorDouble>();
254  auto grad_u_ptr = boost::make_shared<MatrixDouble>();
255  post_proc_fe->getOpPtrVector().push_back(
257 
258  post_proc_fe->getOpPtrVector().push_back(
260 
262 
263  post_proc_fe->getOpPtrVector().push_back(
264 
265  new OpPPMap(post_proc_fe->getPostProcMesh(),
266  post_proc_fe->getMapGaussPts(),
267 
268  OpPPMap::DataMapVec{{"U", u_ptr}},
269 
270  OpPPMap::DataMapMat{{"GRAD_U", grad_u_ptr}},
271 
273 
275 
276  )
277 
278  );
279 
280  pipeline_mng->getDomainRhsFE() = post_proc_fe;
281  CHKERR pipeline_mng->loopFiniteElements();
282  CHKERR post_proc_fe->writeFile("out_result.h5m");
283 
285 }
286 //! [Output results]
287 
288 //! [Check]
291 
292  auto check_result_fe_ptr = boost::make_shared<DomainEle>(mField);
293  auto petscVec =
295  (mField.get_comm_rank() == 0) ? LAST_NORM : 0, LAST_NORM);
296 
298  check_result_fe_ptr->getOpPtrVector(), {H1})),
299  "Apply transform");
300 
301  check_result_fe_ptr->getRuleHook = [](int, int, int p) { return p; };
302  auto analyticalFunction = [&](double x, double y, double z) {
303  return sin(M_PI * x) * sin(M_PI * y);
304  };
305 
306  auto u_ptr = boost::make_shared<VectorDouble>();
307 
308  check_result_fe_ptr->getOpPtrVector().push_back(
310  auto mValFuncPtr = boost::make_shared<VectorDouble>();
311  check_result_fe_ptr->getOpPtrVector().push_back(
312  new OpGetTensor0fromFunc(mValFuncPtr, analyticalFunction));
313  check_result_fe_ptr->getOpPtrVector().push_back(
314  new OpCalcNormL2Tensor0(u_ptr, petscVec, NORM, mValFuncPtr));
315  CHKERR VecZeroEntries(petscVec);
318  check_result_fe_ptr);
319  CHKERR VecAssemblyBegin(petscVec);
320  CHKERR VecAssemblyEnd(petscVec);
321  MOFEM_LOG_CHANNEL("SELF"); // Clear channel from old tags
322  // print norm in general log
323  if (mField.get_comm_rank() == 0) {
324  const double *norms;
325  CHKERR VecGetArrayRead(petscVec, &norms);
326  MOFEM_TAG_AND_LOG("SELF", Sev::inform, "Errors")
327  << "NORM: " << std::sqrt(norms[NORM]);
328  CHKERR VecRestoreArrayRead(petscVec, &norms);
329  }
330  // compare norm for ctest
331  if (atom_test && !mField.get_comm_rank()) {
332  const double *t_ptr;
333  CHKERR VecGetArrayRead(petscVec, &t_ptr);
334  double ref_norm = 2.2e-04;
335  double cal_norm;
336  switch (atom_test) {
337  case 1: // 2D
338  cal_norm = sqrt(t_ptr[0]);
339  break;
340  default:
341  SETERRQ1(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
342  "atom test %d does not exist", atom_test);
343  }
344  if (cal_norm > ref_norm) {
345  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
346  "atom test %d failed! Calculated Norm %3.16e is greater than "
347  "reference Norm %3.16e",
348  atom_test, cal_norm, ref_norm);
349  }
350  CHKERR VecRestoreArrayRead(petscVec, &t_ptr);
351  }
353 }
354 //! [Check]
355 
356 //! [Run program]
359 
360  CHKERR readMesh();
368 
370 }
371 //! [Run program]
372 
373 //! [Main]
374 int main(int argc, char *argv[]) {
375 
376  // Initialisation of MoFEM/PETSc and MOAB data structures
377  const char param_file[] = "param_file.petsc";
378  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
379 
380  // Error handling
381  try {
382  // Register MoFEM discrete manager in PETSc
383  DMType dm_name = "DMMOFEM";
384  CHKERR DMRegister_MoFEM(dm_name);
385 
386  // Create MOAB instance
387  moab::Core mb_instance; // mesh database
388  moab::Interface &moab = mb_instance; // mesh database interface
389 
390  // Create MoFEM instance
391  MoFEM::Core core(moab); // finite element database
392  MoFEM::Interface &m_field = core; // finite element interface
393 
394  // Run the main analysis
395  Poisson2DHomogeneous poisson_problem(m_field);
396  CHKERR poisson_problem.runProgram();
397  }
398  CATCH_ERRORS;
399 
400  // Finish work: cleaning memory, getting statistics, etc.
402 
403  return 0;
404 }
405 //! [Main]
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
Poisson2DHomogeneous::solveSystem
MoFEMErrorCode solveSystem()
[Set integration rules]
Definition: poisson_2d_homogeneous.cpp:216
EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION
Definition: plastic.cpp:13
Poisson2DHomogeneous::NORMS
NORMS
Definition: poisson_2d_homogeneous.cpp:59
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
Poisson2DHomogeneous::sourceTermFunction
static double sourceTermFunction(const double x, const double y, const double z)
Definition: poisson_2d_homogeneous.cpp:52
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
Poisson2DHomogeneousOperators::OpDomainRhsVectorF
Definition: poisson_2d_homogeneous.hpp:94
Poisson2DHomogeneous::oRder
int oRder
Definition: poisson_2d_homogeneous.cpp:50
MoFEM::PipelineManager::ElementsAndOpsByDim
Definition: PipelineManager.hpp:38
Poisson2DHomogeneous::setupProblem
MoFEMErrorCode setupProblem()
[Read mesh]
Definition: poisson_2d_homogeneous.cpp:78
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
NOBASE
@ NOBASE
Definition: definitions.h:59
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
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::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
Poisson2DHomogeneous::assembleSystem
MoFEMErrorCode assembleSystem()
[Boundary condition]
Definition: poisson_2d_homogeneous.cpp:148
Poisson2DHomogeneous::outputResults
MoFEMErrorCode outputResults()
[Solve system]
Definition: poisson_2d_homogeneous.cpp:243
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
SPACE_DIM
constexpr int SPACE_DIM
Definition: poisson_2d_homogeneous.cpp:15
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Poisson2DHomogeneous::Poisson2DHomogeneous
Poisson2DHomogeneous(MoFEM::Interface &m_field)
Definition: poisson_2d_homogeneous.cpp:62
Poisson2DHomogeneous::simpleInterface
Simple * simpleInterface
Definition: poisson_2d_homogeneous.cpp:46
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
atom_test
int atom_test
Definition: contact.cpp:97
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::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
convert.type
type
Definition: convert.py:64
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::BcScalarMeshsetType
Definition: BcManager.hpp:15
OpDomainLhsMatrixK
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< BASE_DIM, FIELD_DIM, SPACE_DIM > OpDomainLhsMatrixK
Definition: electrostatics.hpp:40
Poisson2DHomogeneous::petscVec
SmartPetscObj< Vec > petscVec
Definition: poisson_2d_homogeneous.cpp:57
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
Poisson2DHomogeneous::atom_test
int atom_test
Definition: poisson_2d_homogeneous.cpp:58
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::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::Simple::getDomainFEName
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition: Simple.hpp:368
Poisson2DHomogeneous
Definition: poisson_2d_homogeneous.cpp:27
MoFEM::type_from_handle
auto type_from_handle(const EntityHandle h)
get type from entity handle
Definition: Templates.hpp:1898
Poisson2DHomogeneous::mField
MoFEM::Interface & mField
Definition: poisson_2d_homogeneous.cpp:48
Poisson2DHomogeneous::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
[Assemble system]
Definition: poisson_2d_homogeneous.cpp:201
Poisson2DHomogeneous::LAST_NORM
@ LAST_NORM
Definition: poisson_2d_homogeneous.cpp:59
MoFEM::OpGetTensor0fromFunc
Get values from scalar function at integration points and save them to VectorDouble for Tensor0.
Definition: NormsOperators.hpp:105
MoFEM::OpPostProcMapInMoab::DataMapMat
std::map< std::string, boost::shared_ptr< MatrixDouble > > DataMapMat
Definition: PostProcBrokenMeshInMoabBase.hpp:701
main
int main(int argc, char *argv[])
[Run program]
Definition: poisson_2d_homogeneous.cpp:374
MOFEM_NOT_FOUND
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
help
static char help[]
Definition: poisson_2d_homogeneous.cpp:25
Poisson2DHomogeneous::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Setup problem]
Definition: poisson_2d_homogeneous.cpp:132
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:413
Poisson2DHomogeneous::checkResults
MoFEMErrorCode checkResults()
[Output results]
Definition: poisson_2d_homogeneous.cpp:289
MoFEM::EssentialPreProc< TemperatureCubitBcData >
Specialization for TemperatureCubitBcData.
Definition: EssentialTemperatureCubitBcData.hpp:91
Range
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
DomainEleOp
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
Poisson2DHomogeneous::runProgram
MoFEMErrorCode runProgram()
[Check]
Definition: poisson_2d_homogeneous.cpp:357
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
Poisson2DHomogeneousOperators
Definition: poisson_2d_homogeneous.hpp:31
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
poisson_2d_homogeneous.hpp
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::OpCalcNormL2Tensor0
Get norm of input VectorDouble for Tensor0.
Definition: NormsOperators.hpp:15
MoFEM::createVectorMPI
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
Definition: PetscSmartObj.hpp:202
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
Poisson2DHomogeneous::readMesh
MoFEMErrorCode readMesh()
[Read mesh]
Definition: poisson_2d_homogeneous.cpp:66
MoFEM::SmartPetscObj< Vec >
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:683
MoFEM::Simple::getProblemName
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:389
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:586
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
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
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
Poisson2DHomogeneous::NORM
@ NORM
Definition: poisson_2d_homogeneous.cpp:59
F
@ F
Definition: free_surface.cpp:394
MoFEM::PostProcBrokenMeshInMoab< FaceElementForcesAndSourcesCore >
Definition: PostProcBrokenMeshInMoabBase.hpp:677
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698