v0.14.0
helmholtz.cpp
Go to the documentation of this file.
1 /**
2  * \file helmholtz.cpp
3  * \example helmholtz.cpp
4  *
5  * Using PipelineManager interface calculate Helmholtz problem. Example show how
6  * to manage complex variable fields.
7  */
8 
9 #include <MoFEM.hpp>
10 
11 using namespace MoFEM;
12 
13 static char help[] = "...\n\n";
14 
15 
16 
21 
22 constexpr int SPACE_DIM = 2;
23 
31  PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
32 
33 struct Example {
34 
35  Example(MoFEM::Interface &m_field) : mField(m_field) {}
36 
37  MoFEMErrorCode runProblem();
38 
39 private:
40  MoFEM::Interface &mField;
42  boost::shared_ptr<std::vector<unsigned char>> boundaryMarker;
43 
44  MoFEMErrorCode readMesh();
45  MoFEMErrorCode setupProblem();
46  MoFEMErrorCode boundaryCondition();
47  MoFEMErrorCode assembleSystem();
48  MoFEMErrorCode solveSystem();
49  MoFEMErrorCode outputResults();
50  MoFEMErrorCode checkResults();
51 };
52 
53 //! [run problem]
56  CHKERR readMesh();
57  CHKERR setupProblem();
58  CHKERR boundaryCondition();
59  CHKERR assembleSystem();
60  CHKERR solveSystem();
61  CHKERR outputResults();
62  CHKERR checkResults();
64 }
65 //! [run problem]
66 
67 //! [Read mesh]
70 
71  CHKERR mField.getInterface(simpleInterface);
72  CHKERR simpleInterface->getOptions();
73  CHKERR simpleInterface->loadFile();
74 
76 }
77 //! [Read mesh]
78 
79 //! [Set up problem]
82  // Add field
83  CHKERR simpleInterface->addDomainField("P_REAL", H1,
85  CHKERR simpleInterface->addDomainField("P_IMAG", H1,
87  CHKERR simpleInterface->addBoundaryField("P_REAL", H1,
89  CHKERR simpleInterface->addBoundaryField("P_IMAG", H1,
91  int order = 6;
92  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
93  CHKERR simpleInterface->setFieldOrder("P_REAL", order);
94  CHKERR simpleInterface->setFieldOrder("P_IMAG", order);
95  CHKERR simpleInterface->setUp();
97 }
98 //! [Set up problem]
99 
100 //! [Applying essential BC]
103 
104  auto get_ents_on_mesh_skin = [&]() {
105  Range boundary_entities;
107  std::string entity_name = it->getName();
108  if (entity_name.compare(0, 2, "BC") == 0) {
109  CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 1,
110  boundary_entities, true);
111  }
112  }
113  // Add vertices to boundary entities
114  Range boundary_vertices;
115  CHKERR mField.get_moab().get_connectivity(boundary_entities,
116  boundary_vertices, true);
117  boundary_entities.merge(boundary_vertices);
118 
119  return boundary_entities;
120  };
121 
122  auto mark_boundary_dofs = [&](Range &&skin_edges) {
123  auto problem_manager = mField.getInterface<ProblemsManager>();
124  auto marker_ptr = boost::make_shared<std::vector<unsigned char>>();
125  problem_manager->markDofs(simpleInterface->getProblemName(), ROW,
126  skin_edges, *marker_ptr);
127  return marker_ptr;
128  };
129 
130  auto remove_dofs_from_problem = [&](Range &&skin_edges) {
132  auto problem_manager = mField.getInterface<ProblemsManager>();
133  CHKERR problem_manager->removeDofsOnEntities(
134  simpleInterface->getProblemName(), "P_IMAG", skin_edges, 0, 1);
136  };
137 
138  // Get global local vector of marked DOFs. Is global, since is set for all
139  // DOFs on processor. Is local since only DOFs on processor are in the
140  // vector. To access DOFs use local indices.
141  boundaryMarker = mark_boundary_dofs(get_ents_on_mesh_skin());
142  CHKERR remove_dofs_from_problem(get_ents_on_mesh_skin());
143 
145 }
146 //! [Applying essential BC]
147 
148 //! [Push operators to pipeline]
151  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
152 
153  double k = 90;
154  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-k", &k, PETSC_NULL);
155 
156  auto beta = [](const double, const double, const double) { return -1; };
157  auto k2 = [k](const double, const double, const double) { return pow(k, 2); };
158  auto kp = [k](const double, const double, const double) { return k; };
159  auto km = [k](const double, const double, const double) { return -k; };
160  auto integration_rule = [](int, int, int p_data) { return 2 * p_data; };
161 
162  auto set_domain = [&]() {
165  pipeline_mng->getOpDomainLhsPipeline(), {H1});
166 
167  pipeline_mng->getOpDomainLhsPipeline().push_back(
168  new OpSetBc("P_REAL", true, boundaryMarker));
169 
170  pipeline_mng->getOpDomainLhsPipeline().push_back(
171  new OpDomainGradGrad("P_REAL", "P_REAL", beta));
172  pipeline_mng->getOpDomainLhsPipeline().push_back(
173  new OpDomainGradGrad("P_IMAG", "P_IMAG", beta));
174 
175  pipeline_mng->getOpDomainLhsPipeline().push_back(
176  new OpDomainMass("P_REAL", "P_REAL", k2));
177  pipeline_mng->getOpDomainLhsPipeline().push_back(
178  new OpDomainMass("P_IMAG", "P_IMAG", k2));
179 
180  pipeline_mng->getOpDomainLhsPipeline().push_back(new OpUnSetBc("P_REAL"));
181 
184  };
185 
186  auto set_boundary = [&]() {
189  pipeline_mng->getOpBoundaryLhsPipeline(), {});
190 
191  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
192  new OpSetBc("P_REAL", true, boundaryMarker));
193  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
194  new OpBoundaryMass("P_IMAG", "P_REAL", kp));
195  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
196  new OpBoundaryMass("P_REAL", "P_IMAG", km));
197  pipeline_mng->getOpBoundaryLhsPipeline().push_back(new OpUnSetBc("P_REAL"));
198 
199  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
200  new OpSetBc("P_REAL", false, boundaryMarker));
201  pipeline_mng->getOpBoundaryLhsPipeline().push_back(
202  new OpBoundaryMass("P_REAL", "P_REAL", beta));
203  pipeline_mng->getOpBoundaryLhsPipeline().push_back(new OpUnSetBc("P_REAL"));
204 
205  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
206  new OpSetBc("P_REAL", false, boundaryMarker));
207  pipeline_mng->getOpBoundaryRhsPipeline().push_back(
208  new OpBoundarySource("P_REAL", beta));
209  pipeline_mng->getOpBoundaryRhsPipeline().push_back(new OpUnSetBc("P_REAL"));
210 
214  };
215 
216  CHKERR set_domain();
217  CHKERR set_boundary();
218 
220 }
221 //! [Push operators to pipeline]
222 
223 //! [Solve]
226  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
227  auto solver = pipeline_mng->createKSP();
228  CHKERR KSPSetFromOptions(solver);
229  CHKERR KSPSetUp(solver);
230 
231  auto dm = simpleInterface->getDM();
232  auto D = createDMVector(dm);
233  auto F = vectorDuplicate(D);
234 
235  CHKERR KSPSolve(solver, F, D);
236  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
237  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
238  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
240 }
241 //! [Solve]
242 
243 //! [Postprocess results]
246  PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
247  pipeline_mng->getDomainLhsFE().reset();
248  pipeline_mng->getDomainRhsFE().reset();
249  pipeline_mng->getBoundaryLhsFE().reset();
250  pipeline_mng->getBoundaryRhsFE().reset();
251 
253 
254  auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
255 
256  auto p_real_ptr = boost::make_shared<VectorDouble>();
257  auto p_imag_ptr = boost::make_shared<VectorDouble>();
258 
259  post_proc_fe->getOpPtrVector().push_back(
260  new OpCalculateScalarFieldValues("P_REAL", p_real_ptr));
261  post_proc_fe->getOpPtrVector().push_back(
262  new OpCalculateScalarFieldValues("P_IMAG", p_imag_ptr));
263 
265 
266  post_proc_fe->getOpPtrVector().push_back(
267 
268  new OpPPMap(
269 
270  post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
271 
272  {{"P_REAL", p_real_ptr}, {"P_IMAG", p_imag_ptr}},
273 
274  {}, {}, {}
275 
276  )
277 
278  );
279 
280  pipeline_mng->getDomainRhsFE() = post_proc_fe;
281  CHKERR pipeline_mng->loopFiniteElements();
282  CHKERR post_proc_fe->writeFile("out_helmholtz.h5m");
284 }
285 //! [Postprocess results]
286 
287 //! [Check results]
291 
292  auto dm = simpleInterface->getDM();
293  auto D = createDMVector(dm);
294  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
295  double nrm2;
296  CHKERR VecNorm(D, NORM_2, &nrm2);
297  MOFEM_LOG("WORLD", Sev::inform)
298  << std::setprecision(9) << "Solution norm " << nrm2;
299 
300  PetscBool test_flg = PETSC_FALSE;
301  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-test", &test_flg, PETSC_NULL);
302  if (test_flg) {
303  constexpr double regression_test = 97.261672;
304  constexpr double eps = 1e-6;
305  if (std::abs(nrm2 - regression_test) / regression_test > eps)
306  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
307  "Not converged solution");
308  }
310 }
311 //! [Check results]
312 
313 int main(int argc, char *argv[]) {
314 
315  // Initialisation of MoFEM/PETSc and MOAB data structures
316  const char param_file[] = "param_file.petsc";
317  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
318 
319  try {
320 
321  //! [Register MoFEM discrete manager in PETSc]
322  DMType dm_name = "DMMOFEM";
323  CHKERR DMRegister_MoFEM(dm_name);
324  //! [Register MoFEM discrete manager in PETSc
325 
326  //! [Create MoAB]
327  moab::Core mb_instance; ///< mesh database
328  moab::Interface &moab = mb_instance; ///< mesh database interface
329  //! [Create MoAB]
330 
331  //! [Create MoFEM]
332  MoFEM::Core core(moab); ///< finite element database
333  MoFEM::Interface &m_field = core; ///< finite element database insterface
334  //! [Create MoFEM]
335 
336  //! [Example]
337  Example ex(m_field);
338  CHKERR ex.runProblem();
339  //! [Example]
340  }
341  CATCH_ERRORS;
342 
344 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
Example::checkResults
MoFEMErrorCode checkResults()
[Postprocess results]
Definition: dynamic_first_order_con_law.cpp:1205
SPACE_DIM
constexpr int SPACE_DIM
Definition: helmholtz.cpp:22
main
int main(int argc, char *argv[])
[Check results]
Definition: helmholtz.cpp:313
MoFEM::PipelineManager::getDomainRhsFE
boost::shared_ptr< FEMethod > & getDomainRhsFE()
Definition: PipelineManager.hpp:405
Example::assembleSystem
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
Definition: dynamic_first_order_con_law.cpp:647
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
Example::boundaryMarker
boost::shared_ptr< std::vector< unsigned char > > boundaryMarker
Definition: helmholtz.cpp:42
H1
@ H1
continuous field
Definition: definitions.h:85
Example::Example
Example(MoFEM::Interface &m_field)
Definition: helmholtz.cpp:35
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::PipelineManager::loopFiniteElements
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
Definition: PipelineManager.cpp:19
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::PipelineManager::getBoundaryLhsFE
boost::shared_ptr< FEMethod > & getBoundaryLhsFE()
Definition: PipelineManager.hpp:409
Example::simpleInterface
Simple * simpleInterface
Definition: helmholtz.cpp:41
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
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
MoFEM::PipelineManager::getBoundaryRhsFE
boost::shared_ptr< FEMethod > & getBoundaryRhsFE()
Definition: PipelineManager.hpp:413
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Definition: child_and_parent.cpp:53
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
OpBoundaryMass
FormsIntegrators< BoundaryEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpBoundaryMass
[Only used with Hencky/nonlinear material]
Definition: seepage.cpp:71
ROW
@ ROW
Definition: definitions.h:136
MoFEM::PipelineManager::getOpBoundaryLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpBoundaryLhsPipeline()
Get the Op Boundary Lhs Pipeline object.
Definition: PipelineManager.hpp:797
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
Example::readMesh
MoFEMErrorCode readMesh()
[Run problem]
Definition: dynamic_first_order_con_law.cpp:463
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
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
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
Example
[Example]
Definition: plastic.cpp:177
double
OpBoundarySource
FormsIntegrators< EdgeEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, 1 > OpBoundarySource
Definition: helmholtz.cpp:31
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
Example::solveSystem
MoFEMErrorCode solveSystem()
[Solve]
Definition: dynamic_first_order_con_law.cpp:893
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
help
static char help[]
Definition: helmholtz.cpp:13
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::ProblemsManager::markDofs
MoFEMErrorCode markDofs(const std::string problem_name, RowColData rc, const enum MarkOP op, const Range ents, std::vector< unsigned char > &marker) const
Create vector with marked indices.
Definition: ProblemsManager.cpp:3548
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::PipelineManager::setBoundaryLhsIntegrationRule
MoFEMErrorCode setBoundaryLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:557
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
BiLinearForm
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
MoFEM::EdgeElementForcesAndSourcesCore::UserDataOperator
default operator for EDGE element
Definition: EdgeElementForcesAndSourcesCore.hpp:68
MoFEM::PipelineManager::setDomainLhsIntegrationRule
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:503
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
MoFEM::OpUnSetBc
Definition: FormsIntegrators.hpp:49
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_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
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: helmholtz.cpp:25
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:225
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Example::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Set up problem]
Definition: dynamic_first_order_con_law.cpp:536
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
MoFEM::OpSetBc
Set indices on entities on finite element.
Definition: FormsIntegrators.hpp:38
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::PipelineManager::getOpBoundaryRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpBoundaryRhsPipeline()
Get the Op Boundary Rhs Pipeline object.
Definition: PipelineManager.hpp:821
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:213
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::PipelineManager::setBoundaryRhsIntegrationRule
MoFEMErrorCode setBoundaryRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:584
Example::mField
MoFEM::Interface & mField
Definition: plastic.cpp:184
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
Example::outputResults
MoFEMErrorCode outputResults()
[Solve]
Definition: dynamic_first_order_con_law.cpp:1183
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
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpDomainMass
Definition: helmholtz.cpp:27
OpBoundaryMass
FormsIntegrators< EdgeEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpBoundaryMass
Definition: helmholtz.cpp:29
F
@ F
Definition: free_surface.cpp:394
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698