v0.15.5
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 <PostProc.hpp>
#include <MeshSmoothingOpFactory.hpp>

Go to the source code of this file.

Functions

static MoFEMErrorCode extract_mesh_edges (MoFEM::Interface &m_field, Range &output_edges, Range &output_vertices)
 
int main (int argc, char *argv[])
 

Variables

static char help [] = "mesh smoothing\n\n"
 
int edges_block_set = 1
 
int vertex_block_set = 2
 
PetscBool output_vtk = PETSC_TRUE
 
double tol = 0.1
 
double gamma_factor = 0.8
 
constexpr bool fix_edge_blocks = true
 

Detailed Description

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

Definition in file mesh_smoothing.cpp.

Function Documentation

◆ extract_mesh_edges()

static MoFEMErrorCode extract_mesh_edges ( MoFEM::Interface m_field,
Range output_edges,
Range output_vertices 
)
static
Examples
mofem/tools/mesh_smoothing.cpp.

Definition at line 349 of file mesh_smoothing.cpp.

351 {
353
354 Range fineFeatureTris;
355 FTENSOR_INDEX(3, i);
356 Range vols;
357 CHKERR m_field.get_moab().get_entities_by_dimension(0, 3, vols);
358 if (vols.empty()) {
359 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
360 "no 3d entities found in the mesh");
361 }
362
363 Skinner skinner(&m_field.get_moab());
364 Range skin_tris;
365 CHKERR skinner.find_skin(0, vols, false, skin_tris);
366 if (skin_tris.empty()) {
367 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
368 "no skin triangles found in the mesh");
369 }
370
371 auto tri_normal = [&](const EntityHandle tri,
374 CHKERR m_field.getInterface<Tools>()->getTriNormal(tri, &normal(0));
375 normal.normalize();
377 };
378
379 const double fine_planar_cos = std::cos(1e-6); //FIXME: parameter?
380
381 Range skin_edges;
382 CHKERR m_field.get_moab().get_adjacencies(skin_tris, 1, true, skin_edges,
383 moab::Interface::UNION);
384
385 for (auto edge : skin_edges) {
386 Range adj_tris;
387 CHKERR m_field.get_moab().get_adjacencies(&edge, 1, 2, false, adj_tris,
388 moab::Interface::UNION);
389 adj_tris = intersect(adj_tris, skin_tris);
390 if (adj_tris.size() != 2) {
391 continue;
392 }
393
396 CHKERR tri_normal(adj_tris[0], n0);
397 CHKERR tri_normal(adj_tris[1], n1);
398 const double dot = n0(i) * n1(i);
399 if (std::abs(dot) <= fine_planar_cos) {
400 fineFeatureTris.insert(adj_tris[0]);
401 fineFeatureTris.insert(adj_tris[1]);
402 output_edges.insert(edge);
403 }
404 }
405
406 std::map<EntityHandle, std::vector<EntityHandle>> vert_edges;
407 for (auto edge : output_edges) {
408 const EntityHandle *conn = nullptr;
409 int nconn = 0;
410 CHKERR m_field.get_moab().get_connectivity(edge, conn, nconn, false);
411 if (nconn != 2)
412 continue;
413 vert_edges[conn[0]].push_back(edge);
414 vert_edges[conn[1]].push_back(edge);
415 }
416
417 Range straight_edges;
418 for (auto edge : output_edges) {
419 const EntityHandle *conn = nullptr;
420 int nconn = 0;
421 CHKERR m_field.get_moab().get_connectivity(edge, conn, nconn, false);
422 if (nconn != 2)
423 continue;
424
425 bool keep_edge = true;
426 for (int vv = 0; vv < 2; ++vv) {
427 const auto v = conn[vv];
428 const auto it = vert_edges.find(v);
429 if (it == vert_edges.end()) {
430 keep_edge = false;
431 break;
432 }
433 const auto &inc_edges = it->second;
434 if (inc_edges.size() > 2) {
435 keep_edge = false;
436 break;
437 }
438 if (inc_edges.size() == 2) {
439 const auto other_edge =
440 (inc_edges[0] == edge) ? inc_edges[1] : inc_edges[0];
441 const EntityHandle *conn_other = nullptr;
442 int nconn_other = 0;
443 CHKERR m_field.get_moab().get_connectivity(other_edge, conn_other,
444 nconn_other, false);
445 if (nconn_other != 2) {
446 keep_edge = false;
447 break;
448 }
449 const EntityHandle v_other =
450 (conn_other[0] == v) ? conn_other[1] : conn_other[0];
451
452 double pv[3], pa[3], pb[3];
453 CHKERR m_field.get_moab().get_coords(&v, 1, pv);
454 EntityHandle va = conn[1 - vv];
455 CHKERR m_field.get_moab().get_coords(&va, 1, pa);
456 CHKERR m_field.get_moab().get_coords(&v_other, 1, pb);
457
461 &pv[2]);
463 &pa[2]);
465 &pb[2]);
466 a(i) = t_pa(i) - t_pv(i);
467 b(i) = t_pb(i) - t_pv(i);
468 const double na = a.l2();
469 const double nb = b.l2();
470 if (na == 0 || nb == 0) {
471 keep_edge = false;
472 break;
473 }
474 const double dot = (a(i) * b(i)) / (na * nb);
475 if (std::abs(dot) < fine_planar_cos) {
476 keep_edge = false;
477 break;
478 }
479 }
480 }
481
482 if (keep_edge) {
483 straight_edges.insert(edge);
484 }
485 }
486
487 output_edges = straight_edges;
488 output_vertices.clear();
489 if (!output_edges.empty()) {
490 std::map<EntityHandle, int> vert_count;
491 for (auto edge : output_edges) {
492 const EntityHandle *conn = nullptr;
493 int nconn = 0;
494 CHKERR m_field.get_moab().get_connectivity(edge, conn, nconn, false);
495 if (nconn != 2)
496 continue;
497 vert_count[conn[0]]++;
498 vert_count[conn[1]]++;
499 }
500 Range filtered_edges;
501 for (auto edge : output_edges) {
502 const EntityHandle *conn = nullptr;
503 int nconn = 0;
504 CHKERR m_field.get_moab().get_connectivity(edge, conn, nconn, false);
505 if (nconn != 2)
506 continue;
507 const int d0 = vert_count[conn[0]];
508 const int d1 = vert_count[conn[1]];
509 if (!(d0 == 1 && d1 == 1)) {
510 filtered_edges.insert(edge);
511 for (int vv = 0; vv < 2; ++vv) {
512 const auto v = conn[vv];
513 auto it = vert_count.find(v);
514 if (it != vert_count.end() && it->second == 1) {
515 output_vertices.insert(v);
516 }
517 }
518 }
519 }
520 output_edges = filtered_edges;
521 }
522
523 if (!fineFeatureTris.empty())
525 "tri_features_test.vtk", fineFeatureTris);
526 if (!output_edges.empty())
527 CHKERR CutMeshInterface::SaveData(m_field.get_moab())("edges_test.vtk",
528 output_edges);
529 if (!output_vertices.empty())
530 CHKERR CutMeshInterface::SaveData(m_field.get_moab())("vertices_test.vtk",
531 output_vertices);
532
534}
#define FTENSOR_INDEX(DIM, I)
constexpr double a
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#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 ...
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
virtual moab::Interface & get_moab()=0
Auxiliary tools.
Definition Tools.hpp:19
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ main()

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

