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
12#include <PostProc.hpp>
13
15
16using namespace MoFEM;
17
18static char help[] = "mesh smoothing\n\n";
19
22PetscBool output_vtk = PETSC_TRUE;
23double tol = 0.1;
24double gamma_factor = 0.8;
25
26constexpr bool fix_edge_blocks = true;
27
28int main(int argc, char *argv[]) {
29
30 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
31
32 try {
33
34 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Mesh cut options", "none");
35 CHKERR PetscOptionsInt("-edges_block_set", "edges side set", "",
36 edges_block_set, &edges_block_set, PETSC_NULLPTR);
37 CHKERR PetscOptionsInt("-vertex_block_set", "vertex side set", "",
38 vertex_block_set, &vertex_block_set, PETSC_NULLPTR);
39 CHKERR PetscOptionsBool("-output_vtk", "if true outout vtk file", "",
40 output_vtk, &output_vtk, PETSC_NULLPTR);
41 CHKERR PetscOptionsScalar("-quality_reduction_tol", "",
42 "Tolerance of quality reduction", tol, &tol,
43 PETSC_NULLPTR);
44 CHKERR PetscOptionsScalar("-gamma_factor", "",
45 "Gamma factor", gamma_factor, &gamma_factor,
46 PETSC_NULLPTR);
47 PetscOptionsEnd();
48
49 // Create MoAB database
50 moab::Core moab_core; // create database
51 moab::Interface &moab = moab_core; // create interface to database
52 // Create MoFEM database and link it to MoAB
53 MoFEM::Core mofem_core(moab); // create database
54 MoFEM::Interface &m_field = mofem_core; // create interface to database
55 // Register DM Manager
56 CHKERR DMRegister_MoFEM("DMMOFEM"); // register MoFEM DM in PETSc
57
58 // Get simple interface is simplified version enabling quick and
59 // easy construction of problem.
60 Simple *simple_interface = m_field.getInterface<Simple>();
61 Range edges, vertices, fixed_vertex;
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
91 ->getEntitiesByDimension(edges_block_set, BLOCKSET, 1, edges,
92 true);
93
94 } else {
95 CHKERR Smoother::extractMeshEdges(m_field, edges, fixed_vertex);
96 }
100 ->getEntitiesByDimension(vertex_block_set, BLOCKSET, 0,
101 fixed_vertex, true);
102 }
103
104 CHKERR m_field.get_moab().get_connectivity(edges, vertices, true);
105
106
107 if(!edges.empty() && !fix_edge_blocks) {
108 // Declare approximation fields
109 CHKERR m_field.add_field("LAMBDA_EDGE", H1, AINSWORTH_LOBATTO_BASE, 2,
110 MB_TAG_SPARSE, MF_ZERO);
111 CHKERR m_field.add_ents_to_field_by_type(edges, MBEDGE,
112 "LAMBDA_EDGE");
114 ->synchroniseFieldEntities("LAMBDA_EDGE");
115 CHKERR m_field.set_field_order(0, MBVERTEX, "LAMBDA_EDGE", 1);
116
117 CHKERR m_field.add_finite_element("EDGE_SLIDING");
119 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
121 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
123 "EDGE_SLIDING", "MESH_NODE_POSITIONS");
124 CHKERR m_field.modify_finite_element_add_field_row("EDGE_SLIDING",
125 "LAMBDA_EDGE");
126 CHKERR m_field.modify_finite_element_add_field_col("EDGE_SLIDING",
127 "LAMBDA_EDGE");
128 CHKERR m_field.modify_finite_element_add_field_data("EDGE_SLIDING",
129 "LAMBDA_EDGE");
130 CHKERR m_field.add_ents_to_finite_element_by_type(edges, MBEDGE,
131 "EDGE_SLIDING");
132 simple_interface->getOtherFiniteElements().push_back("EDGE_SLIDING");
133 }
134
135 CHKERR simple_interface->defineFiniteElements();
136 CHKERR simple_interface->defineProblem();
137 CHKERR simple_interface->buildFields();
138 CHKERR simple_interface->buildFiniteElements();
139 CHKERR simple_interface->buildProblem();
140
141 // Remove vertices form LAMBDA_SURFACE which are on the edges
142 if (!edges.empty()) {
143 Range edges;
144 CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
145 edges_block_set, BLOCKSET, 1, edges, true);
146 Range verts;
147 CHKERR m_field.get_moab().get_connectivity(edges, verts, true);
148 Simple *simple_interface = m_field.getInterface<Simple>();
149 auto prb_mng = m_field.getInterface<ProblemsManager>();
150 CHKERR prb_mng->removeDofsOnEntities(simple_interface->getProblemName(),
151 "LAMBDA_SURFACE", verts, 0, 1);
152 if (fix_edge_blocks) {
153 CHKERR prb_mng->removeDofsOnEntities(
154 simple_interface->getProblemName(), "MESH_NODE_POSITIONS", verts,
155 0, 3);
156 }
157 }
158 }
159
160 MeshSmoothingOps::Config smoothing_cfg;
161 smoothing_cfg.geometryField = "MESH_NODE_POSITIONS";
162 smoothing_cfg.meshReferenceField = "NONE";
163 smoothing_cfg.lambdaSurfaceField = "LAMBDA_SURFACE";
164 smoothing_cfg.lambdaEdgeField = "LAMBDA_EDGE";
165 smoothing_cfg.enableSurfaceSliding = true;
166 smoothing_cfg.enableEdgeSliding = !fix_edge_blocks;
167 smoothing_cfg.fixEdgeBlocks = fix_edge_blocks;
168 smoothing_cfg.qualityType = BARRIER_AND_QUALITY;
169 smoothing_cfg.alpha = 1.0;
170 smoothing_cfg.gamma = 0.0;
171
172 MeshSmoothingOps::Handles smoothing_handles;
173 CHKERR MeshSmoothingOps::createOperators(m_field, smoothing_cfg, fixed_vertex,
174 edges, smoothing_handles);
175
176 DM dm;
177 CHKERR simple_interface->getDM(&dm);
178 CHKERR MeshSmoothingOps::addToSnesDM(dm, smoothing_cfg, smoothing_handles);
179
180 struct Solve {
181
182 MoFEMErrorCode operator()(DM dm) const {
184 auto solver = SmartPetscObj<SNES>();
186
187 auto hook = [&](SNES snes, Vec x, Vec F,
188 boost::shared_ptr<CacheTuple> cache_ptr,
189 void *ctx) -> MoFEMErrorCode {
191
194
195 auto *m_field_ptr = getInterfacePtr(dm);
196 auto post_proc_fe =
197 boost::make_shared<PostProcEleDomain>(*m_field_ptr);
198 auto &pip = post_proc_fe->getOpPtrVector();
199 auto X_ptr = boost::make_shared<MatrixDouble>();
200 pip.push_back(new OpCalculateVectorFieldValues<3>(
201 "MESH_NODE_POSITIONS", X_ptr));
202 auto res_X_ptr = boost::make_shared<MatrixDouble>();
203 pip.push_back(new OpCalculateVectorFieldValues<3>(
204 "MESH_NODE_POSITIONS", res_X_ptr, SmartPetscObj<Vec>(F, true)));
205 auto lambda_ptr = boost::make_shared<VectorDouble>();
206 pip.push_back(
207 new OpCalculateScalarFieldValues("LAMBDA_SURFACE", lambda_ptr));
208 auto res_lambda_ptr = boost::make_shared<VectorDouble>();
209 pip.push_back(new OpCalculateScalarFieldValues(
210 "LAMBDA_SURFACE", res_lambda_ptr, SmartPetscObj<Vec>(F, true)));
211
213 pip.push_back(
214
215 new OpPPMap(
216
217 post_proc_fe->getPostProcMesh(),
218 post_proc_fe->getMapGaussPts(),
219
220 {{"LAMBDA_SURFACE", lambda_ptr},
221 {"RES_LAMBDA_SURFACE", res_lambda_ptr}},
222
223 {{"MESH_NODE_POSITIONS", X_ptr},
224
225 {"RES_MESH_NODE_POSITIONS", res_X_ptr}},
226
227 {},
228
229 {}
230
231 )
232
233 );
234
235 CHKERR DMoFEMLoopFiniteElements(dm, "SURFACE_SLIDING", post_proc_fe);
236 CHKERR post_proc_fe->writeFile("debug_smoothing.h5m");
238 };
239
240 auto snes_ctx_ptr = getDMSnesCtx(dm);
241 snes_ctx_ptr->getRhsHook() = hook;
242
243 CHKERR MeshSmoothingOps::solveSnes(dm, solver, {"LAMBDA_SURFACE"});
244
246 }
247
248 MoFEMErrorCode setCoordsFromField(DM dm) const {
250 MoFEM::Interface *m_field_ptr;
251 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
252 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
253 "MESH_NODE_POSITIONS", it)) {
254 if (it->get()->getEntType() != MBVERTEX)
255 continue;
256 VectorDouble3 coords(3);
257 for(int dd = 0;dd!=3;++dd)
258 coords[dd] = it->get()->getEntFieldData()[dd];
259 EntityHandle ent = it->get()->getEnt();
260 CHKERR m_field_ptr->get_moab().set_coords(&ent, 1, &*coords.begin());
261 }
263 }
264
265 MoFEMErrorCode setFieldFromCoords(DM dm) const {
267 MoFEM::Interface *m_field_ptr;
268 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
269 for (_IT_GET_ENT_FIELD_BY_NAME_FOR_LOOP_(*m_field_ptr,
270 "MESH_NODE_POSITIONS", it)) {
271 if (it->get()->getEntType() != MBVERTEX)
272 continue;
273 EntityHandle ent = it->get()->getEnt();
274 VectorDouble3 coords(3);
275 CHKERR m_field_ptr->get_moab().get_coords(&ent, 1, &*coords.begin());
276 for(int dd = 0;dd!=3;++dd)
277 it->get()->getEntFieldData()[dd] = coords[dd];
278 }
280 }
281
282 };
283
284 Solve solve;
285 CHKERR solve.setFieldFromCoords(dm);
286
287 double min_quality = 0;
289 dm, simple_interface->getDomainFEName(), smoothing_handles, min_quality);
290 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f\n", min_quality);
291
292 double gamma = min_quality > 0 ? gamma_factor * min_quality
293 : min_quality / gamma_factor;
294 CHKERR MeshSmoothingOps::updateBarrierGamma(smoothing_handles, gamma);
295
296 double min_quality_p, eps;
297 do {
298
299 min_quality_p = min_quality;
300
301 CHKERR solve(dm);
302
303 CHKERR solve.setCoordsFromField(dm);
305 dm, simple_interface->getDomainFEName(), smoothing_handles,
306 min_quality);
307
308 eps = (min_quality - min_quality_p) / min_quality;
309 PetscPrintf(PETSC_COMM_WORLD, "Min quality = %4.3f eps = %4.3f\n",
310 min_quality, eps);
311
312 double gamma = min_quality > 0 ? gamma_factor * min_quality
313 : min_quality / gamma_factor;
314 CHKERR MeshSmoothingOps::updateBarrierGamma(smoothing_handles, gamma);
315
316 } while (eps > tol);
317
318 // if (m_field.getInterface<MeshsetsManager>()->checkMeshset(edges_block_set,
319 // BLOCKSET)) {
320 // Range edges;
321 // CHKERR m_field.getInterface<MeshsetsManager>()->getEntitiesByDimension(
322 // edges_block_set, BLOCKSET, 1, edges, true);
323
324 // Range tets;
325 // CHKERR moab.get_entities_by_type(0,MBTET,tets);
326 // Skinner skin(&moab);
327 // Range skin_faces; // skin faces from 3d ents
328 // CHKERR skin.find_skin(0, tets, false, skin_faces);
329
330 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::setTags(moab, edges,
331 // skin_faces);
332 // CHKERR EdgeSlidingConstrains::CalculateEdgeBase::saveEdges(
333 // moab, "out_edges.vtk", edges);
334 // }
335
336 if (output_vtk)
337 CHKERR m_field.getInterface<BitRefManager>()->writeBitLevelByType(
338 BitRefLevel().set(0), BitRefLevel().set(), MBTET, "out.vtk", "VTK",
339 "");
340 }
342
344}
Post-process fields on refined mesh.
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
@ 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.
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
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
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
VolumeLengthQualityType qualityType
Managing BitRefLevels.
Managing BitRefLevels.
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.
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
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
static MoFEMErrorCode extractMeshEdges(MoFEM::Interface &m_field, Range &output_edges, Range &output_vertices)
Definition Smoother.hpp:10
double tol
PetscBool output_vtk
static char help[]
constexpr bool fix_edge_blocks
double gamma_factor
int vertex_block_set
int edges_block_set
@ BARRIER_AND_QUALITY