gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
AnalyseMeshQuality.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #include "GmshConfig.h"
7 
8 #if defined(HAVE_MESH)
9 
10 #include "AnalyseMeshQuality.h"
12 #include "OS.h"
13 #include "Context.h"
14 #include "PView.h"
15 #include "GModel.h"
16 #include "MElement.h"
17 #include "bezierBasis.h"
18 #include <sstream>
19 #include <fstream>
20 #if defined(HAVE_OPENGL)
21 #include "drawContext.h"
22 #endif
23 #if defined(HAVE_VISUDEV)
24 #include "BasisFactory.h"
25 #endif
26 
27 StringXNumber CurvedMeshOptions_Number[] = {
28  {GMSH_FULLRC, "JacobianDeterminant", nullptr, 0},
29  {GMSH_FULLRC, "IGEMeasure", nullptr, 0},
30  {GMSH_FULLRC, "ICNMeasure", nullptr, 0},
31  {GMSH_FULLRC, "HidingThreshold", nullptr, 99},
32  {GMSH_FULLRC, "ThresholdGreater", nullptr, 1},
33  {GMSH_FULLRC, "CreateView", nullptr, 0},
34  {GMSH_FULLRC, "Recompute", nullptr, 0},
35  {GMSH_FULLRC, "DimensionOfElements", nullptr, -1}
36 #if defined(HAVE_VISUDEV)
37  ,
38  {GMSH_FULLRC, "Element to draw quality", nullptr, 0}
39 #endif
40 };
41 
42 extern "C" {
44 {
45  return new GMSH_AnalyseMeshQualityPlugin();
46 }
47 }
48 
50 {
51  return sizeof(CurvedMeshOptions_Number) / sizeof(StringXNumber);
52 }
53 
55 {
56  return &CurvedMeshOptions_Number[iopt];
57 }
58 
60 {
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"
69  "\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 "
78  "details.)\n"
79  "\n"
80  "Parameters:\n"
81  "\n"
82  "- `JacobianDeterminant': compute J?\n"
83  "\n"
84  "- `IGEMeasure': compute IGE?\n"
85  "\n"
86  "- `ICNMeasure': compute ICN?\n"
87  "\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"
93  "\n"
94  "- `CreateView': create a model-based view of min(J)/max(J), "
95  "min(IGE) and/or min(ICN)?\n"
96  "\n"
97  "- `Recompute': force recomputation (set to 1 if the mesh has "
98  "changed).\n"
99  "\n"
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.";
103 }
104 
106 {
107  _m = GModel::current();
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);
116 
117 #if defined(HAVE_VISUDEV)
118  _pwJac = computeJac / 2;
119  _pwIGE = computeIGE / 2;
120  _pwICN = computeICN / 2;
121 
122  _numElementToScan = static_cast<int>(CurvedMeshOptions_Number[8].def);
123  _viewOrder = 0;
124  _dataPViewJac.clear();
125  _dataPViewIGE.clear();
126  _dataPViewICN.clear();
127 #endif
128 
129  if(askedDim < 0 || askedDim > 4) askedDim = _m->getDim();
130 
131  if(recompute) _clear(askedDim);
132 
133  // Compute what have to
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) {
139  if(!_computedJac[dim - 1]) {
140  double time = Cpu(), w = TimeOfDay();
141  Msg::StatusBar(true, "Computing Jacobian for %dD elements...", dim);
143  Msg::StatusBar(true,
144  "Done computing Jacobian for %dD elements (Wall %gs, "
145  "CPU %gs)",
146  dim, TimeOfDay() - w, Cpu() - time);
147  printStatJ = true;
148  }
149  if(computeIGE && !_computedIGE[dim - 1]) {
150  double time = Cpu(), w = TimeOfDay();
151  Msg::StatusBar(true, "Computing IGE for %dD elements...", dim);
152  _computeMinIGE(dim);
153  Msg::StatusBar(true,
154  "Done computing IGE for %dD elements (Wall %gs, "
155  "CPU %gs)",
156  dim, TimeOfDay() - w, Cpu() - time);
157  printStatS = true;
158  }
159  if(computeICN && !_computedICN[dim - 1]) {
160  double time = Cpu(), w = TimeOfDay();
161  Msg::StatusBar(true, "Computing ICN for %dD elements...", dim);
162  _computeMinICN(dim);
163  Msg::StatusBar(true,
164  "Done computing ICN for %dD elements (Wall %gs, "
165  "CPU %gs)",
166  dim, TimeOfDay() - w, Cpu() - time);
167  printStatI = true;
168  }
169  }
170  }
171  if(printStatJ) _printStatJacobian();
172  if(printStatS) _printStatIGE();
173  if(printStatI) _printStatICN();
174 
175 #if defined(HAVE_VISUDEV)
176  _createPViewPointwise();
177 #endif
178 
179  // Create PView
180  PView *view = nullptr;
181  if(createView) {
182  for(int dim = 1; dim <= 3; ++dim) {
183  if((askedDim == 4 && dim > 1) || dim == askedDim) {
184  if(!_pviewJac[dim - 1] && computeJac) {
185  _pviewJac[dim - 1] = true;
186  std::map<int, std::vector<double> > dataPV;
187  for(std::size_t i = 0; i < _data.size(); ++i) {
188  MElement *const el = _data[i].element();
189  if(el->getDim() == dim) {
190  double q = 0;
191  if(_data[i].maxJ() > 0)
192  q = _data[i].minJ() / _data[i].maxJ();
193  else if(_data[i].maxJ() < 0)
194  q = _data[i].maxJ() / _data[i].minJ();
195  dataPV[el->getNum()].push_back(q);
196  }
197  }
198  if(dataPV.size()) {
199  std::stringstream name;
200  name << "minJ/maxJ " << dim << "D";
201  view = new PView(name.str().c_str(), "ElementData", _m, dataPV);
202  }
203  }
204  if(!_pviewIGE[dim - 1] && computeIGE) {
205  _pviewIGE[dim - 1] = true;
206  std::map<int, std::vector<double> > dataPV;
207  for(std::size_t i = 0; i < _data.size(); ++i) {
208  MElement *const el = _data[i].element();
209  if(el->getDim() == dim)
210  dataPV[el->getNum()].push_back(_data[i].minS());
211  }
212  if(dataPV.size()) {
213  std::stringstream name;
214  name << "IGE " << dim << "D";
215  view = new PView(name.str().c_str(), "ElementData", _m, dataPV);
216  }
217  }
218  if(!_pviewICN[dim - 1] && computeICN) {
219  _pviewICN[dim - 1] = true;
220  std::map<int, std::vector<double> > dataPV;
221  for(std::size_t i = 0; i < _data.size(); ++i) {
222  MElement *const el = _data[i].element();
223  if(el->getDim() == dim)
224  dataPV[el->getNum()].push_back(_data[i].minI());
225  }
226  if(dataPV.size()) {
227  std::stringstream name;
228  name << "ICN " << dim << "D";
229  view = new PView(name.str().c_str(), "ElementData", _m, dataPV);
230  }
231  }
232  }
233  }
234  }
235 
236  // hide elements
237  int whichMeasure = computeICN ? 2 : computeIGE ? 1 : computeJac ? 0 : -1;
238  if(threshold < 99 && whichMeasure >= 0) {
239  _hideWithThreshold(askedDim, whichMeasure, threshold, thresholdGreater);
241 #if defined(HAVE_OPENGL)
243 #endif
244  }
245 
246  return view;
247 }
248 
250 {
251  if(_computedJac[dim - 1]) return;
252 
253  std::set<GEntity *, GEntityPtrFullLessThan> entities;
254  switch(dim) {
255  case 3:
256  for(auto it = _m->firstRegion(); it != _m->lastRegion(); it++)
257  entities.insert(*it);
258  break;
259  case 2:
260  for(auto it = _m->firstFace(); it != _m->lastFace(); it++)
261  entities.insert(*it);
262  break;
263  case 1:
264  for(auto it = _m->firstEdge(); it != _m->lastEdge(); it++)
265  entities.insert(*it);
266  break;
267  default: return;
268  }
269 
270  int cntInverted = 0;
271  for(auto it = entities.begin(); it != entities.end(); ++it) {
272  GEntity *entity = *it;
273  unsigned num = entity->getNumMeshElements();
274  fullMatrix<double> *normals = nullptr;
275  switch(dim) {
276  case 3:
277  Msg::StatusBar(true, "Volume %d: checking the Jacobian of %d elements",
278  entity->tag(), num);
279  break;
280  case 2:
281  Msg::StatusBar(true, "Surface %d: checking the Jacobian of %d elements",
282  entity->tag(), num);
283  // check the classical case of 2D planar meshes in the z=0 plane to issue
284  // a warning of the mesh is oriented along -z
285  if(entity->geomType() == GEntity::DiscreteSurface) {
286  SBoundingBox3d bb = entity->bounds();
287  if(!bb.empty() && bb.max().z() - bb.min().z() == .0) {
288  normals = new fullMatrix<double>(1, 3);
289  normals->set(0, 0, 0);
290  normals->set(0, 1, 0);
291  normals->set(0, 2, 1);
292  }
293  }
294  break;
295  case 1:
296  Msg::StatusBar(true, "Line %d: checking the Jacobian of %d elements",
297  entity->tag(), num);
298  break;
299  default: break;
300  }
301 
302  MsgProgressStatus progress(num);
303 
304  _data.reserve(_data.size() + num);
305  for(unsigned i = 0; i < num; ++i) {
306  MElement *el = entity->getMeshElement(i);
307  double min, max;
309  _data.push_back(data_elementMinMax(el, min, max));
310  if(min < 0 && max < 0) ++cntInverted;
311  progress.next();
312 
313 #if defined(HAVE_VISUDEV)
314  _computePointwiseQuantities(el, normals);
315 #endif
316  }
317  if(normals) delete normals;
318  }
319  if(cntInverted) {
320  Msg::Warning("%d element%s completely inverted",
321  (cntInverted > 1) ? "s are" : " is");
322  }
323  _computedJac[dim - 1] = true;
325 }
326 
328 {
329  if(_computedIGE[dim - 1]) return;
330 
331  MsgProgressStatus progress(_data.size());
332 
333  for(std::size_t i = 0; i < _data.size(); ++i) {
334  MElement *const el = _data[i].element();
335  if(el->getDim() != dim) continue;
336  if(_data[i].minJ() <= 0 && _data[i].maxJ() >= 0) { _data[i].setMinS(0); }
337  else {
338  _data[i].setMinS(jacobianBasedQuality::minIGEMeasure(el, true));
339  }
340  progress.next();
341  }
342 
343  _computedIGE[dim - 1] = true;
344 }
345 
347 {
348  if(_computedICN[dim - 1]) return;
349 
350  MsgProgressStatus progress(_data.size());
351 
352  for(std::size_t i = 0; i < _data.size(); ++i) {
353  MElement *const el = _data[i].element();
354  if(el->getDim() != dim) continue;
355  if(_data[i].minJ() <= 0 && _data[i].maxJ() >= 0) { _data[i].setMinI(0); }
356  else {
357  _data[i].setMinI(jacobianBasedQuality::minICNMeasure(el, true));
358  }
359  progress.next();
360  }
361 
362  _computedICN[dim - 1] = true;
363 }
364 
366  int whichMeasure,
367  double threshold,
368  bool greater)
369 {
370  int nHidden = 0;
371  for(std::size_t i = 0; i < _data.size(); ++i) {
372  MElement *const el = _data[i].element();
373  const int dim = el->getDim();
374  if((askedDim == 4 && dim > 1) || dim == askedDim) {
375  double q = 1.;
376  switch(whichMeasure) {
377  case 2: q = _data[i].minI(); break;
378  case 1: q = _data[i].minS(); break;
379  case 0:
380  if(_data[i].maxJ() > 0)
381  q = _data[i].minJ() / _data[i].maxJ();
382  else if(_data[i].maxJ() < 0)
383  q = _data[i].maxJ() / _data[i].minJ();
384  }
385  if((greater && q > threshold) || (!greater && q < threshold)) {
386  el->setVisibility(0);
387  ++nHidden;
388  }
389  else {
390  el->setVisibility(1);
391  }
392  }
393  }
394  return nHidden;
395 }
396 
398 {
399  if(_data.empty()) {
400  Msg::Info("No stat to print");
401  return;
402  }
403  double infminJ, supminJ, avgminJ;
404  double infratJ, supratJ, avgratJ;
405  double avgratJc;
406  int count = 0, countc = 0;
407 
408  infminJ = infratJ = 1e10;
409  supminJ = supratJ = -1e10;
410  avgminJ = avgratJ = avgratJc = 0;
411 
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();
416 
417  double q = 0;
418  if(_data[i].maxJ() > 0)
419  q = _data[i].minJ() / _data[i].maxJ();
420  else if(_data[i].maxJ() < 0)
421  q = _data[i].maxJ() / _data[i].minJ();
422  infratJ = std::min(infratJ, q);
423  supratJ = std::max(supratJ, q);
424  avgratJ += q;
425  if(q < 1 - 1e-5) {
426  avgratJc += _data[i].minJ() / _data[i].maxJ();
427  ++countc;
428  }
429 
430  ++count;
431  }
432  avgminJ /= count;
433  avgratJ /= count;
434  avgratJc /= countc;
435 
436  Msg::Info("minJ = %8.3g, %8.3g, %8.3g (min, avg, max)", infminJ, avgminJ,
437  supminJ);
438  if(countc && countc < count)
439  Msg::Info("minJ/maxJ = %8.3g (avg on the %d "
440  "non-constant elements)",
441  avgratJc, countc);
442  Msg::Info("minJ/maxJ = %8.3g, %8.3g, %8.3g (worst, avg, best)", infratJ,
443  avgratJ, supratJ);
444 }
445 
447 {
448  if(_data.empty()) {
449  Msg::Info("No stat to print");
450  return;
451  }
452  double infminS, supminS, avgminS;
453  infminS = supminS = avgminS = _data[0].minS();
454 
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();
459  }
460  avgminS /= _data.size();
461 
462  Msg::Info("IGE = %8.3g, %8.3g, %8.3g (worst, avg, best)", infminS,
463  avgminS, supminS);
464 }
465 
467 {
468  if(_data.empty()) {
469  Msg::Info("No stat to print");
470  return;
471  }
472  double infminI, supminI, avgminI;
473  infminI = supminI = avgminI = _data[0].minI();
474 
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();
479  }
480  avgminI /= _data.size();
481 
482  Msg::Info("ICN = %8.3g, %8.3g, %8.3g (worst, avg, best)", infminI,
483  avgminI, supminI);
484 }
485 
486 void GMSH_AnalyseMeshQualityPlugin::_clear(int askedDim)
487 {
488  _data.clear();
489  if(askedDim < 4) {
490  _computedJac[askedDim - 1] = false;
491  _computedIGE[askedDim - 1] = false;
492  _computedICN[askedDim - 1] = false;
493  _pviewJac[askedDim - 1] = false;
494  _pviewIGE[askedDim - 1] = false;
495  _pviewICN[askedDim - 1] = false;
496  }
497  else {
498  _computedJac[1] = _computedJac[2] = false;
499  _computedIGE[1] = _computedIGE[2] = false;
500  _computedICN[1] = _computedICN[2] = false;
501  _pviewJac[1] = _pviewJac[2] = false;
502  _pviewIGE[1] = _pviewIGE[2] = false;
503  _pviewICN[1] = _pviewICN[2] = false;
504  }
505 }
506 
507 #if defined(HAVE_VISUDEV)
508 
509 void GMSH_AnalyseMeshQualityPlugin::_computePointwiseQuantities(
510  MElement *el, const fullMatrix<double> *normals)
511 {
512  if(_numElementToScan != -1 && el->getNum() != _numElementToScan) return;
513 
514  if(!_type2tag[el->getType()])
515  _type2tag[el->getType()] = el->getTypeForMSH();
516  else
517  // Skip if element tag is different from previous elements of same type
518  if(_type2tag[el->getType()] != el->getTypeForMSH())
519  return;
520 
521  static fullVector<double> tmpVector;
522  int num = el->getNum();
523 
524  if(!_viewOrder) {
525  // _viewOrder = std::min(10, 2 * el->getPolynomialOrder());
526  _viewOrder = 9;
527  }
528 
529  if(_pwJac) {
530  jacobianBasedQuality::sampleJacobianDeterminant(el, _viewOrder, tmpVector,
531  normals);
532  std::vector<double> &vec = _dataPViewJac[num];
533  for(int j = 0; j < tmpVector.size(); ++j) vec.push_back(tmpVector(j));
534  }
535 
536  if(_pwIGE) {
537  jacobianBasedQuality::sampleIGEMeasure(el, _viewOrder, tmpVector);
538  std::vector<double> &vec = _dataPViewIGE[num];
539  for(int j = 0; j < tmpVector.size(); ++j) vec.push_back(tmpVector(j));
540  }
541 
542  if(_pwICN) {
543  // jacobianBasedQuality::sampleICNMeasure(el, _viewOrder, tmpVector);
544  // std::vector<double> &vec = _dataPViewICN[num];
545  // for (int j = 0; j < tmpVector.size(); ++j)
546  // vec.push_back(tmpVector(j));
547  }
548 }
549 
550 void GMSH_AnalyseMeshQualityPlugin::_setInterpolationMatrices(PView *view)
551 {
552  PViewData *viewData = view->getData();
553  for(int type = 0; type < 20; ++type) {
554  if(!_type2tag[type]) continue;
555  viewData->deleteInterpolationMatrices(type);
556  const nodalBasis *fsE = BasisFactory::getNodalBasis(_type2tag[type]);
557  const polynomialBasis *pfsE = dynamic_cast<const polynomialBasis *>(fsE);
558  const int hoTag = ElementType::getType(type, _viewOrder);
559  const nodalBasis *fsV = BasisFactory::getNodalBasis(hoTag);
560  const polynomialBasis *pfsV = dynamic_cast<const polynomialBasis *>(fsV);
561  viewData->setInterpolationMatrices(type, pfsV->coefficients,
562  pfsV->monomials, pfsE->coefficients,
563  pfsE->monomials);
564  }
565 }
566 
567 void GMSH_AnalyseMeshQualityPlugin::_createPViewPointwise()
568 {
569  if(_pwJac) {
570  _setInterpolationMatrices(new PView("Pointwise Jacobian", "ElementNodeData",
571  _m, _dataPViewJac, 0, 1));
572  }
573 
574  if(_pwIGE) {
575  _setInterpolationMatrices(
576  new PView("Pointwise IGE", "ElementNodeData", _m, _dataPViewIGE, 0, 1));
577  }
578 
579  if(_pwICN) {
580  _setInterpolationMatrices(
581  new PView("Pointwise ICN", "ElementNodeData", _m, _dataPViewICN, 0, 1));
582  }
583 }
584 
585 #endif
586 
587 #endif
MElement::getNum
virtual std::size_t getNum() const
Definition: MElement.h:68
TimeOfDay
double TimeOfDay()
Definition: OS.cpp:399
GModel::firstEdge
eiter firstEdge()
Definition: GModel.h:356
PView
Definition: PView.h:27
GMSH_AnalyseMeshQualityPlugin
Definition: AnalyseMeshQuality.h:37
fullMatrix::set
void set(int r, int c, scalar v)
Definition: fullMatrix.h:289
GMSH_AnalyseMeshQualityPlugin::_printStatIGE
void _printStatIGE()
polynomialBasis::monomials
fullMatrix< double > monomials
Definition: polynomialBasis.h:20
fullVector< double >
MElement::getTypeForMSH
virtual int getTypeForMSH() const
Definition: MElement.h:488
GMSH_Plugin
Definition: Plugin.h:26
Msg::Info
static void Info(const char *fmt,...)
Definition: GmshMessage.cpp:587
StringXNumber::def
double def
Definition: Options.h:922
OS.h
MElement::getDim
virtual int getDim() const =0
GMSH_AnalyseMeshQualityPlugin::_data
std::vector< data_elementMinMax > _data
Definition: AnalyseMeshQuality.h:56
Msg::Warning
static void Warning(const char *fmt,...)
Definition: GmshMessage.cpp:543
Msg::StatusBar
static void StatusBar(bool log, const char *fmt,...)
Definition: GmshMessage.cpp:686
data_elementMinMax
Definition: AnalyseMeshQuality.h:17
StringXNumber
Definition: Options.h:918
SBoundingBox3d::min
SPoint3 min() const
Definition: SBoundingBox3d.h:90
jacobianBasedQuality::minMaxJacobianDeterminant
void minMaxJacobianDeterminant(MElement *el, double &min, double &max, const fullMatrix< double > *normals, bool debug)
Definition: qualityMeasuresJacobian.cpp:242
fullVector::size
int size() const
Definition: fullMatrix.h:69
GMSH_AnalyseMeshQualityPlugin::_computedJac
bool _computedJac[3]
Definition: AnalyseMeshQuality.h:53
PView.h
jacobianBasedQuality::sampleJacobianDeterminant
void sampleJacobianDeterminant(MElement *el, int deg, double &min, double &max, const fullMatrix< double > *normals)
Definition: qualityMeasuresJacobian.cpp:372
GMSH_AnalyseMeshQualityPlugin::execute
PView * execute(PView *)
GEntity
Definition: GEntity.h:31
GMSH_AnalyseMeshQualityPlugin::_hideWithThreshold
int _hideWithThreshold(int askedDim, int whichMeasure, double threshold, bool greater)
GMSH_RegisterAnalyseMeshQualityPlugin
GMSH_Plugin * GMSH_RegisterAnalyseMeshQualityPlugin()
drawContext::global
static drawContextGlobal * global()
Definition: drawContext.cpp:85
drawContextGlobal::draw
virtual void draw(bool rateLimited=true)
Definition: drawContext.h:97
jacobianBasedQuality::minICNMeasure
double minICNMeasure(MElement *el, bool knownValid, bool reversedOk, const fullMatrix< double > *normals, bool debug)
Definition: qualityMeasuresJacobian.cpp:326
CTX::instance
static CTX * instance()
Definition: Context.cpp:122
fullMatrix< double >
GModel::lastFace
fiter lastFace()
Definition: GModel.h:359
PViewData::setInterpolationMatrices
void setInterpolationMatrices(int type, const fullMatrix< double > &coefVal, const fullMatrix< double > &expVal)
Definition: PViewData.cpp:154
BasisFactory::getNodalBasis
static const nodalBasis * getNodalBasis(int tag)
Definition: BasisFactory.cpp:23
MElement::getType
virtual int getType() const =0
AnalyseMeshQuality.h
GMSH_AnalyseMeshQualityPlugin::_printStatJacobian
void _printStatJacobian()
GEntity::getNumMeshElements
virtual std::size_t getNumMeshElements() const
Definition: GEntity.h:348
GMSH_AnalyseMeshQualityPlugin::_printStatICN
void _printStatICN()
SBoundingBox3d::empty
bool empty()
Definition: SBoundingBox3d.h:36
GMSH_AnalyseMeshQualityPlugin::_computeMinIGE
void _computeMinIGE(int dim)
ElementType::getType
int getType(int parentType, int order, bool serendip=false)
Definition: ElementType.cpp:757
GMSH_FULLRC
#define GMSH_FULLRC
Definition: Options.h:20
GModel::getDim
int getDim() const
Definition: GModel.cpp:989
qualityMeasuresJacobian.h
MsgProgressStatus
Definition: GmshMessage.h:181
Cpu
double Cpu()
Definition: OS.cpp:366
PView::getData
PViewData * getData(bool useAdaptiveIfAvailable=false)
Definition: PView.cpp:233
GEntity::geomType
virtual GeomType geomType() const
Definition: GEntity.h:238
GMSH_AnalyseMeshQualityPlugin::_computedIGE
bool _computedIGE[3]
Definition: AnalyseMeshQuality.h:53
GModel::firstFace
fiter firstFace()
Definition: GModel.h:355
PViewData::deleteInterpolationMatrices
void deleteInterpolationMatrices(int type=0)
Definition: PViewData.cpp:194
GMSH_AnalyseMeshQualityPlugin::getNbOptions
int getNbOptions() const
polynomialBasis::coefficients
fullMatrix< double > coefficients
Definition: polynomialBasis.h:21
jacobianBasedQuality::sampleIGEMeasure
void sampleIGEMeasure(MElement *el, int deg, double &min, double &max)
Definition: qualityMeasuresJacobian.cpp:386
CTX::mesh
contextMeshOptions mesh
Definition: Context.h:313
MElement
Definition: MElement.h:30
jacobianBasedQuality::minIGEMeasure
double minIGEMeasure(MElement *el, bool knownValid, bool reversedOk, const fullMatrix< double > *normals, bool debug)
Definition: qualityMeasuresJacobian.cpp:280
GEntity::tag
int tag() const
Definition: GEntity.h:280
GMSH_AnalyseMeshQualityPlugin::_pviewICN
bool _pviewICN[3]
Definition: AnalyseMeshQuality.h:54
GMSH_AnalyseMeshQualityPlugin::_pviewIGE
bool _pviewIGE[3]
Definition: AnalyseMeshQuality.h:54
PViewData
Definition: PViewData.h:29
GEntity::DiscreteSurface
@ DiscreteSurface
Definition: GEntity.h:117
contextMeshOptions::changed
int changed
Definition: Context.h:82
GModel::firstRegion
riter firstRegion()
Definition: GModel.h:354
GModel::lastRegion
riter lastRegion()
Definition: GModel.h:358
GMSH_AnalyseMeshQualityPlugin::_pviewJac
bool _pviewJac[3]
Definition: AnalyseMeshQuality.h:54
ENT_ALL
#define ENT_ALL
Definition: GmshDefines.h:235
Context.h
nodalBasis
Definition: nodalBasis.h:12
GEntity::getMeshElement
virtual MElement * getMeshElement(std::size_t index) const
Definition: GEntity.h:363
GMSH_AnalyseMeshQualityPlugin::getHelp
std::string getHelp() const
MElement.h
bezierCoeff::releasePools
static void releasePools()
Definition: bezierBasis.cpp:830
GMSH_AnalyseMeshQualityPlugin::_computedICN
bool _computedICN[3]
Definition: AnalyseMeshQuality.h:53
GMSH_AnalyseMeshQualityPlugin::_computeMinMaxJandValidity
void _computeMinMaxJandValidity(int dim)
SPoint3::z
double z(void) const
Definition: SPoint3.h:129
polynomialBasis
Definition: polynomialBasis.h:16
MElement::setVisibility
virtual void setVisibility(char val)
Definition: MElement.h:97
GMSH_AnalyseMeshQualityPlugin::getOption
StringXNumber * getOption(int)
GModel::lastEdge
eiter lastEdge()
Definition: GModel.h:360
GModel.h
SBoundingBox3d::max
SPoint3 max() const
Definition: SBoundingBox3d.h:91
GEntity::bounds
virtual SBoundingBox3d bounds(bool fast=false)
Definition: GEntity.h:301
bezierBasis.h
SBoundingBox3d
Definition: SBoundingBox3d.h:21
GMSH_AnalyseMeshQualityPlugin::_m
GModel * _m
Definition: AnalyseMeshQuality.h:39
GModel::current
static GModel * current(int index=-1)
Definition: GModel.cpp:136
GMSH_AnalyseMeshQualityPlugin::_clear
void _clear(int askedDim)
GMSH_AnalyseMeshQualityPlugin::_computeMinICN
void _computeMinICN(int dim)
BasisFactory.h
drawContext.h