Definition at line 32 of file mesh_smoothing.cpp.

32 {
33
34 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
35
36 try {
37
38 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Mesh cut options", "none");
39 CHKERR PetscOptionsInt("-edges_block_set", "edges side set", "",
40 edges_block_set, &edges_block_set, PETSC_NULLPTR);
41 CHKERR PetscOptionsInt("-vertex_block_set", "vertex side set", "",
42 vertex_block_set, &vertex_block_set, PETSC_NULLPTR);
43 CHKERR PetscOptionsBool("-output_vtk", "if true outout vtk file", "",
44 output_vtk, &output_vtk, PETSC_NULLPTR);
45 CHKERR PetscOptionsScalar("-quality_reduction_tol", "",
46 "Tolerance of quality reduction", tol, &tol,
47 PETSC_NULLPTR);
48 CHKERR PetscOptionsScalar("-gamma_factor", "",
49 "Gamma factor", gamma_factor, &gamma_factor,
50 PETSC_NULLPTR);
51 PetscOptionsEnd();
52
53 // Create MoAB database
54 moab::Core moab_core; // create database
55 moab::Interface &moab = moab_core; // create interface to database
56 // Create MoFEM database and link it to MoAB
57 MoFEM::Core mofem_core(moab); // create database
58 MoFEM::Interface &m_field = mofem_core; // create interface to database
59 // Register DM Manager
60 CHKERR DMRegister_MoFEM("DMMOFEM"); // register MoFEM DM in PETSc
61
62 // Get simple interface is simplified version enabling quick and
63 // easy construction of problem.
64 Simple *simple_interface = m_field.getInterface<Simple>();
65 Range edges, vertices, fixed_vertex;
66
67 // Build problem with simple interface
68 {
69 // Get options for simple interface from command line
70 CHKERR simple_interface->getOptions();
71 // Load mesh file to database
72 CHKERR simple_interface->loadFile();
73
74 // Add domain field "U" in space H1 and Legendre base, Ainsworth recipe is
75 // used to construct base functions.
76 CHKERR simple_interface->addDomainField("MESH_NODE_POSITIONS", H1,
78 // Add Lagrange multiplier field on body boundary
79 CHKERR simple_interface->addBoundaryField("LAMBDA_SURFACE", H1,
81
82 // Set fields order domain and boundary fields.
83 CHKERR simple_interface->setFieldOrder("MESH_NODE_POSITIONS",
84 1); // to approximate function
85 CHKERR simple_interface->setFieldOrder("LAMBDA_SURFACE",
86 1); // to Lagrange multipliers
87
88 simple_interface->getDomainFEName() = "SMOOTHING";
89 simple_interface->getBoundaryFEName() = "SURFACE_SLIDING";
90
91 // Other fields and finite elements added to database directly
95 ->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
96 true);
97
98 } else {
99 CHKERR extract_mesh_edges(m_field, edges, fixed_vertex);
100 }
104 ->getEntitiesByDimension(vertex_block_set, BLOCKSET, 0,
105 fixed_vertex, true);
106 }
107
108 CHKERR m_field.get_moab().get_connectivity(edges, vertices, true);
109
110
111 if(!edges.empty() && !fix_edge_blocks) {
112 // Declare approximation fields
113 CHKERR m_field.add_field("LAMBDA_EDGE", H1, AINSWORTH_LOBATTO_BASE, 2,
114 MB_TAG_SPARSE, MF_ZERO);
115 CHKERR m_field.add_ents_to_field_by_type(edges, MBEDGE,
116 "LAMBDA_EDGE");
118 ->synchroniseFieldEntities("LAMBDA_EDGE");
119 CHKERR m_field.set_field_order(0, MBVERTEX, "LAMBDA_EDGE", 1);
120
121 CHKERR m_field.add_finite_element("EDGE_SLIDING");
123 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
125 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
127 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
128 CHKERR m_field.modify_finite_element_add_field_row("EDGE_SLIDING",
129 "LAMBDA_EDGE");
130 CHKERR m_field.modify_finite_element_add_field_col("EDGE_SLIDING",
131 "LAMBDA_EDGE");
132 CHKERR m_field.modify_finite_element_add_field_data("EDGE_SLIDING",
133 "LAMBDA_EDGE");
134 CHKERR m_field.add_ents_to_finite_element_by_type(edges, MBEDGE,
135 "EDGE_SLIDING");
136 simple_interface->getOtherFiniteElements().push_back("EDGE_SLIDING");
137 }
138
139 CHKERR simple_interface->defineFiniteElements();
140 CHKERR simple_interface->defineProblem();
141 CHKERR simple_interface->buildFields();
142 CHKERR simple_interface->buildFiniteElements();
143 CHKERR simple_interface->buildProblem();
144
145 // Remove vertices form LAMBDA_SURFACE which are on the edges
146 if (!edges.empty()) {
147 Range edges;
148 CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
149 edges_block_set, BLOCKSET, 1, edges, true);
150 Range verts;
151 CHKERR m_field.get_moab().get_connectivity(edges, verts, true);
152 Simple *simple_interface = m_field.getInterface<Simple>();
153 auto prb_mng = m_field.getInterface<ProblemsManager>();
154 CHKERR prb_mng->removeDofsOnEntities(simple_interface->getProblemName(),
155 "LAMBDA_SURFACE", verts, 0, 1);
156 if (fix_edge_blocks) {
157 CHKERR prb_mng->removeDofsOnEntities(
158 simple_interface->getProblemName(), "MESH_NODE_POSITIONS", verts,
159 0, 3);
160 }
161 }
162 }
163
164 MeshSmoothingOps::Config smoothing_cfg;
165 smoothing_cfg.geometryField = "MESH_NODE_POSITIONS";
166 smoothing_cfg.meshReferenceField = "NONE";
167 smoothing_cfg.lambdaSurfaceField = "LAMBDA_SURFACE";
168 smoothing_cfg.lambdaEdgeField = "LAMBDA_EDGE";
169 smoothing_cfg.enableSurfaceSliding = true;
170 smoothing_cfg.enableEdgeSliding = !fix_edge_blocks;
171 smoothing_cfg.fixEdgeBlocks = fix_edge_blocks;
173 smoothing_cfg.gamma = 0.0;
174
175 MeshSmoothingOps::Handles smoothing_handles;
176 CHKERR MeshSmoothingOps::createOperators(m_field, smoothing_cfg, fixed_vertex,
177 edges, smoothing_handles);
178
179 DM dm;
180 CHKERR simple_interface->getDM(&dm);
181 CHKERR MeshSmoothingOps::addToSnesDM(dm, smoothing_cfg, smoothing_handles);
182
183 struct Solve {
184
185 MoFEMErrorCode operator()(DM dm) const {
187 auto solver = SmartPetscObj<SNES>();
189
190 auto hook = [&](SNES snes, Vec x, Vec F,
191 boost::shared_ptr<CacheTuple> cache_ptr,
192 void *ctx) -> MoFEMErrorCode {
194
197
198 auto *m_field_ptr = getInterfacePtr(dm);
199 auto post_proc_fe =
200 boost::make_shared<PostProcEleDomain>(*m_field_ptr);
201 auto &pip = post_proc_fe->getOpPtrVector();
202 auto X_ptr = boost::make_shared<MatrixDouble>();
203 pip.push_back(new OpCalculateVectorFieldValues<3>(
204 "MESH_NODE_POSITIONS", X_ptr));
205 auto res_X_ptr = boost::make_shared<MatrixDouble>();
206 pip.push_back(new OpCalculateVectorFieldValues<3>(
207 "MESH_NODE_POSITIONS", res_X_ptr, SmartPetscObj<Vec>(F, true)));
208 auto lambda_ptr = boost::make_shared<VectorDouble>();
209 pip.push_back(
210 new OpCalculateScalarFieldValues("LAMBDA_SURFACE", lambda_ptr));
211 auto res_lambda_ptr = boost::make_shared<VectorDouble>();
212 pip.push_back(new OpCalculateScalarFieldValues(
213 "LAMBDA_SURFACE", res_lambda_ptr, SmartPetscObj<Vec>(F, true)));
214
216 pip.push_back(
217
218 new OpPPMap(
219
220 post_proc_fe->getPostProcMesh(),
221 post_proc_fe->getMapGaussPts(),
222
223 {{"LAMBDA_SURFACE", lambda_ptr},
224 {"RES_LAMBDA_SURFACE", res_lambda_ptr}},
225
226 {{"MESH_NODE_POSITIONS", X_ptr},
227
228 {"RES_MESH_NODE_POSITIONS", res_X_ptr}},
229
230 {},
231
232 {}
233
234 )
235
236 );
237
238 CHKERR DMoFEMLoopFiniteElements(dm, "SURFACE_SLIDING", post_proc_fe);
239 CHKERR post_proc_fe->writeFile("debug_smoothing.h5m");
241 };
242
243 auto snes_ctx_ptr = getDMSnesCtx(dm);
244 snes_ctx_ptr->getRhsHook() = hook;
245
246 CHKERR MeshSmoothingOps::solveSnes(dm, solver, {"LAMBDA_SURFACE"});
247
249 }
250
251 MoFEMErrorCode setCoordsFromField(DM dm) const {
253 MoFEM::Interface *m_field_ptr;
254 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
255 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
256 "MESH_NODE_POSITIONS", it)) {
257 if (it->get()->getEntType() != MBVERTEX)
258 continue;
259 VectorDouble3 coords(3);
260 for(int dd = 0;dd!=3;++dd)
261 coords[dd] = it->get()->getEntFieldData()[dd];
262 EntityHandle ent = it->get()->getEnt();
263 CHKERR m_field_ptr->get_moab().set_coords(&ent, 1, &*coords.begin());
264 }
266 }
267
268 MoFEMErrorCode setFieldFromCoords(DM dm) const {
270 MoFEM::Interface *m_field_ptr;
271 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
272 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
273 "MESH_NODE_POSITIONS", it)) {
274 if (it->get()->getEntType() != MBVERTEX)
275 continue;
276 EntityHandle ent = it->get()->getEnt();
277 VectorDouble3 coords(3);
278 CHKERR m_field_ptr->get_moab().get_coords(&ent, 1, &*coords.begin());
279 for(int dd = 0;dd!=3;++dd)
280 it->get()->getEntFieldData()[dd] = coords[dd];
281 }
283 }
284
285 };
286
287 Solve solve;
288 CHKERR solve.setFieldFromCoords(dm);
289
290 double min_quality = 0;
292 dm, simple_interface->getDomainFEName(), smoothing_handles, min_quality);
293 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f\n", min_quality);
294
295 double gamma = min_quality > 0 ? gamma_factor * min_quality
296 : min_quality / gamma_factor;
297 CHKERR MeshSmoothingOps::updateBarrierGamma(smoothing_handles, gamma);
298
299 double min_quality_p, eps;
300 do {
301
302 min_quality_p = min_quality;
303
304 CHKERR solve(dm);
305
306 CHKERR solve.setCoordsFromField(dm);
308 dm, simple_interface->getDomainFEName(), smoothing_handles,
309 min_quality);
310
311 eps = (min_quality - min_quality_p) / min_quality;
312 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f eps = %4.3f\n",
313 min_quality, eps);
314
315 double gamma = min_quality > 0 ? gamma_factor * min_quality
316 : min_quality / gamma_factor;
317 CHKERR MeshSmoothingOps::updateBarrierGamma(smoothing_handles, gamma);
318
319 } while (eps > tol);
320
321 // if (m_field.getInterface<MeshsetsManager>()->checkMeshset(edges_block_set,
322 // BLOCKSET)) {
323 // Range edges;
324 // CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
325 // edges_block_set, BLOCKSET, 1, edges, true);
326
327 // Range tets;
328 // CHKERR moab.get_entities_by_type(0,MBTET,tets);
329 // Skinner skin(&moab);
330 // Range skin_faces; // skin faces from 3d ents
331 // CHKERR skin.find_skin(0, tets, false, skin_faces);
332
333 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::setTags(moab, edges,
334 // skin_faces);
335 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
336 // moab, "out_edges.vtk", edges);
337 // }
338
339 if (output_vtk)
340 CHKERR m_field.getInterface<BitRefManager>()->writeBitLevelByType(
341 BitRefLevel().set(0), BitRefLevel().set(), MBTET, "out.vtk", "VTK",
342 "");
343 }
345
347}
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
@ H1
continuous field
Definition definitions.h:85
@ BLOCKSET
PostProcEleByDim< SPACE_DIM >::PostProcEleDomain PostProcEleDomain
@ F
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
#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
const FTensor::Tensor2< T, Dim, Dim > Vec
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
MoFEMErrorCode updateBarrierGamma(Handles &handles, double gamma)
MoFEMErrorCode evaluateMinQuality(DM dm, const std::string &domain_fe_name, const Handles &handles, double &min_quality)
MoFEMErrorCode createSnes(DM dm, SmartPetscObj< SNES > &snes)
MoFEMErrorCode solveSnes(DM dm, SNES snes, const std::vector< std::string > &zero_vertex_fields={})
MoFEMErrorCode createOperators(MoFEM::Interface &m_field, const Config &cfg, const Range &fixed_vertices, const Range &edge_entities, Handles &handles)
MoFEMErrorCode addToSnesDM(DM dm, const Config &cfg, const Handles &handles)
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
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:1171
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1265
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
AdolCOps::VolumeLengthQualityType qualityType
Managing BitRefLevels.
Managing BitRefLevels.
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.
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
double tol
PetscBool output_vtk
static char help[]
constexpr bool fix_edge_blocks
static MoFEMErrorCode extract_mesh_edges(MoFEM::Interface &m_field, Range &output_edges, Range &output_vertices)
double gamma_factor
int vertex_block_set
int edges_block_set

Variable Documentation

◆ edges_block_set

int edges_block_set = 1

◆ fix_edge_blocks

constexpr bool fix_edge_blocks = true
constexpr
Examples
mofem/tools/mesh_smoothing.cpp.

Definition at line 26 of file mesh_smoothing.cpp.

◆ gamma_factor

double gamma_factor = 0.8
Examples
mofem/tools/mesh_smoothing.cpp.

Definition at line 24 of file mesh_smoothing.cpp.

◆ help

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

Definition at line 18 of file mesh_smoothing.cpp.

◆ output_vtk

PetscBool output_vtk = PETSC_TRUE

◆ tol

double tol = 0.1

◆ vertex_block_set

int vertex_block_set = 2