v0.15.5
Loading...
Searching...
No Matches
mesh_smoothing.cpp
Go to the documentation of this file.
1/** \file mesh_smoothing.cpp
2 * \brief Improve mesh quality using Volume-length quality measure with barrier
3 * \example mofem/tools/mesh_smoothing.cpp
4 *
5 */
6
7
8
9#include <MoFEM.hpp>
10using namespace MoFEM;
11
13#include <Smoother.hpp>
14#include <SurfaceSlidingConstrains.hpp>
15#include <VolumeLengthQuality.hpp>
16
17using namespace MoFEM;
18
19static char help[] = "mesh smoothing\n\n";
20
23PetscBool output_vtk = PETSC_TRUE;
24double tol = 0.1;
25double gamma_factor = 0.8;
26
27constexpr bool fix_edge_blocks = true;
28
29int main(int argc, char *argv[]) {
30
31 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
32
33 try {
34
35 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Mesh cut options", "none");
36 CHKERR PetscOptionsInt("-edges_block_set", "edges side set", "",
37 edges_block_set, &edges_block_set, PETSC_NULLPTR);
38 CHKERR PetscOptionsInt("-vertex_block_set", "vertex side set", "",
39 vertex_block_set, &vertex_block_set, PETSC_NULLPTR);
40 CHKERR PetscOptionsBool("-output_vtk", "if true outout vtk file", "",
41 output_vtk, &output_vtk, PETSC_NULLPTR);
42 CHKERR PetscOptionsScalar("-quality_reduction_tol", "",
43 "Tolerance of quality reduction", tol, &tol,
44 PETSC_NULLPTR);
45 CHKERR PetscOptionsScalar("-gamma_factor", "",
46 "Gamma factor", gamma_factor, &gamma_factor,
47 PETSC_NULLPTR);
48 PetscOptionsEnd();
49
50 // Create MoAB database
51 moab::Core moab_core; // create database
52 moab::Interface &moab = moab_core; // create interface to database
53 // Create MoFEM database and link it to MoAB
54 MoFEM::Core mofem_core(moab); // create database
55 MoFEM::Interface &m_field = mofem_core; // create interface to database
56 // Register DM Manager
57 CHKERR DMRegister_MoFEM("DMMOFEM"); // register MoFEM DM in PETSc
58
59 // Get simple interface is simplified version enabling quick and
60 // easy construction of problem.
61 Simple *simple_interface = m_field.getInterface<Simple>();
62
63 // Build problem with simple interface
64 {
65 // Get options for simple interface from command line
66 CHKERR simple_interface->getOptions();
67 // Load mesh file to database
68 CHKERR simple_interface->loadFile();
69
70 // Add domain field "U" in space H1 and Legendre base, Ainsworth recipe is
71 // used to construct base functions.
72 CHKERR simple_interface->addDomainField("MESH_NODE_POSITIONS", H1,
74 // Add Lagrange multiplier field on body boundary
75 CHKERR simple_interface->addBoundaryField("LAMBDA_SURFACE", H1,
77
78 // Set fields order domain and boundary fields.
79 CHKERR simple_interface->setFieldOrder("MESH_NODE_POSITIONS",
80 1); // to approximate function
81 CHKERR simple_interface->setFieldOrder("LAMBDA_SURFACE",
82 1); // to Lagrange multipliers
83
84 simple_interface->getDomainFEName() = "SMOOTHING";
85 simple_interface->getBoundaryFEName() = "SURFACE_SLIDING";
86
87 // Other fields and finite elements added to database directly
88 if (fix_edge_blocks) {
89 // do no add field
90 } else {
93 // Declare approximation fields
94 CHKERR m_field.add_field("LAMBDA_EDGE", H1, AINSWORTH_LOBATTO_BASE, 2,
95 MB_TAG_SPARSE, MF_ZERO);
96 Range edges;
98 ->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
99 true);
100 CHKERR m_field.add_ents_to_field_by_type(edges, MBEDGE,
101 "LAMBDA_EDGE");
103 ->synchroniseFieldEntities("LAMBDA_EDGE");
104 CHKERR m_field.set_field_order(0, MBVERTEX, "LAMBDA_EDGE", 1);
105
106 CHKERR m_field.add_finite_element("EDGE_SLIDING");
108 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
110 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
112 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
113 CHKERR m_field.modify_finite_element_add_field_row("EDGE_SLIDING",
114 "LAMBDA_EDGE");
115 CHKERR m_field.modify_finite_element_add_field_col("EDGE_SLIDING",
116 "LAMBDA_EDGE");
117 CHKERR m_field.modify_finite_element_add_field_data("EDGE_SLIDING",
118 "LAMBDA_EDGE");
119 CHKERR m_field.add_ents_to_finite_element_by_type(edges, MBEDGE,
120 "EDGE_SLIDING");
121 simple_interface->getOtherFiniteElements().push_back("EDGE_SLIDING");
122 }
123 }
124
125 CHKERR simple_interface->defineFiniteElements();
126 CHKERR simple_interface->defineProblem();
127 CHKERR simple_interface->buildFields();
128 CHKERR simple_interface->buildFiniteElements();
129 CHKERR simple_interface->buildProblem();
130
131 // Remove vertices form LAMBDA_SURFACE which are on the edges
134 Range edges;
135 CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
136 edges_block_set, BLOCKSET, 1, edges, true);
137 Range verts;
138 CHKERR m_field.get_moab().get_connectivity(edges, verts, true);
139 Simple *simple_interface = m_field.getInterface<Simple>();
140 auto prb_mng = m_field.getInterface<ProblemsManager>();
141 CHKERR prb_mng->removeDofsOnEntities(simple_interface->getProblemName(),
142 "LAMBDA_SURFACE", verts, 0, 1);
143 if (fix_edge_blocks) {
144 CHKERR prb_mng->removeDofsOnEntities(
145 simple_interface->getProblemName(), "MESH_NODE_POSITIONS", verts,
146 0, 3);
147 }
148 }
149 }
150
151 struct ElementsAndOperators {
152
153 MoFEM::Interface &mField;
154 Vec minQualityVec;
155 double *minQualityPtr;
156
157 ElementsAndOperators(MoFEM::Interface &m_field) : mField(m_field) {
158 ierr = VecCreateMPI(PETSC_COMM_WORLD, 1, m_field.get_comm_size(),
159 &minQualityVec);
160 CHKERRABORT(PETSC_COMM_WORLD, ierr);
161 ierr = VecGetArray(minQualityVec, &minQualityPtr);
162 CHKERRABORT(PETSC_COMM_WORLD, ierr);
163 }
164
165 virtual ~ElementsAndOperators() {
166 ierr = VecRestoreArray(minQualityVec, &minQualityPtr);
167 CHKERRABORT(PETSC_COMM_WORLD, ierr);
168 ierr = VecDestroy(&minQualityVec);
169 CHKERRABORT(PETSC_COMM_WORLD, ierr);
170 }
171
172 double getMinQuality() const { return *minQualityPtr; }
173
174 enum Tags {
175 SMOOTHING_TAG = 1,
176 SURFACE_CONSTRAINS_TAG,
177 EDGE_CONSTRAINS_TAG
178 };
179
180 boost::shared_ptr<Smoother> smootherFe;
181 boost::shared_ptr<Smoother::MyVolumeFE>
182 feSmootherRhs; ///< Integrate smoothing operators
183 boost::shared_ptr<Smoother::MyVolumeFE>
184 feSmootherLhs; ///< Integrate smoothing operators
185 boost::shared_ptr<VolumeLengthQuality<double> > volumeLengthDouble;
186 boost::shared_ptr<VolumeLengthQuality<adouble> > volumeLengthAdouble;
187
188 boost::shared_ptr<SurfaceSlidingConstrains> surfaceConstrain;
189 boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
190 skinOrientation;
191 boost::shared_ptr<EdgeSlidingConstrains> edgeConstrain;
192
193 boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore> minQualityFe;
194
196 double *minQualityPtr;
197 MinQualityOp(double *min_quality_ptr)
199 "MESH_NODE_POSITIONS", UserDataOperator::OPROW),
200 minQualityPtr(min_quality_ptr) {}
201 MoFEMErrorCode doWork(int side, EntityType type,
204 if (type != MBVERTEX)
206 double q = Tools::volumeLengthQuality(&*data.getFieldData().begin());
207 *minQualityPtr = fmin(*minQualityPtr, q);
209 }
210 };
211
212 MoFEMErrorCode createSmoothingFE() {
214 smootherFe = boost::shared_ptr<Smoother>(new Smoother(mField));
215 volumeLengthAdouble = boost::shared_ptr<VolumeLengthQuality<adouble> >(
217 volumeLengthDouble = boost::shared_ptr<VolumeLengthQuality<double> >(
219
220 Range tets;
221 CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
222 smootherFe->setOfBlocks[0].tEts.merge(tets);
223
224 smootherFe->setOfBlocks[0].materialDoublePtr = volumeLengthDouble;
225 smootherFe->setOfBlocks[0].materialAdoublePtr = volumeLengthAdouble;
226
227 // set element data
228 smootherFe->commonData.spatialPositions = "MESH_NODE_POSITIONS";
229 smootherFe->commonData.meshPositions = "NONE";
230
231 smootherFe->feRhs.meshPositionsFieldName = "NONE";
232 smootherFe->feLhs.meshPositionsFieldName = "NONE";
233 smootherFe->feRhs.addToRule = 0;
234 smootherFe->feLhs.addToRule = 0;
235
236 feSmootherRhs = smootherFe->feRhsPtr;
237 feSmootherLhs = smootherFe->feLhsPtr;
238
239 feSmootherRhs->getOpPtrVector().push_back(
241 "MESH_NODE_POSITIONS", smootherFe->commonData));
242 feSmootherRhs->getOpPtrVector().push_back(
244 "MESH_NODE_POSITIONS", smootherFe->setOfBlocks.at(0),
245 smootherFe->commonData, SMOOTHING_TAG, false));
246 feSmootherRhs->getOpPtrVector().push_back(new Smoother::OpRhsSmoother(
247 "MESH_NODE_POSITIONS", smootherFe->setOfBlocks[0],
248 smootherFe->commonData, smootherFe->smootherData));
249
250 feSmootherLhs->getOpPtrVector().push_back(
252 "MESH_NODE_POSITIONS", smootherFe->commonData));
253 feSmootherLhs->getOpPtrVector().push_back(
255 "MESH_NODE_POSITIONS", smootherFe->setOfBlocks.at(0),
256 smootherFe->commonData, SMOOTHING_TAG, true));
257 feSmootherLhs->getOpPtrVector().push_back(new Smoother::OpLhsSmoother(
258 "MESH_NODE_POSITIONS", "MESH_NODE_POSITIONS",
259 smootherFe->setOfBlocks.at(0), smootherFe->commonData,
260 smootherFe->smootherData, "LAMBDA_CRACKFRONT_AREA_TANGENT"));
261
262 minQualityFe =
263 boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore>(
265 minQualityFe->getOpPtrVector().push_back(
266 new MinQualityOp(minQualityPtr));
267
268 Range fixed_vertex;
271 CHKERR mField.getInterface<MeshsetsManager>()->getEntitiesByDimension(
272 vertex_block_set, BLOCKSET, 0, fixed_vertex, true);
273 }
274 auto prb_mng = mField.getInterface<ProblemsManager>();
275 Simple *simple_interface = mField.getInterface<Simple>();
276 CHKERR prb_mng->removeDofsOnEntities(simple_interface->getProblemName(),
277 "MESH_NODE_POSITIONS",
278 fixed_vertex, 0, 3);
279 CHKERR prb_mng->removeDofsOnEntities(simple_interface->getProblemName(),
280 "LAMBDA_SURFACE",
281 fixed_vertex, 0, 1);
282 if (!fix_edge_blocks)
285 CHKERR prb_mng->removeDofsOnEntities(
286 simple_interface->getProblemName(), "LAMBDA_EDGE", fixed_vertex,
287 0, 1);
288 }
289
291 }
292
293 MoFEMErrorCode createConstrians() {
295 skinOrientation = boost::shared_ptr<
298 surfaceConstrain = boost::shared_ptr<SurfaceSlidingConstrains>(
299 skinOrientation,
300 new SurfaceSlidingConstrains(mField, *skinOrientation));
301 surfaceConstrain->setOperators(SURFACE_CONSTRAINS_TAG, "LAMBDA_SURFACE",
302 "MESH_NODE_POSITIONS");
303
304 if (!fix_edge_blocks)
307 Range edges;
309 ->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
310 true);
311
312 Range tets;
313 CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
314 Skinner skin(&mField.get_moab());
315 Range skin_faces; // skin faces from 3d ents
316 CHKERR skin.find_skin(0, tets, false, skin_faces);
317
318 edgeConstrain = boost::shared_ptr<EdgeSlidingConstrains>(
319 new EdgeSlidingConstrains(mField));
320 CHKERR edgeConstrain->setOperators(EDGE_CONSTRAINS_TAG, edges,
321 skin_faces, "LAMBDA_EDGE",
322 "MESH_NODE_POSITIONS");
323
324 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
325 // mField.get_moab(), "out_edges.vtk", edges);
326 }
328 }
329
330 MoFEMErrorCode addFEtoDM(DM dm) {
332 boost::shared_ptr<ForcesAndSourcesCore> null;
333
334 CHKERR DMMoFEMSNESSetFunction(dm, "SMOOTHING", feSmootherRhs, null,
335 null);
336 CHKERR DMMoFEMSNESSetFunction(dm, "SURFACE_SLIDING",
337 surfaceConstrain->feRhsPtr, null, null);
338 if (!fix_edge_blocks)
339 CHKERR DMMoFEMSNESSetFunction(dm, "EDGE_SLIDING",
340 edgeConstrain->feRhsPtr, null, null);
341
342 CHKERR DMMoFEMSNESSetJacobian(dm, "SMOOTHING", feSmootherLhs, null,
343 null);
344 CHKERR DMMoFEMSNESSetJacobian(dm, "SURFACE_SLIDING",
345 surfaceConstrain->feLhsPtr, null, null);
346 if (!fix_edge_blocks)
347 CHKERR DMMoFEMSNESSetJacobian(dm, "EDGE_SLIDING",
348 edgeConstrain->feLhsPtr, null, null);
349
350 // MoFEM::SnesCtx *snes_ctx;
351 // DMMoFEMGetSnesCtx(dm,&snes_ctx);
352 // snes_ctx->vErify = true;
353
355 }
356
357 MoFEMErrorCode calcuteMinQuality(DM dm) {
359 *minQualityPtr = 1;
360 CHKERR DMoFEMLoopFiniteElements(dm, "SMOOTHING", minQualityFe.get());
361 CHKERR VecMin(minQualityVec, PETSC_NULLPTR, minQualityPtr);
363 }
364 };
365
366 ElementsAndOperators elements_and_operators(m_field);
367 CHKERR elements_and_operators.createSmoothingFE();
368 CHKERR elements_and_operators.createConstrians();
369
370 DM dm;
371 CHKERR simple_interface->getDM(&dm);
372 CHKERR elements_and_operators.addFEtoDM(dm);
373
374 struct Solve {
375
376 MoFEMErrorCode operator()(DM dm) const {
378
379 // Create the right hand side vector and vector of unknowns
380 Vec F, D;
381 CHKERR DMCreateGlobalVector(dm, &F);
382 // Create unknown vector by creating duplicate copy of F vector. only
383 // structure is duplicated no values.
384 CHKERR VecDuplicate(F, &D);
385
386 CHKERR zeroLambdaFields(dm);
387 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
388
389 // Create solver and link it to DM
390 SNES solver;
391 CHKERR SNESCreate(PETSC_COMM_WORLD, &solver);
392 CHKERR SNESSetFromOptions(solver);
393 CHKERR SNESSetDM(solver, dm);
394 // Set-up solver, is type of solver and pre-conditioners
395 CHKERR SNESSetUp(solver);
396
397 auto hook = [&](SNES snes, Vec x, Vec F,
398 boost::shared_ptr<CacheTuple> cache_ptr,
399 void *ctx) -> MoFEMErrorCode {
401
404
405 auto *m_field_ptr = getInterfacePtr(dm);
406 auto post_proc_fe =
407 boost::make_shared<PostProcEleDomain>(*m_field_ptr);
408 auto &pip = post_proc_fe->getOpPtrVector();
409 auto X_ptr = boost::make_shared<MatrixDouble>();
410 pip.push_back(new OpCalculateVectorFieldValues<3>(
411 "MESH_NODE_POSITIONS", X_ptr));
412 auto res_X_ptr = boost::make_shared<MatrixDouble>();
413 pip.push_back(new OpCalculateVectorFieldValues<3>(
414 "MESH_NODE_POSITIONS", res_X_ptr, SmartPetscObj<Vec>(F, true)));
415 auto lambda_ptr = boost::make_shared<VectorDouble>();
416 pip.push_back(
417 new OpCalculateScalarFieldValues("LAMBDA_SURFACE", lambda_ptr));
418 auto res_lambda_ptr = boost::make_shared<VectorDouble>();
419 pip.push_back(new OpCalculateScalarFieldValues(
420 "LAMBDA_SURFACE", res_lambda_ptr, SmartPetscObj<Vec>(F, true)));
421
423 pip.push_back(
424
425 new OpPPMap(
426
427 post_proc_fe->getPostProcMesh(),
428 post_proc_fe->getMapGaussPts(),
429
430 {{"LAMBDA_SURFACE", lambda_ptr},
431 {"RES_LAMBDA_SURFACE", res_lambda_ptr}},
432
433 {{"MESH_NODE_POSITIONS", X_ptr},
434
435 {"RES_MESH_NODE_POSITIONS", res_X_ptr}},
436
437 {},
438
439 {}
440
441 )
442
443 );
444
445 CHKERR DMoFEMLoopFiniteElements(dm, "SURFACE_SLIDING", post_proc_fe);
446 CHKERR post_proc_fe->writeFile("debug_smoothing.h5m");
448 };
449
450 auto snes_ctx_ptr = getDMSnesCtx(dm);
451 snes_ctx_ptr->getRhsHook() = hook;
452
453 // At solution process, KSP solver using DM creates matrices, Calculate
454 // values of the left hand side and the right hand side vector. then
455 // solves system of equations. Results are stored in vector D.
456 CHKERR SNESSolve(solver, F, D);
457
458 // Scatter solution on the mesh. Stores unknown vector on field on the
459 // mesh.
460 CHKERR DMoFEMMeshToGlobalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
461 // Clean data. Solver and vector are not needed any more.
462 CHKERR SNESDestroy(&solver);
463 CHKERR VecDestroy(&D);
464 CHKERR VecDestroy(&F);
465
467 }
468
469 MoFEMErrorCode setCoordsFromField(DM dm) const {
471 MoFEM::Interface *m_field_ptr;
472 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
473 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
474 "MESH_NODE_POSITIONS", it)) {
475 if (it->get()->getEntType() != MBVERTEX)
476 continue;
477 VectorDouble3 coords(3);
478 for(int dd = 0;dd!=3;++dd)
479 coords[dd] = it->get()->getEntFieldData()[dd];
480 EntityHandle ent = it->get()->getEnt();
481 CHKERR m_field_ptr->get_moab().set_coords(&ent, 1, &*coords.begin());
482 }
484 }
485
486 MoFEMErrorCode setFieldFromCoords(DM dm) const {
488 MoFEM::Interface *m_field_ptr;
489 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
490 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
491 "MESH_NODE_POSITIONS", it)) {
492 if (it->get()->getEntType() != MBVERTEX)
493 continue;
494 EntityHandle ent = it->get()->getEnt();
495 VectorDouble3 coords(3);
496 CHKERR m_field_ptr->get_moab().get_coords(&ent, 1, &*coords.begin());
497 for(int dd = 0;dd!=3;++dd)
498 it->get()->getEntFieldData()[dd] = coords[dd];
499 }
501 }
502
503 private:
504 MoFEMErrorCode zeroLambdaFields(DM dm) const {
506 MoFEM::Interface *m_field_ptr;
507 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
508 CHKERR m_field_ptr->getInterface<FieldBlas>()->setField(
509 0, MBVERTEX, "LAMBDA_SURFACE");
511 }
512
513 };
514
515 Solve solve;
516 CHKERR solve.setFieldFromCoords(dm);
517
518 CHKERR elements_and_operators.calcuteMinQuality(dm);
519 double min_quality = elements_and_operators.getMinQuality();
520 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f\n", min_quality);
521
522 double gamma = min_quality > 0 ? gamma_factor * min_quality
523 : min_quality / gamma_factor;
524 elements_and_operators.volumeLengthDouble->gAmma = gamma;
525 elements_and_operators.volumeLengthAdouble->gAmma = gamma;
526
527 double min_quality_p, eps;
528 do {
529
530 min_quality_p = min_quality;
531
532 CHKERR solve(dm);
533
534 CHKERR solve.setCoordsFromField(dm);
535 CHKERR elements_and_operators.calcuteMinQuality(dm);
536 min_quality = elements_and_operators.getMinQuality();
537
538 eps = (min_quality - min_quality_p) / min_quality;
539 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f eps = %4.3f\n",
540 min_quality, eps);
541
542 double gamma = min_quality > 0 ? gamma_factor * min_quality
543 : min_quality / gamma_factor;
544 elements_and_operators.volumeLengthDouble->gAmma = gamma;
545 elements_and_operators.volumeLengthAdouble->gAmma = gamma;
546
547 } while (eps > tol);
548
549 // if (m_field.getInterface<MeshsetsManager>()->checkMeshset(edges_block_set,
550 // BLOCKSET)) {
551 // Range edges;
552 // CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
553 // edges_block_set, BLOCKSET, 1, edges, true);
554
555 // Range tets;
556 // CHKERR moab.get_entities_by_type(0,MBTET,tets);
557 // Skinner skin(&moab);
558 // Range skin_faces; // skin faces from 3d ents
559 // CHKERR skin.find_skin(0, tets, false, skin_faces);
560
561 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::setTags(moab, edges,
562 // skin_faces);
563 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
564 // moab, "out_edges.vtk", edges);
565 // }
566
567 if (output_vtk)
568 CHKERR m_field.getInterface<BitRefManager>()->writeBitLevelByType(
569 BitRefLevel().set(0), BitRefLevel().set(), MBTET, "out.vtk", "VTK",
570 "");
571 }
573
575}
576
int main()
static const double eps
#define CATCH_ERRORS
Catch errors.
@ MF_ZERO
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ AINSWORTH_LOBATTO_BASE
Definition definitions.h:62
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ BLOCKSET
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PostProcEleByDim< SPACE_DIM >::PostProcEleDomain PostProcEleDomain
@ F
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition DMMoFEM.cpp:708
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:514
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition DMMoFEM.cpp:749
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
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:576
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition DMMoFEM.cpp:410
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition DMMoFEM.cpp:525
#define _IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(MFIELD, NAME, IT)
loop over all dofs from a moFEM field and particular field
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
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.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
bool checkMeshset(const int ms_id, const CubitBCType cubit_bc_type) const
Check for CUBIT meshset by ID and type.
MoFEM::PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle Ele
double D
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorBoundedArray< double, 3 > VectorDouble3
Definition Types.hpp:92
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
static PetscErrorCode solve(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1322
auto getInterfacePtr(DM dm)
Get the Interface Ptr object.
Definition DMMoFEM.hpp:1168
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1262
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Managing BitRefLevels.
Managing BitRefLevels.
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode add_field(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_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorDouble & getFieldData() const
Get DOF values on entity.
Basic algebra on fields.
Definition FieldBlas.hpp:21
@ OPROW
operator doWork function is executed on FE rows
Interface for managing meshsets containing materials and boundary conditions.
Specialization for double precision scalar field values calculation.
Specialization for double precision vector field values calculation.
Post post-proc data at points from hash maps.
Template struct for dimension-specific finite element types.
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode buildProblem()
Build problem.
Definition Simple.cpp:721
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:261
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
Definition Simple.hpp:436
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition Simple.cpp:471
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:191
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition Simple.cpp:659
MoFEMErrorCode addBoundaryField(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 boundary.
Definition Simple.cpp:355
std::vector< std::string > & getOtherFiniteElements()
Get the Other Finite Elements.
Definition Simple.hpp:508
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:799
MoFEMErrorCode buildFields()
Build fields.
Definition Simple.cpp:584
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:575
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition Simple.cpp:551
const std::string getProblemName() const
Get the Problem Name.
Definition Simple.hpp:450
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition Simple.hpp:429
intrusive_ptr for managing petsc objects
static double volumeLengthQuality(const double *coords)
Calculate tetrahedron volume length quality.
Definition Tools.cpp:15
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Class implemented by user to detect face orientation.
Shape preserving constrains, i.e. nodes sliding on body surface.
Volume Length Quality.
constexpr bool fix_edge_blocks
int vertex_block_set
int edges_block_set
@ BARRIER_AND_QUALITY
double tol
PetscBool output_vtk
static char help[]
constexpr bool fix_edge_blocks
double gamma_factor
int vertex_block_set
int edges_block_set