v0.15.5
Loading...
Searching...
No Matches
MeshSmoothingOpFactory.hpp
Go to the documentation of this file.
1#ifndef __MESH_SMOOTHING_OP_FACTORY_HPP__
2#define __MESH_SMOOTHING_OP_FACTORY_HPP__
3
4#ifndef WITH_ADOL_C
5#error "MoFEM need to be compiled with ADOL-C"
6#endif
7
8#include <MoFEM.hpp>
9
11#include <AdolCOps.hpp>
12#include <AdolCElastic.hpp>
13
14#include <algorithm>
15#include <limits>
16#include <vector>
17
19
20struct Tags {
21 int smoothingTag = 1;
22 int surfaceTag = 2;
23 int edgeTag = 3;
24};
25
26struct MinQualityOp;
27
28struct Config {
29 std::string geometryField = "MESH_NODE_POSITIONS";
30 std::string meshReferenceField = "NONE";
31 std::string lambdaSurfaceField = "LAMBDA_SURFACE";
32 std::string lambdaEdgeField = "LAMBDA_EDGE";
34 bool enableEdgeSliding = false;
35 bool fixEdgeBlocks = true;
37 double gamma = 0.0;
39};
40
41struct Handles {
43 boost::shared_ptr<AdolCOps::PhysicalEquations> volQualityPhysicalPtr;
44 boost::shared_ptr<ForcesAndSourcesCore> feSmootherRhs;
45 boost::shared_ptr<ForcesAndSourcesCore> feSmootherLhs;
46 boost::shared_ptr<SurfaceSlidingConstrains> surfaceConstraint;
47 boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
49 boost::shared_ptr<EdgeSlidingConstrains> edgeConstraint;
50 boost::shared_ptr<MoFEM::VolumeElementForcesAndSourcesCore> minQualityFe;
51 boost::shared_ptr<double> minQualityValue;
52};
53
55
57
58 MinQualityOp(const std::string &field_name,
59 double *min_quality_ptr)
62 minQualityPtr(min_quality_ptr) {}
63
64 MoFEMErrorCode doWork(int, EntityType type, EntitiesFieldData::EntData &data) {
66 if (type != MBVERTEX)
68 const double q = MoFEM::Tools::volumeLengthQuality(&*data.getFieldData().begin());
69 *minQualityPtr = std::min(*minQualityPtr, q);
71 }
72};
73
74inline MoFEMErrorCode createOperators(MoFEM::Interface &m_field,
75 const Config &cfg,
76 const Range &fixed_vertices,
77 const Range &edge_entities,
78 Handles &handles) {
80
81 handles.m_field_ptr = &m_field;
82
84 using DomainEleOp = DomainEle::UserDataOperator;
85 handles.feSmootherRhs = boost::make_shared<DomainEle>(m_field);
86 handles.feSmootherLhs = boost::make_shared<DomainEle>(m_field);
87
92 cfg.tags.smoothingTag));
93 CHKERR handles.volQualityPhysicalPtr->getOptions(&m_field);
94 CHKERR handles.volQualityPhysicalPtr->recordTape();
95
96 using DomainEleOp =
98 auto &rhs_ops = handles.feSmootherRhs->getOpPtrVector();
99 rhs_ops.clear();
101 template opRhsFactory<PETSC, GAUSS, DomainEleOp>(
102 m_field, rhs_ops, cfg.geometryField, handles.volQualityPhysicalPtr);
103
104 auto &lhs_ops = handles.feSmootherLhs->getOpPtrVector();
105 lhs_ops.clear();
107 template opLhsFactory<PETSC, GAUSS, DomainEleOp>(
108 m_field, lhs_ops, cfg.geometryField, handles.volQualityPhysicalPtr);
109
110 handles.minQualityFe =
111 boost::make_shared<MoFEM::VolumeElementForcesAndSourcesCore>(m_field);
112 handles.minQualityValue =
113 boost::make_shared<double>(std::numeric_limits<double>::max());
114
115 handles.minQualityFe->getOpPtrVector().push_back(
116 new MinQualityOp(cfg.geometryField, handles.minQualityValue.get()));
117
118 if (!fixed_vertices.empty()) {
119 auto *prb_mng = m_field.getInterface<ProblemsManager>();
120 auto *simple = m_field.getInterface<Simple>();
121
122 CHKERR prb_mng->removeDofsOnEntities(simple->getProblemName(),
123 cfg.geometryField, fixed_vertices, 0,
124 3);
125
126 if (cfg.enableSurfaceSliding && m_field.check_field(cfg.lambdaSurfaceField)) {
127 CHKERR prb_mng->removeDofsOnEntities(simple->getProblemName(),
128 cfg.lambdaSurfaceField, fixed_vertices,
129 0, 1);
130 }
131
132 if (cfg.enableEdgeSliding && !cfg.fixEdgeBlocks &&
133 m_field.check_field(cfg.lambdaEdgeField)) {
134 CHKERR prb_mng->removeDofsOnEntities(simple->getProblemName(),
135 cfg.lambdaEdgeField, fixed_vertices,
136 0, 1);
137 }
138 }
139
140 if (cfg.enableSurfaceSliding) {
141 handles.surfaceOrientation = boost::make_shared<
143 handles.surfaceConstraint = boost::shared_ptr<SurfaceSlidingConstrains>(
144 handles.surfaceOrientation,
145 new SurfaceSlidingConstrains(m_field, *handles.surfaceOrientation));
146 CHKERR handles.surfaceConstraint->setOperators(cfg.tags.surfaceTag,
148 cfg.geometryField);
149 }
150
151 if (cfg.enableEdgeSliding && !edge_entities.empty() && !cfg.fixEdgeBlocks) {
152 Range tets;
153 CHKERR m_field.get_moab().get_entities_by_type(0, MBTET, tets);
154 Range skin_faces;
155 Skinner skin(&m_field.get_moab());
156 CHKERR skin.find_skin(0, tets, false, skin_faces);
157 if (skin_faces.empty()) {
158 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
159 "No skin faces found for edge constrains");
160 }
161 handles.edgeConstraint = boost::make_shared<EdgeSlidingConstrains>(m_field);
162 CHKERR handles.edgeConstraint->setOperators(
163 cfg.tags.edgeTag, edge_entities, skin_faces, cfg.lambdaEdgeField,
164 cfg.geometryField);
165 }
166
168}
169
170inline MoFEMErrorCode addToSnesDM(DM dm, const Config &cfg,
171 const Handles &handles) {
173 boost::shared_ptr<ForcesAndSourcesCore> null;
174
175 CHKERR DMMoFEMSNESSetFunction(dm, "SMOOTHING", handles.feSmootherRhs, null,
176 null);
177 CHKERR DMMoFEMSNESSetJacobian(dm, "SMOOTHING", handles.feSmootherLhs, null,
178 null);
179
180 if (cfg.enableSurfaceSliding && handles.surfaceConstraint) {
181 CHKERR DMMoFEMSNESSetFunction(dm, "SURFACE_SLIDING",
182 handles.surfaceConstraint->feRhsPtr, null,
183 null);
184 CHKERR DMMoFEMSNESSetJacobian(dm, "SURFACE_SLIDING",
185 handles.surfaceConstraint->feLhsPtr, null,
186 null);
187 }
188
189 if (cfg.enableEdgeSliding && !cfg.fixEdgeBlocks && handles.edgeConstraint) {
190 CHKERR DMMoFEMSNESSetFunction(dm, "EDGE_SLIDING",
191 handles.edgeConstraint->feRhsPtr, null, null);
192 CHKERR DMMoFEMSNESSetJacobian(dm, "EDGE_SLIDING",
193 handles.edgeConstraint->feLhsPtr, null, null);
194 }
195
197}
198
199inline MoFEMErrorCode evaluateMinQuality(DM dm,
200 const std::string &domain_fe_name,
201 const Handles &handles,
202 double &min_quality) {
204
205 if (!handles.minQualityFe || !handles.minQualityValue) {
206 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
207 "Min quality operators are not initialised");
208 }
209
210 *handles.minQualityValue = std::numeric_limits<double>::max();
211 CHKERR DMoFEMLoopFiniteElements(dm, domain_fe_name.c_str(),
212 handles.minQualityFe);
213
214 MoFEM::Interface *m_field_ptr = nullptr;
215 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
216 double local_min = *handles.minQualityValue;
217 CHKERR MPI_Allreduce(&local_min, &min_quality, 1, MPI_DOUBLE, MPI_MIN,
218 m_field_ptr->get_comm());
219
221}
222
223inline MoFEMErrorCode updateBarrierGamma(Handles &handles, double gamma) {
225 char gamma_value[64];
226 CHKERR PetscSNPrintf(gamma_value, sizeof(gamma_value), "%.17g", gamma);
227 CHKERR PetscOptionsSetValue(NULL, "-volume_length_gamma", gamma_value);
228 CHKERR handles.volQualityPhysicalPtr->getOptions(handles.m_field_ptr);
229 CHKERR handles.volQualityPhysicalPtr->recordTape();
231}
232
233inline MoFEMErrorCode createSnes(DM dm, SmartPetscObj<SNES> &snes) {
235 MoFEM::Interface *m_field_ptr = nullptr;
236 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
237 snes = createSNES(m_field_ptr->get_comm());
238 CHKERR SNESSetFromOptions(snes);
239 CHKERR SNESSetDM(snes, dm);
240 CHKERR SNESSetUp(snes);
242}
243
244inline MoFEMErrorCode solveSnes(
245 DM dm, SNES snes,
246 const std::vector<std::string> &zero_vertex_fields = {}) {
248 if (!zero_vertex_fields.empty()) {
249 MoFEM::Interface *m_field_ptr = nullptr;
250 CHKERR DMoFEMGetInterfacePtr(dm, &m_field_ptr);
251 auto field_blas = m_field_ptr->getInterface<FieldBlas>();
252 for (const auto &field_name : zero_vertex_fields) {
253 CHKERR field_blas->setField(0, MBVERTEX, field_name);
254 }
255 }
256
257 auto f = createDMVector(dm);
258 auto d = vectorDuplicate(f);
259 CHKERR DMoFEMMeshToLocalVector(dm, d, INSERT_VALUES, SCATTER_FORWARD);
260 CHKERR SNESSolve(snes, f, d);
261 CHKERR DMoFEMMeshToGlobalVector(dm, d, INSERT_VALUES, SCATTER_REVERSE);
263}
264
265} // namespace MeshSmoothingOps
266
267#endif // __MESH_SMOOTHING_OP_FACTORY_HPP__
ForcesAndSourcesCore::UserDataOperator UserDataOperator
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
#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 ...
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode, RowColData rc=RowColData::COL)
set ghosted vector values on all existing mesh entities
Definition DMMoFEM.cpp:525
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode, RowColData rc=RowColData::COL)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
auto createDMVector(DM dm, RowColData rc=RowColData::COL)
Get smart vector from DM.
Definition DMMoFEM.hpp:1237
virtual bool check_field(const std::string &name) const =0
check if field is in database
boost::shared_ptr< ADolCData > createADolCDataPtr()
Definition AdolOps.cpp:245
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
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)
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
constexpr auto field_name
static auto setTagName(std::string name, int tag=-1)
Definition AdolCOps.hpp:16
AdolCOps::VolumeLengthQualityType qualityType
boost::shared_ptr< double > minQualityValue
boost::shared_ptr< ForcesAndSourcesCore > feSmootherLhs
boost::shared_ptr< AdolCOps::PhysicalEquations > volQualityPhysicalPtr
boost::shared_ptr< EdgeSlidingConstrains > edgeConstraint
boost::shared_ptr< MoFEM::VolumeElementForcesAndSourcesCore > minQualityFe
boost::shared_ptr< SurfaceSlidingConstrains > surfaceConstraint
boost::shared_ptr< ForcesAndSourcesCore > feSmootherRhs
boost::shared_ptr< SurfaceSlidingConstrains::DriverElementOrientation > surfaceOrientation
MinQualityOp(const std::string &field_name, double *min_quality_ptr)
MoFEMErrorCode doWork(int, EntityType type, EntitiesFieldData::EntData &data)
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
Deprecated interface functions.
@ OPROW
operator doWork function is executed on FE rows
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.
Implementing surface sliding constrains.