v0.14.0
Loading...
Searching...
No Matches
Functions | Variables
mesh_smoothing.cpp File Reference

Improve mesh quality using Volume-length quality measure with barrier. More...

#include <MoFEM.hpp>
#include <BasicFiniteElements.hpp>
#include <Smoother.hpp>
#include <SurfaceSlidingConstrains.hpp>
#include <VolumeLengthQuality.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "mesh smoothing\n\n"
 
PetscBool flg_myfile = PETSC_TRUE
 
char mesh_file_name [255]
 
int edges_block_set = 1
 
int vertex_block_set = 2
 
PetscBool output_vtk = PETSC_TRUE
 
double tol = 0.1
 
double gamma_factor = 0.8
 

Detailed Description

Improve mesh quality using Volume-length quality measure with barrier.

Definition in file mesh_smoothing.cpp.

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

< Integrate smoothing operators

< Integrate smoothing operators

Definition at line 30 of file mesh_smoothing.cpp.

30 {
31
32 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
33
34 try {
35
36 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Mesh cut options", "none");
37 CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
39 CHKERR PetscOptionsInt("-edges_block_set", "edges side set", "",
40 edges_block_set, &edges_block_set, PETSC_NULL);
41 CHKERR PetscOptionsInt("-vertex_block_set", "vertex side set", "",
43 CHKERR PetscOptionsBool("-output_vtk", "if true outout vtk file", "",
44 output_vtk, &output_vtk, PETSC_NULL);
45 CHKERR PetscOptionsScalar("-quality_reduction_tol", "",
46 "Tolerance of quality reduction", tol, &tol,
47 PETSC_NULL);
48 CHKERR PetscOptionsScalar("-gamma_factor", "",
49 "Gamma factor", gamma_factor, &gamma_factor,
50 PETSC_NULL);
51 ierr = PetscOptionsEnd();
53
54 // Create MoAB database
55 moab::Core moab_core; // create database
56 moab::Interface &moab = moab_core; // create interface to database
57 // Create MoFEM database and link it to MoAB
58 MoFEM::Core mofem_core(moab); // create database
59 MoFEM::Interface &m_field = mofem_core; // create interface to database
60 // Register DM Manager
61 CHKERR DMRegister_MoFEM("DMMOFEM"); // register MoFEM DM in PETSc
62
63 // Get simple interface is simplified version enabling quick and
64 // easy construction of problem.
65 Simple *simple_interface;
66 // Query interface and get pointer to Simple interface
67 CHKERR m_field.getInterface(simple_interface);
68
69 // Build problem with simple interface
70 {
71 // Get options for simple interface from command line
72 CHKERR simple_interface->getOptions();
73 // Load mesh file to database
74 CHKERR simple_interface->loadFile();
75
76 // Add domain filed "U" in space H1 and Legendre base, Ainsworth recipe is
77 // used to construct base functions.
78 CHKERR simple_interface->addDomainField("MESH_NODE_POSITIONS", H1,
80 // Add Lagrange multiplier field on body boundary
81 CHKERR simple_interface->addBoundaryField("LAMBDA_SURFACE", H1,
83
84 // Set fields order domain and boundary fields.
85 CHKERR simple_interface->setFieldOrder("MESH_NODE_POSITIONS",
86 1); // to approximate function
87 CHKERR simple_interface->setFieldOrder("LAMBDA_SURFACE",
88 1); // to Lagrange multipliers
89
90 simple_interface->getDomainFEName() = "SMOOTHING";
91 simple_interface->getBoundaryFEName() = "SURFACE_SLIDING";
92
93 // Other fields and finite elements added to database directly
94 {
97 // Declare approximation fields
98 CHKERR m_field.add_field("LAMBDA_EDGE", H1, AINSWORTH_LOBATTO_BASE, 2,
99 MB_TAG_SPARSE, MF_ZERO);
100 Range edges;
102 ->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
103 true);
104 CHKERR m_field.add_ents_to_field_by_type(edges, MBEDGE,
105 "LAMBDA_EDGE");
107 ->synchroniseFieldEntities("LAMBDA_EDGE");
108 CHKERR m_field.set_field_order(0, MBVERTEX, "LAMBDA_EDGE", 1);
109
110 CHKERR m_field.add_finite_element("EDGE_SLIDING");
112 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
114 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
116 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
117 CHKERR m_field.modify_finite_element_add_field_row("EDGE_SLIDING",
118 "LAMBDA_EDGE");
119 CHKERR m_field.modify_finite_element_add_field_col("EDGE_SLIDING",
120 "LAMBDA_EDGE");
121 CHKERR m_field.modify_finite_element_add_field_data("EDGE_SLIDING",
122 "LAMBDA_EDGE");
123 CHKERR m_field.add_ents_to_finite_element_by_type(edges, MBEDGE,
124 "EDGE_SLIDING");
125 simple_interface->getOtherFiniteElements().push_back("EDGE_SLIDING");
126 }
127 }
128
129 CHKERR simple_interface->defineFiniteElements();
130 CHKERR simple_interface->defineProblem();
131 CHKERR simple_interface->buildFields();
132 // Remove vertices form LAMBDA_SURFACE which are on the edges
135 Range edges;
136 CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
137 edges_block_set, BLOCKSET, 1, edges, true);
138 Range verts;
139 CHKERR m_field.get_moab().get_connectivity(edges, verts, true);
140 CHKERR m_field.remove_ents_from_field("LAMBDA_SURFACE", verts);
141 }
142 CHKERR simple_interface->buildFiniteElements();
143 CHKERR simple_interface->buildProblem();
144 }
145
146 struct ElementsAndOperators {
147
148 MoFEM::Interface &mField;
149 Vec minQualityVec;
150 double *minQualityPtr;
151
152 ElementsAndOperators(MoFEM::Interface &m_field) : mField(m_field) {
153 ierr = VecCreateMPI(PETSC_COMM_WORLD, 1, m_field.get_comm_size(),
154 &minQualityVec);
155 CHKERRABORT(PETSC_COMM_WORLD, ierr);
156 ierr = VecGetArray(minQualityVec, &minQualityPtr);
157 CHKERRABORT(PETSC_COMM_WORLD, ierr);
158 }
159
160 virtual ~ElementsAndOperators() {
161 ierr = VecRestoreArray(minQualityVec, &minQualityPtr);
162 CHKERRABORT(PETSC_COMM_WORLD, ierr);
163 ierr = VecDestroy(&minQualityVec);
164 CHKERRABORT(PETSC_COMM_WORLD, ierr);
165 }
166
167 double getMinQuality() const { return *minQualityPtr; }
168
169 enum Tags {
170 SMOOTHING_TAG = 1,
171 SURFACE_CONSTRAINS_TAG,
172 EDGE_CONSTRAINS_TAG
173 };
174
175 boost::shared_ptr<Smoother> smootherFe;
176 boost::shared_ptr<Smoother::MyVolumeFE>
177 feSmootherRhs; ///< Integrate smoothing operators
178 boost::shared_ptr<Smoother::MyVolumeFE>
179 feSmootherLhs; ///< Integrate smoothing operators
180 boost::shared_ptr<VolumeLengthQuality<double> > volumeLengthDouble;
181 boost::shared_ptr<VolumeLengthQuality<adouble> > volumeLengthAdouble;
182
183 boost::shared_ptr<SurfaceSlidingConstrains> surfaceConstrain;
184 boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
185 skinOrientation;
186 boost::shared_ptr<EdgeSlidingConstrains> edgeConstrain;
187
188 boost::shared_ptr<DirichletFixFieldAtEntitiesBc> fixMaterialEnts;
189
190 boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore> minQualityFe;
191
193 double *minQualityPtr;
194 MinQualityOp(double *min_quality_ptr)
196 "MESH_NODE_POSITIONS", UserDataOperator::OPROW),
197 minQualityPtr(min_quality_ptr) {}
198 MoFEMErrorCode doWork(int side, EntityType type,
201 if (type != MBVERTEX)
203 double q = Tools::volumeLengthQuality(&*data.getFieldData().begin());
204 *minQualityPtr = fmin(*minQualityPtr, q);
206 }
207 };
208
209 MoFEMErrorCode createSmoothingFE() {
211 smootherFe = boost::shared_ptr<Smoother>(new Smoother(mField));
212 volumeLengthAdouble = boost::shared_ptr<VolumeLengthQuality<adouble> >(
214 volumeLengthDouble = boost::shared_ptr<VolumeLengthQuality<double> >(
216
217 Range tets;
218 CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
219 smootherFe->setOfBlocks[0].tEts.merge(tets);
220
221 smootherFe->setOfBlocks[0].materialDoublePtr = volumeLengthDouble;
222 smootherFe->setOfBlocks[0].materialAdoublePtr = volumeLengthAdouble;
223
224 // set element data
225 smootherFe->commonData.spatialPositions = "MESH_NODE_POSITIONS";
226 smootherFe->commonData.meshPositions = "NONE";
227
228 smootherFe->feRhs.meshPositionsFieldName = "NONE";
229 smootherFe->feLhs.meshPositionsFieldName = "NONE";
230 smootherFe->feRhs.addToRule = 0;
231 smootherFe->feLhs.addToRule = 0;
232
233 feSmootherRhs = smootherFe->feRhsPtr;
234 feSmootherLhs = smootherFe->feLhsPtr;
235
236 // Smoother right hand side
237 feSmootherRhs->getOpPtrVector().push_back(
239 "MESH_NODE_POSITIONS", smootherFe->commonData));
240 feSmootherRhs->getOpPtrVector().push_back(
242 "MESH_NODE_POSITIONS", smootherFe->setOfBlocks.at(0),
243 smootherFe->commonData, SMOOTHING_TAG, false));
244 feSmootherRhs->getOpPtrVector().push_back(new Smoother::OpRhsSmoother(
245 "MESH_NODE_POSITIONS", smootherFe->setOfBlocks[0],
246 smootherFe->commonData, smootherFe->smootherData));
247
248 // Smoother left hand side
249 feSmootherLhs->getOpPtrVector().push_back(
251 "MESH_NODE_POSITIONS", smootherFe->commonData));
252 feSmootherLhs->getOpPtrVector().push_back(
254 "MESH_NODE_POSITIONS", smootherFe->setOfBlocks.at(0),
255 smootherFe->commonData, SMOOTHING_TAG, true));
256 feSmootherLhs->getOpPtrVector().push_back(new Smoother::OpLhsSmoother(
257 "MESH_NODE_POSITIONS", "MESH_NODE_POSITIONS",
258 smootherFe->setOfBlocks.at(0), smootherFe->commonData,
259 smootherFe->smootherData, "LAMBDA_CRACKFRONT_AREA_TANGENT"));
260
261 minQualityFe =
262 boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore>(
264 minQualityFe->getOpPtrVector().push_back(
265 new MinQualityOp(minQualityPtr));
266
267 Range fixed_vertex;
270 CHKERR mField.getInterface<MeshsetsManager>()->getEntitiesByDimension(
271 vertex_block_set, BLOCKSET, 0, fixed_vertex, true);
272 }
273 fixMaterialEnts = boost::shared_ptr<DirichletFixFieldAtEntitiesBc>(
274 new DirichletFixFieldAtEntitiesBc(mField, "MESH_NODE_POSITIONS",
275 fixed_vertex));
276 fixMaterialEnts->fieldNames.push_back("LAMBDA_SURFACE");
277 fixMaterialEnts->fieldNames.push_back("LAMBDA_EDGE");
278
280 }
281
282 MoFEMErrorCode createConstrians() {
284 skinOrientation = boost::shared_ptr<
287 surfaceConstrain = boost::shared_ptr<SurfaceSlidingConstrains>(
288 skinOrientation,
289 new SurfaceSlidingConstrains(mField, *skinOrientation));
290 surfaceConstrain->setOperators(SURFACE_CONSTRAINS_TAG, "LAMBDA_SURFACE",
291 "MESH_NODE_POSITIONS");
292
295 Range edges;
297 ->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
298 true);
299
300 Range tets;
301 CHKERR mField.get_moab().get_entities_by_type(0, MBTET, tets);
302 Skinner skin(&mField.get_moab());
303 Range skin_faces; // skin faces from 3d ents
304 CHKERR skin.find_skin(0, tets, false, skin_faces);
305
306 edgeConstrain = boost::shared_ptr<EdgeSlidingConstrains>(
307 new EdgeSlidingConstrains(mField));
308 CHKERR edgeConstrain->setOperators(EDGE_CONSTRAINS_TAG, edges,
309 skin_faces, "LAMBDA_EDGE",
310 "MESH_NODE_POSITIONS");
311
312 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
313 // mField.get_moab(), "out_edges.vtk", edges);
314
315 }
317 }
318
319 MoFEMErrorCode addFEtoDM(DM dm) {
321 boost::shared_ptr<ForcesAndSourcesCore> null;
322
323 CHKERR DMMoFEMSNESSetFunction(dm, DM_NO_ELEMENT, null, fixMaterialEnts,
324 null);
325 CHKERR DMMoFEMSNESSetFunction(dm, "SMOOTHING", feSmootherRhs, null,
326 null);
327 CHKERR DMMoFEMSNESSetFunction(dm, "SURFACE_SLIDING",
328 surfaceConstrain->feRhsPtr, null, null);
329 CHKERR DMMoFEMSNESSetFunction(dm, "EDGE_SLIDING",
330 edgeConstrain->feRhsPtr, null, null);
332 fixMaterialEnts);
333
334 CHKERR DMMoFEMSNESSetJacobian(dm, DM_NO_ELEMENT, null, fixMaterialEnts,
335 null);
336 CHKERR DMMoFEMSNESSetJacobian(dm, "SMOOTHING", feSmootherLhs, null,
337 null);
338 CHKERR DMMoFEMSNESSetJacobian(dm, "SURFACE_SLIDING",
339 surfaceConstrain->feLhsPtr, null, null);
340 CHKERR DMMoFEMSNESSetJacobian(dm, "EDGE_SLIDING",
341 edgeConstrain->feLhsPtr, null, null);
343 fixMaterialEnts);
344
345 // MoFEM::SnesCtx *snes_ctx;
346 // DMMoFEMGetSnesCtx(dm,&snes_ctx);
347 // snes_ctx->vErify = true;
348
350 }
351
352 MoFEMErrorCode calcuteMinQuality(DM dm) {
354 *minQualityPtr = 1;
355 CHKERR DMoFEMLoopFiniteElements(dm, "SMOOTHING", minQualityFe.get());
356 CHKERR VecMin(minQualityVec, PETSC_NULL, minQualityPtr);
358 }
359 };
360
361 ElementsAndOperators elements_and_operators(m_field);
362 CHKERR elements_and_operators.createSmoothingFE();
363 CHKERR elements_and_operators.createConstrians();
364
365 DM dm;
366 CHKERR simple_interface->getDM(&dm);
367 CHKERR elements_and_operators.addFEtoDM(dm);
368
369 struct Solve {
370
371 MoFEMErrorCode operator()(DM dm) const {
373
374 // Create the right hand side vector and vector of unknowns
375 Vec F, D;
376 CHKERR DMCreateGlobalVector(dm, &F);
377 // Create unknown vector by creating duplicate copy of F vector. only
378 // structure is duplicated no values.
379 CHKERR VecDuplicate(F, &D);
380
381 CHKERR zeroLambdaFields(dm);
382 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_FORWARD);
383
384 // Create solver and link it to DM
385 SNES solver;
386 CHKERR SNESCreate(PETSC_COMM_WORLD, &solver);
387 CHKERR SNESSetFromOptions(solver);
388 CHKERR SNESSetDM(solver, dm);
389 // Set-up solver, is type of solver and pre-conditioners
390 CHKERR SNESSetUp(solver);
391 // At solution process, KSP solver using DM creates matrices, Calculate
392 // values of the left hand side and the right hand side vector. then
393 // solves system of equations. Results are stored in vector D.
394 CHKERR SNESSolve(solver, F, D);
395
396 // Scatter solution on the mesh. Stores unknown vector on field on the
397 // mesh.
398 CHKERR DMoFEMMeshToGlobalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
399 // Clean data. Solver and vector are not needed any more.
400 CHKERR SNESDestroy(&solver);
401 CHKERR VecDestroy(&D);
402 CHKERR VecDestroy(&F);
403
405 }
406
407 MoFEMErrorCode setCoordsFromField(DM dm) const {
409 MoFEM::Interface *m_field_ptr;
410 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
411 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
412 "MESH_NODE_POSITIONS", it)) {
413 if (it->get()->getEntType() != MBVERTEX)
414 continue;
415 VectorDouble3 coords(3);
416 for(int dd = 0;dd!=3;++dd)
417 coords[dd] = it->get()->getEntFieldData()[dd];
418 EntityHandle ent = it->get()->getEnt();
419 CHKERR m_field_ptr->get_moab().set_coords(&ent, 1, &*coords.begin());
420 }
422 }
423
424 MoFEMErrorCode setFieldFromCoords(DM dm) const {
426 MoFEM::Interface *m_field_ptr;
427 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
428 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
429 "MESH_NODE_POSITIONS", it)) {
430 if (it->get()->getEntType() != MBVERTEX)
431 continue;
432 EntityHandle ent = it->get()->getEnt();
433 VectorDouble3 coords(3);
434 CHKERR m_field_ptr->get_moab().get_coords(&ent, 1, &*coords.begin());
435 for(int dd = 0;dd!=3;++dd)
436 it->get()->getEntFieldData()[dd] = coords[dd];
437 }
439 }
440
441 private:
442 MoFEMErrorCode zeroLambdaFields(DM dm) const {
444 MoFEM::Interface *m_field_ptr;
445 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
446 CHKERR m_field_ptr->getInterface<FieldBlas>()->setField(
447 0, MBVERTEX, "LAMBDA_SURFACE");
449 }
450
451 };
452
453 Solve solve;
454 CHKERR solve.setFieldFromCoords(dm);
455
456 CHKERR elements_and_operators.calcuteMinQuality(dm);
457 double min_quality = elements_and_operators.getMinQuality();
458 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f\n", min_quality);
459
460 double gamma = min_quality > 0 ? gamma_factor * min_quality
461 : min_quality / gamma_factor;
462 elements_and_operators.volumeLengthDouble->gAmma = gamma;
463 elements_and_operators.volumeLengthAdouble->gAmma = gamma;
464
465 double min_quality_p, eps;
466 do {
467
468 min_quality_p = min_quality;
469
470 CHKERR solve(dm);
471
472 CHKERR solve.setCoordsFromField(dm);
473 CHKERR elements_and_operators.calcuteMinQuality(dm);
474 min_quality = elements_and_operators.getMinQuality();
475
476 eps = (min_quality - min_quality_p) / min_quality;
477 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f eps = %4.3f\n",
478 min_quality, eps);
479
480 double gamma = min_quality > 0 ? gamma_factor * min_quality
481 : min_quality / gamma_factor;
482 elements_and_operators.volumeLengthDouble->gAmma = gamma;
483 elements_and_operators.volumeLengthAdouble->gAmma = gamma;
484
485 } while (eps > tol);
486
487 // if (m_field.getInterface<MeshsetsManager>()->checkMeshset(edges_block_set,
488 // BLOCKSET)) {
489 // Range edges;
490 // CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
491 // edges_block_set, BLOCKSET, 1, edges, true);
492
493 // Range tets;
494 // CHKERR moab.get_entities_by_type(0,MBTET,tets);
495 // Skinner skin(&moab);
496 // Range skin_faces; // skin faces from 3d ents
497 // CHKERR skin.find_skin(0, tets, false, skin_faces);
498
499 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::setTags(moab, edges,
500 // skin_faces);
501 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
502 // moab, "out_edges.vtk", edges);
503 // }
504
505 if (output_vtk)
506 CHKERR m_field.getInterface<BitRefManager>()->writeBitLevelByType(
507 BitRefLevel().set(0), BitRefLevel().set(), MBTET, "out.vtk", "VTK",
508 "");
509 }
511
513}
#define DM_NO_ELEMENT
Definition DMMoFEM.hpp:10
ForcesAndSourcesCore::UserDataOperator UserDataOperator
@ BARRIER_AND_QUALITY
static const double eps
#define CATCH_ERRORS
Catch errors.
@ MF_ZERO
Definition definitions.h:98
@ 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 ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ 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 ...
@ 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:704
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:509
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:745
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:47
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:572
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition DMMoFEM.cpp:400
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition DMMoFEM.cpp:521
#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_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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 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.
virtual MoFEMErrorCode remove_ents_from_field(const std::string name, const EntityHandle meshset, const EntityType type, int verb=DEFAULT_VERBOSITY)=0
remove entities from field
bool checkMeshset(const int ms_id, const CubitBCType cubit_bc_type) const
check for CUBIT Id and CUBIT type
double D
double tol
PetscBool output_vtk
static char help[]
PetscBool flg_myfile
double gamma_factor
int vertex_block_set
int edges_block_set
char mesh_file_name[255]
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
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
Fix dofs on entities.
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:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorDouble & getFieldData() const
get dofs values
Basic algebra on fields.
Definition FieldBlas.hpp:21
structure to get information form mofem into EntitiesFieldData
Interface for managing meshsets containing materials and boundary conditions.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode buildProblem()
Build problem.
Definition Simple.cpp:597
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
Definition Simple.hpp:348
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition Simple.cpp:380
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:194
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition Simple.cpp:571
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
std::vector< std::string > & getOtherFiniteElements()
Get the Other Finite Elements.
Definition Simple.hpp:400
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:667
MoFEMErrorCode buildFields()
Build fields.
Definition Simple.cpp:482
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:473
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:282
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition Simple.cpp:452
const std::string getDomainFEName() const
Get the Domain FE Name.
Definition Simple.hpp:341
static double volumeLengthQuality(const double *coords)
Calculate tetrahedron volume length quality.
Definition Tools.cpp:15
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce 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.

Variable Documentation

◆ edges_block_set

int edges_block_set = 1
Examples
mesh_cut.cpp, and mesh_smoothing.cpp.

Definition at line 24 of file mesh_smoothing.cpp.

◆ flg_myfile

PetscBool flg_myfile = PETSC_TRUE
Examples
mesh_cut.cpp, and mesh_smoothing.cpp.

Definition at line 22 of file mesh_smoothing.cpp.

◆ gamma_factor

double gamma_factor = 0.8
Examples
mesh_smoothing.cpp.

Definition at line 28 of file mesh_smoothing.cpp.

◆ help

char help[] = "mesh smoothing\n\n"
static

Definition at line 20 of file mesh_smoothing.cpp.

◆ mesh_file_name

char mesh_file_name[255]

◆ output_vtk

PetscBool output_vtk = PETSC_TRUE
Examples
mesh_cut.cpp, mesh_smoothing.cpp, and split_sideset.cpp.

Definition at line 26 of file mesh_smoothing.cpp.

◆ tol

double tol = 0.1

◆ vertex_block_set

int vertex_block_set = 2
Examples
mesh_cut.cpp, and mesh_smoothing.cpp.

Definition at line 25 of file mesh_smoothing.cpp.