6 #include "GmshConfig.h"
20 #if defined(HAVE_OPENGL)
23 #if defined(HAVE_VISUDEV)
36 #if defined(HAVE_VISUDEV)
38 {
GMSH_FULLRC,
"Element to draw quality",
nullptr, 0}
51 return sizeof(CurvedMeshOptions_Number) /
sizeof(
StringXNumber);
56 return &CurvedMeshOptions_Number[iopt];
61 return "Plugin(AnalyseMeshQuality) analyses the quality of the elements "
62 "of a given dimension in the current model. Depending on the input "
63 "parameters it computes the minimum of the Jacobian determinant (J), "
64 "the IGE quality measure (Inverse Gradient Error) and/or the ICN "
65 "quality measure (Condition Number). Statistics are printed and, "
66 "if requested, a model-based post-processing view is created for each "
67 "quality measure. The plugin can optionally hide elements by "
68 "comparing the measure to a prescribed threshold.\n"
70 "J is faster to compute but gives information only on element "
71 "validity while the other measures also give information on element "
72 "quality. The IGE measure is related to the error on the gradient of "
73 "the finite element solution. It is the scaled Jacobian for quads and "
74 "hexes and a new measure for triangles and tetrahedra. "
75 "The ICN measure is related to the condition number of the stiffness "
76 "matrix. (See the article \"Efficient computation of the minimum of "
77 "shape quality measures on curvilinear finite elements\" for "
82 "- `JacobianDeterminant': compute J?\n"
84 "- `IGEMeasure': compute IGE?\n"
86 "- `ICNMeasure': compute ICN?\n"
88 "- `HidingThreshold': hide all elements for which min(mu) is "
89 "strictly greater than (if `ThresholdGreater' == 1) or less than "
90 "(if `ThresholdGreater' == 0) the threshold, where mu is ICN if "
91 "`ICNMeasure' == 1, IGE if `IGEMeasure' == 1 or min(J)/max(J) if "
92 "`JacobianDeterminant' == 1.\n"
94 "- `CreateView': create a model-based view of min(J)/max(J), "
95 "min(IGE) and/or min(ICN)?\n"
97 "- `Recompute': force recomputation (set to 1 if the mesh has "
100 "- `DimensionOfElements': analyse elements of the given dimension if "
101 "equal to 1, 2 or 3; analyse 2D and 3D elements if equal to 4; "
102 "or analyse elements of the highest dimension if equal to -1.";
108 int computeJac =
static_cast<int>(CurvedMeshOptions_Number[0].
def);
109 int computeIGE =
static_cast<int>(CurvedMeshOptions_Number[1].
def);
110 int computeICN =
static_cast<int>(CurvedMeshOptions_Number[2].
def);
111 double threshold = CurvedMeshOptions_Number[3].
def;
112 bool thresholdGreater =
static_cast<bool>(CurvedMeshOptions_Number[4].
def);
113 bool createView =
static_cast<bool>(CurvedMeshOptions_Number[5].
def);
114 bool recompute =
static_cast<bool>(CurvedMeshOptions_Number[6].
def);
115 int askedDim =
static_cast<int>(CurvedMeshOptions_Number[7].
def);
117 #if defined(HAVE_VISUDEV)
118 _pwJac = computeJac / 2;
119 _pwIGE = computeIGE / 2;
120 _pwICN = computeICN / 2;
122 _numElementToScan =
static_cast<int>(CurvedMeshOptions_Number[8].
def);
124 _dataPViewJac.clear();
125 _dataPViewIGE.clear();
126 _dataPViewICN.clear();
129 if(askedDim < 0 || askedDim > 4) askedDim =
_m->
getDim();
131 if(recompute)
_clear(askedDim);
134 bool printStatJ =
false;
135 bool printStatS =
false;
136 bool printStatI =
false;
137 for(
int dim = 1; dim <= 3 && dim <=
_m->
getDim(); ++dim) {
138 if((askedDim == 4 && dim > 1) || dim == askedDim) {
141 Msg::StatusBar(
true,
"Computing Jacobian for %dD elements...", dim);
144 "Done computing Jacobian for %dD elements (Wall %gs, "
154 "Done computing IGE for %dD elements (Wall %gs, "
164 "Done computing ICN for %dD elements (Wall %gs, "
175 #if defined(HAVE_VISUDEV)
176 _createPViewPointwise();
180 PView *view =
nullptr;
182 for(
int dim = 1; dim <= 3; ++dim) {
183 if((askedDim == 4 && dim > 1) || dim == askedDim) {
186 std::map<int, std::vector<double> > dataPV;
187 for(std::size_t i = 0; i <
_data.size(); ++i) {
191 if(
_data[i].maxJ() > 0)
193 else if(
_data[i].maxJ() < 0)
195 dataPV[el->
getNum()].push_back(q);
199 std::stringstream name;
200 name <<
"minJ/maxJ " << dim <<
"D";
201 view =
new PView(name.str().c_str(),
"ElementData",
_m, dataPV);
206 std::map<int, std::vector<double> > dataPV;
207 for(std::size_t i = 0; i <
_data.size(); ++i) {
213 std::stringstream name;
214 name <<
"IGE " << dim <<
"D";
215 view =
new PView(name.str().c_str(),
"ElementData",
_m, dataPV);
220 std::map<int, std::vector<double> > dataPV;
221 for(std::size_t i = 0; i <
_data.size(); ++i) {
227 std::stringstream name;
228 name <<
"ICN " << dim <<
"D";
229 view =
new PView(name.str().c_str(),
"ElementData",
_m, dataPV);
237 int whichMeasure = computeICN ? 2 : computeIGE ? 1 : computeJac ? 0 : -1;
238 if(threshold < 99 && whichMeasure >= 0) {
241 #if defined(HAVE_OPENGL)
253 std::set<GEntity *, GEntityPtrFullLessThan> entities;
257 entities.insert(*it);
261 entities.insert(*it);
265 entities.insert(*it);
271 for(
auto it = entities.begin(); it != entities.end(); ++it) {
277 Msg::StatusBar(
true,
"Volume %d: checking the Jacobian of %d elements",
281 Msg::StatusBar(
true,
"Surface %d: checking the Jacobian of %d elements",
289 normals->
set(0, 0, 0);
290 normals->
set(0, 1, 0);
291 normals->
set(0, 2, 1);
296 Msg::StatusBar(
true,
"Line %d: checking the Jacobian of %d elements",
305 for(
unsigned i = 0; i < num; ++i) {
310 if(min < 0 && max < 0) ++cntInverted;
313 #if defined(HAVE_VISUDEV)
314 _computePointwiseQuantities(el, normals);
317 if(normals)
delete normals;
321 (cntInverted > 1) ?
"s are" :
" is");
333 for(std::size_t i = 0; i <
_data.size(); ++i) {
335 if(el->
getDim() != dim)
continue;
336 if(
_data[i].minJ() <= 0 &&
_data[i].maxJ() >= 0) {
_data[i].setMinS(0); }
352 for(std::size_t i = 0; i <
_data.size(); ++i) {
354 if(el->
getDim() != dim)
continue;
355 if(
_data[i].minJ() <= 0 &&
_data[i].maxJ() >= 0) {
_data[i].setMinI(0); }
371 for(std::size_t i = 0; i <
_data.size(); ++i) {
373 const int dim = el->
getDim();
374 if((askedDim == 4 && dim > 1) || dim == askedDim) {
376 switch(whichMeasure) {
377 case 2: q =
_data[i].minI();
break;
378 case 1: q =
_data[i].minS();
break;
380 if(
_data[i].maxJ() > 0)
382 else if(
_data[i].maxJ() < 0)
385 if((greater && q > threshold) || (!greater && q < threshold)) {
403 double infminJ, supminJ, avgminJ;
404 double infratJ, supratJ, avgratJ;
406 int count = 0, countc = 0;
408 infminJ = infratJ = 1e10;
409 supminJ = supratJ = -1e10;
410 avgminJ = avgratJ = avgratJc = 0;
412 for(std::size_t i = 0; i <
_data.size(); ++i) {
413 infminJ = std::min(infminJ,
_data[i].minJ());
414 supminJ = std::max(supminJ,
_data[i].minJ());
415 avgminJ +=
_data[i].minJ();
418 if(
_data[i].maxJ() > 0)
420 else if(
_data[i].maxJ() < 0)
422 infratJ = std::min(infratJ, q);
423 supratJ = std::max(supratJ, q);
426 avgratJc +=
_data[i].minJ() /
_data[i].maxJ();
436 Msg::Info(
"minJ = %8.3g, %8.3g, %8.3g (min, avg, max)", infminJ, avgminJ,
438 if(countc && countc < count)
439 Msg::Info(
"minJ/maxJ = %8.3g (avg on the %d "
440 "non-constant elements)",
442 Msg::Info(
"minJ/maxJ = %8.3g, %8.3g, %8.3g (worst, avg, best)", infratJ,
452 double infminS, supminS, avgminS;
453 infminS = supminS = avgminS =
_data[0].minS();
455 for(std::size_t i = 1; i <
_data.size(); ++i) {
456 infminS = std::min(infminS,
_data[i].minS());
457 supminS = std::max(supminS,
_data[i].minS());
458 avgminS +=
_data[i].minS();
460 avgminS /=
_data.size();
462 Msg::Info(
"IGE = %8.3g, %8.3g, %8.3g (worst, avg, best)", infminS,
472 double infminI, supminI, avgminI;
473 infminI = supminI = avgminI =
_data[0].minI();
475 for(std::size_t i = 1; i <
_data.size(); ++i) {
476 infminI = std::min(infminI,
_data[i].minI());
477 supminI = std::max(supminI,
_data[i].minI());
478 avgminI +=
_data[i].minI();
480 avgminI /=
_data.size();
482 Msg::Info(
"ICN = %8.3g, %8.3g, %8.3g (worst, avg, best)", infminI,
507 #if defined(HAVE_VISUDEV)
509 void GMSH_AnalyseMeshQualityPlugin::_computePointwiseQuantities(
512 if(_numElementToScan != -1 && el->
getNum() != _numElementToScan)
return;
532 std::vector<double> &vec = _dataPViewJac[num];
533 for(
int j = 0; j < tmpVector.
size(); ++j) vec.push_back(tmpVector(j));
538 std::vector<double> &vec = _dataPViewIGE[num];
539 for(
int j = 0; j < tmpVector.
size(); ++j) vec.push_back(tmpVector(j));
550 void GMSH_AnalyseMeshQualityPlugin::_setInterpolationMatrices(
PView *view)
553 for(
int type = 0; type < 20; ++type) {
554 if(!_type2tag[type])
continue;
567 void GMSH_AnalyseMeshQualityPlugin::_createPViewPointwise()
570 _setInterpolationMatrices(
new PView(
"Pointwise Jacobian",
"ElementNodeData",
571 _m, _dataPViewJac, 0, 1));
575 _setInterpolationMatrices(
576 new PView(
"Pointwise IGE",
"ElementNodeData",
_m, _dataPViewIGE, 0, 1));
580 _setInterpolationMatrices(
581 new PView(
"Pointwise ICN",
"ElementNodeData",
_m, _dataPViewICN, 0, 1));