Newer
Older
// Gmsh - Copyright (C) 1997-2015 C. Geuzaine, J.-F. Remacle
// See the LICENSE.txt file for license information. Please report all
// bugs and problems to the public mailing list <gmsh@geuz.org>.
#include "GmshConfig.h"

Christophe Geuzaine
committed
#include "GmshMessage.h"
#include "Numeric.h"
#include "Context.h"
#include "MLine.h"
#include "MTriangle.h"
#include "MQuadrangle.h"
#include "MTetrahedron.h"
#include "MHexahedron.h"
#include "MPrism.h"
#include "MPyramid.h"
#include "BoundaryLayers.h"
#include "CenterlineField.h"
#include "GFaceCompound.h"
#include "Field.h"
#include "yamakawa.h"
Paul-Emile Bernard
committed
#include "pointInsertion.h"
#if defined(HAVE_OPTHOM)
#include "OptHomRun.h"
#include "OptHomElastic.h"
#endif

Christophe Geuzaine
committed
#include "PView.h"
#include "PViewData.h"

Christophe Geuzaine
committed
#endif

Christophe Geuzaine
committed
template<class T>
static void GetQualityMeasure(std::vector<T*> &ele,
double &gamma, double &gammaMin, double &gammaMax,
double &minSICN, double &minSICNMin, double &minSICNMax,
double quality[3][100])
for(unsigned int i = 0; i < ele.size(); i++){
double g = ele[i]->gammaShapeMeasure();
gamma += g;
gammaMin = std::min(gammaMin, g);
double s = ele[i]->minSICNShapeMeasure();
minSICN += s;
minSICNMin = std::min(minSICNMin, s);
minSICNMax = std::max(minSICNMax, s);
rho += r;
rhoMin = std::min(rhoMin, r);
rhoMax = std::max(rhoMax, r);
for(int j = 0; j < 100; j++){
if(s > (2*j-100) / 100. && s <= (2*j-98) / 100.) quality[0][j]++;
if(g > j / 100. && g <= (j + 1) / 100.) quality[1][j]++;
if(r > j / 100. && r <= (j + 1) / 100.) quality[2][j]++;
}
}
void GetStatistics(double stat[50], double quality[3][100])
for(int i = 0; i < 50; i++) stat[i] = 0.;
stat[0] = m->getNumVertices();
stat[1] = m->getNumEdges();
stat[2] = m->getNumFaces();
stat[3] = m->getNumRegions();
std::map<int, std::vector<GEntity*> > physicals[4];
stat[45] = physicals[0].size() + physicals[1].size() +
physicals[2].size() + physicals[3].size();
for(GModel::eiter it = m->firstEdge(); it != m->lastEdge(); ++it)
for(GModel::fiter it = m->firstFace(); it != m->lastFace(); ++it){
// TODO: make this an option! if((*it)->getVisibility()){
stat[5] += (*it)->mesh_vertices.size();
stat[7] += (*it)->triangles.size();
stat[8] += (*it)->quadrangles.size();
for(GModel::riter it = m->firstRegion(); it != m->lastRegion(); ++it){
stat[6] += (*it)->mesh_vertices.size();
stat[9] += (*it)->tetrahedra.size();
stat[10] += (*it)->hexahedra.size();
stat[11] += (*it)->prisms.size();
stat[12] += (*it)->pyramids.size();
}
stat[13] = CTX::instance()->meshTimer[0];
stat[14] = CTX::instance()->meshTimer[1];
stat[15] = CTX::instance()->meshTimer[2];
if(quality){
for(int i = 0; i < 3; i++)
for(int j = 0; j < 100; j++)
double minSICN = 0., minSICNMin = 1., minSICNMax = -1.;
double gamma = 0., gammaMin = 1., gammaMax = 0.;
double rho = 0., rhoMin = 1., rhoMax = 0.;
double N = stat[9] + stat[10] + stat[11] + stat[12];
if(N){ // if we have 3D elements
for(GModel::riter it = m->firstRegion(); it != m->lastRegion(); ++it){
GetQualityMeasure((*it)->tetrahedra, gamma, gammaMin, gammaMax,
minSICN, minSICNMin, minSICNMax, rho, rhoMin, rhoMax, quality);
GetQualityMeasure((*it)->hexahedra, gamma, gammaMin, gammaMax,
minSICN, minSICNMin, minSICNMax, rho, rhoMin, rhoMax, quality);
GetQualityMeasure((*it)->prisms, gamma, gammaMin, gammaMax,
minSICN, minSICNMin, minSICNMax, rho, rhoMin, rhoMax, quality);
GetQualityMeasure((*it)->pyramids, gamma, gammaMin, gammaMax,
minSICN, minSICNMin, minSICNMax, rho, rhoMin, rhoMax, quality);
else{ // 2D elements
N = stat[7] + stat[8];
for(GModel::fiter it = m->firstFace(); it != m->lastFace(); ++it){
GetQualityMeasure((*it)->quadrangles, gamma, gammaMin, gammaMax,
minSICN, minSICNMin, minSICNMax, rho, rhoMin, rhoMax, quality);
GetQualityMeasure((*it)->triangles, gamma, gammaMin, gammaMax,
minSICN, minSICNMin, minSICNMax, rho, rhoMin, rhoMax, quality);
stat[17] = minSICN / N; stat[18] = minSICNMin; stat[19] = minSICNMax;
stat[20] = gamma / N; stat[21] = gammaMin; stat[22] = gammaMax;
stat[23] = rho / N; stat[24] = rhoMin; stat[25] = rhoMax;
stat[26] = PView::list.size();
for(unsigned int i = 0; i < PView::list.size(); i++) {
PViewData *data = PView::list[i]->getData(true);
stat[27] += data->getNumPoints();
stat[28] += data->getNumLines();
stat[29] += data->getNumTriangles();
stat[30] += data->getNumQuadrangles();
stat[31] += data->getNumTetrahedra();
stat[32] += data->getNumHexahedra();
stat[33] += data->getNumPrisms();
stat[34] += data->getNumPyramids();

Christophe Geuzaine
committed
#endif
if(CTX::instance()->expertMode || !m->getNumVertices()) return false;
// try to detect obvious mistakes in characteristic lenghts (one of
// the most common cause for erroneous bug reports on the mailing
// list)
for(GModel::viter it = m->firstVertex(); it != m->lastVertex(); ++it)
sumAllLc += (*it)->prescribedMeshSizeAtVertex() * CTX::instance()->mesh.lcFactor;
sumAllLc /= (double)m->getNumVertices();
if(!sumAllLc || pow(CTX::instance()->lc / sumAllLc, dim) > 1.e10)

Christophe Geuzaine
committed
("Your choice of mesh element sizes will likely produce a very\n"

Christophe Geuzaine
committed
"large mesh. Do you really want to continue?\n\n"
"(To disable this warning in the future, select `Enable expert mode'\n"

Christophe Geuzaine
committed
return false;
}
static bool CancelDelaunayHybrid(GModel *m)
{
if(CTX::instance()->expertMode) return false;

Christophe Geuzaine
committed
int n = 0;
for(GModel::riter it = m->firstRegion(); it != m->lastRegion(); ++it)
n += (*it)->getNumMeshElements();
if(n)

Christophe Geuzaine
committed
("You are trying to generate a mixed structured/unstructured grid using\n"
"the 3D Delaunay algorithm. This algorithm cannot garantee that the\n"
"final mesh will be conforming. (You should probably use the 3D Frontal\n"
"algorithm instead.) Do you really want to continue with the Delaunay?\n\n"

Christophe Geuzaine
committed
"(To disable this warning in the future, select `Enable expert mode'\n"
}
static void Mesh0D(GModel *m)
{
for(GModel::viter it = m->firstVertex(); it != m->lastVertex(); ++it){
if(gv->mesh_vertices.empty())
gv->mesh_vertices.push_back(new MVertex(gv->x(), gv->y(), gv->z(), gv));
if(gv->points.empty())
gv->points.push_back(new MPoint(gv->mesh_vertices.back()));
}
for(GModel::viter it = m->firstVertex(); it != m->lastVertex(); ++it){
if (gv->meshMaster() != gv->tag()){
if (gv->correspondingVertices.empty()){
Paul-Emile Bernard
committed
GVertex *master = m->getVertexByTag(abs(gv->meshMaster()));
if(master)gv->correspondingVertices[gv->mesh_vertices[0]] = (master->mesh_vertices[0]);

Christophe Geuzaine
committed
Msg::StatusBar(true, "Meshing 1D...");
double t1 = Cpu();
for(GModel::eiter it = m->firstEdge(); it != m->lastEdge(); ++it)
(*it)->meshStatistics.status = GEdge::PENDING;
Msg::ResetProgressMeter();
int nIter = 0, nTot = m->getNumEdges();
for(GModel::eiter it = m->firstEdge(); it != m->lastEdge(); ++it){
if ((*it)->meshStatistics.status == GEdge::PENDING){
Paul-Emile Bernard
committed
mesher(*it);
nPending++;
if(!nIter) Msg::ProgressMeter(nPending, nTot, false, "Meshing 1D...");
if(nIter++ > 10) break;
}
double t2 = Cpu();
CTX::instance()->meshTimer[0] = t2 - t1;

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done meshing 1D (%g s)", CTX::instance()->meshTimer[0]);
}
if(CTX::instance()->createAppendMeshStatReport == 1)
statreport = Fopen(CTX::instance()->meshStatReportFileName.c_str(), "w");
else if(CTX::instance()->createAppendMeshStatReport == 2)
statreport = Fopen(CTX::instance()->meshStatReportFileName.c_str(), "a");
else
return;
if(statreport){
Msg::Error("Could not open file '%s'", CTX::instance()->meshStatReportFileName.c_str());
return;
}
double worst = 1, best = 0, avg = 0;
double e_long = 0, e_short = 1.e22, e_avg = 0;
int nTotT = 0, nTotE = 0, nTotGoodLength = 0, nTotGoodQuality = 0;
int nUnmeshed = 0, numFaces = 0;
if(CTX::instance()->createAppendMeshStatReport == 1){
fprintf(statreport, "2D stats\tname\t\t#faces\t\t#fail\t\t"

Christophe Geuzaine
committed
"#t\t\tQavg\t\tQbest\t\tQworst\t\t#Q>90\t\t#Q>90/#t\t"
"#e\t\ttau\t\t#Egood\t\t#Egood/#e\tCPU\n");
if(m->empty()){
fclose(statreport);
return;
}
for(GModel::fiter it = m->firstFace() ; it != m->lastFace(); ++it){
worst = std::min((*it)->meshStatistics.worst_element_shape, worst);
best = std::max((*it)->meshStatistics.best_element_shape, best);
avg += (*it)->meshStatistics.average_element_shape * (*it)->meshStatistics.nbTriangle;
e_avg += (*it)->meshStatistics.efficiency_index;
e_long = std::max((*it)->meshStatistics.longest_edge_length, e_long);
e_short = std::min((*it)->meshStatistics.smallest_edge_length, e_short);
if ((*it)->meshStatistics.status == GFace::FAILED ||
Paul-Emile Bernard
committed
(*it)->meshStatistics.status == GFace::PENDING) nUnmeshed++;
nTotT += (*it)->meshStatistics.nbTriangle;
nTotE += (*it)->meshStatistics.nbEdge;
nTotGoodLength += (*it)->meshStatistics.nbGoodLength;
nTotGoodQuality += (*it)->meshStatistics.nbGoodQuality;
numFaces++;
}
Msg::Info("*** Efficiency index for surface mesh tau=%g ", 100*exp(e_avg/(double)nTotE));
fprintf(statreport,"\t%16s\t%d\t\t%d\t\t", m->getName().c_str(), numFaces, nUnmeshed);
fprintf(statreport,"%d\t\t%8.7f\t%8.7f\t%8.7f\t%d\t\t%8.7f\t",

Christophe Geuzaine
committed
nTotT, avg / (double)nTotT, best, worst, nTotGoodQuality,
(double)nTotGoodQuality / nTotT);
fprintf(statreport,"%d\t\t%8.7f\t%d\t\t%8.7f\t%8.1f\n",

Christophe Geuzaine
committed
nTotE, exp(e_avg / (double)nTotE), nTotGoodLength,
(double)nTotGoodLength / nTotE, CTX::instance()->meshTimer[1]);

Christophe Geuzaine
committed
Msg::StatusBar(true, "Meshing 2D...");
for(GModel::fiter it = m->firstFace(); it != m->lastFace(); ++it)
(*it)->meshStatistics.status = GFace::PENDING;
// boundary layers are special: their generation (including vertices
// and curve meshes) is global as it depends on a smooth normal
// field generated from the surface mesh of the source surfaces

Christophe Geuzaine
committed
std::set<GFace*, GEntityLessThan> cf, f;

Christophe Geuzaine
committed
for(GModel::fiter it = m->firstFace(); it != m->lastFace(); ++it)
if ((*it)->geomType() == GEntity::CompoundSurface)
cf.insert(*it);
else
f.insert(*it);
Msg::ResetProgressMeter();
int nIter = 0, nTot = m->getNumFaces();

Christophe Geuzaine
committed
std::vector<GFace*> temp;
temp.insert(temp.begin(), f.begin(), f.end());
#pragma omp parallel for schedule (dynamic)

Christophe Geuzaine
committed
for(size_t K = 0 ; K < temp.size() ; K++){
Paul-Emile Bernard
committed
if (temp[K]->meshStatistics.status == GFace::PENDING){
Paul-Emile Bernard
committed
meshGFace mesher(true);
mesher(temp[K]);
Tristan Carrier Baudouin
committed
#if defined(HAVE_BFGS)
if(CTX::instance()->mesh.optimizeLloyd){
if (temp[K]->geomType()==GEntity::CompoundSurface ||

Christophe Geuzaine
committed
temp[K]->geomType()==GEntity::Plane ||
temp[K]->geomType()==GEntity::RuledSurface) {
if (temp[K]->meshAttributes.method != MESH_TRANSFINITE &&
!temp[K]->meshAttributes.extrude) {
smoothing smm(CTX::instance()->mesh.optimizeLloyd,6);
//m->writeMSH("beforeLLoyd.msh");
smm.optimize_face(temp[K]);
int rec = ((CTX::instance()->mesh.recombineAll ||

Christophe Geuzaine
committed
temp[K]->meshAttributes.recombine) &&
!CTX::instance()->mesh.recombine3DAll);
//m->writeMSH("afterLLoyd.msh");
if (rec) recombineIntoQuads(temp[K]);
//m->writeMSH("afterRecombine.msh");
}
}
}
Tristan Carrier Baudouin
committed
#endif
Paul-Emile Bernard
committed
{
nPending++;
}
}
if(!nIter) Msg::ProgressMeter(nPending, nTot, false, "Meshing 2D...");

Christophe Geuzaine
committed
}

Christophe Geuzaine
committed
for(std::set<GFace*, GEntityLessThan>::iterator it = cf.begin();
it != cf.end(); ++it){
meshGFace mesher(true);
Tristan Carrier Baudouin
committed
#if defined(HAVE_BFGS)
if(CTX::instance()->mesh.optimizeLloyd){
if ((*it)->geomType()==GEntity::CompoundSurface ||

Christophe Geuzaine
committed
(*it)->geomType()==GEntity::Plane ||
(*it)->geomType()==GEntity::RuledSurface) {
if ((*it)->meshAttributes.method != MESH_TRANSFINITE &&
!(*it)->meshAttributes.extrude) {
smoothing smm(CTX::instance()->mesh.optimizeLloyd,6);
//m->writeMSH("beforeLLoyd.msh");
smm.optimize_face(*it);
int rec = ((CTX::instance()->mesh.recombineAll ||

Christophe Geuzaine
committed
(*it)->meshAttributes.recombine) &&
!CTX::instance()->mesh.recombine3DAll);
//m->writeMSH("afterLLoyd.msh");
if (rec) recombineIntoQuads(*it);
//m->writeMSH("afterRecombine.msh");
}
}
}
Tristan Carrier Baudouin
committed
#endif
if(!nIter) Msg::ProgressMeter(nPending, nTot, false, "Meshing 2D...");

Christophe Geuzaine
committed
CTX::instance()->meshTimer[1] = t2 - t1;

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done meshing 2D (%g s)", CTX::instance()->meshTimer[1]);
}
static void FindConnectedRegions(std::vector<GRegion*> &delaunay,

Christophe Geuzaine
committed
std::vector<std::vector<GRegion*> > &connected)
const unsigned int nbVolumes = delaunay.size();
if (!nbVolumes)return;
while (delaunay.size()){
std::set<GRegion*> oneDomain;
GRegion *r = delaunay[0];
while(!_stack.empty()){
r = _stack.top();
_stack.pop();
oneDomain.insert(r);
std::list<GFace*> faces = r->faces();
for (std::list<GFace*> :: iterator it = faces.begin(); it != faces.end() ; ++it){
Paul-Emile Bernard
committed
GFace *gf = *it;
GRegion *other = gf->getRegion(0) == r ? gf->getRegion(1) : gf->getRegion(0);
if (other != 0 && oneDomain.find(other) == oneDomain.end())
_stack.push (other);
}
}
std::vector<GRegion*> temp1,temp2;
for (unsigned int i=0;i<delaunay.size();i++){
r = delaunay[i];
if (oneDomain.find(r) == oneDomain.end())temp1.push_back(r);
else temp2.push_back(r);
}
connected.push_back(temp2);
delaunay=temp1;
}
Msg::Info("Delaunay Meshing %d volumes with %d connected components",

Christophe Geuzaine
committed
nbVolumes,connected.size());

Christophe Geuzaine
committed
Msg::StatusBar(true, "Meshing 3D...");
double t1 = Cpu();
std::for_each(m->firstRegion(), m->lastRegion(), meshGRegionExtruded());
// then subdivide if necessary (unfortunately the subdivision is a
// global operation, which can require changing the surface mesh!)
// then mesh all the non-delaunay regions (front3D with netgen)
std::for_each(m->firstRegion(), m->lastRegion(), meshGRegion(delaunay));

Christophe Geuzaine
committed
// warn if attempting to use Delaunay for mixed meshes
if(delaunay.size() && CancelDelaunayHybrid(m)) return;
// and finally mesh the delaunay regions (again, this is global; but
// we mesh each connected part separately for performance and mesh
// quality reasons)
std::vector<std::vector<GRegion*> > connected;
FindConnectedRegions(delaunay, connected);
// remove quads elements for volumes that are recombined
for(unsigned int i = 0; i < connected.size(); i++){
if(CTX::instance()->mesh.recombine3DAll || gr->meshAttributes.recombine3D){
Paul-Emile Bernard
committed
std::list<GFace*> f = gr->faces();
for(std::list<GFace*>::iterator it = f.begin(); it != f.end() ; ++it)
quadsToTriangles(*it, 1000000);
double time_recombination = 0., vol_element_recombination = 0.;
double vol_hexa_recombination = 0.;
int nb_elements_recombination = 0, nb_hexa_recombination = 0;
// pragma OMP here ?
// additional code for experimental hex mesh
for(unsigned j = 0; j < connected[i].size(); j++){
Paul-Emile Bernard
committed
if (old_algo_hexa()){
Filler f;
f.treat_region(gr);
Paul-Emile Bernard
committed
}
else{
Filler3D f;
treat_region_ok = f.treat_region(gr);
}
if(treat_region_ok && (CTX::instance()->mesh.recombine3DAll ||
gr->meshAttributes.recombine3D)){
double a = Cpu();
Paul-Emile Bernard
committed
Recombinator rec;
rec.execute(gr);
Supplementary sup;
sup.execute(gr);
PostOp post;
post.execute(gr,0);
nb_elements_recombination += post.get_nb_elements();
nb_hexa_recombination += post.get_nb_hexahedra();
vol_element_recombination += post.get_vol_elements();
vol_hexa_recombination += post.get_vol_hexahedra();
Paul-Emile Bernard
committed
stringstream ss;
ss << "yamakawa_part_";
ss << gr->tag();
ss << ".msh";
export_gregion_mesh(gr, ss.str().c_str());
time_recombination += (Cpu() - a);
Paul-Emile Bernard
committed
if(CTX::instance()->mesh.recombine3DAll){
Msg::Info("RECOMBINATION timing:");
Msg::Info(" --- CUMULATIVE TIME RECOMBINATION : %g s.", time_recombination);
Msg::Info("RECOMBINATION CUMULATIVE STATISTICS:");
Msg::Info(".... Percentage of hexahedra (#) : %g",
nb_hexa_recombination*100./nb_elements_recombination);
Msg::Info(".... Percentage of hexahedra (Vol) : %g",
vol_hexa_recombination*100./vol_element_recombination);
Paul-Emile Bernard
committed
}
Thomas Toulorge
committed
// Ensure that all volume Jacobians are positive
m->setAllVolumesPositive();
double t2 = Cpu();
CTX::instance()->meshTimer[2] = t2 - t1;

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done meshing 3D (%g s)", CTX::instance()->meshTimer[2]);
}
void OptimizeMeshNetgen(GModel *m)
{

Christophe Geuzaine
committed
Msg::StatusBar(true, "Optimizing 3D mesh with Netgen...");
double t1 = Cpu();
std::for_each(m->firstRegion(), m->lastRegion(), optimizeMeshGRegionNetgen());
Thomas Toulorge
committed
// Ensure that all volume Jacobians are positive
m->setAllVolumesPositive();

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done optimizing 3D mesh with Netgen (%g s)", t2 - t1);

Christophe Geuzaine
committed
Msg::StatusBar(true, "Optimizing 3D mesh...");
std::for_each(m->firstRegion(), m->lastRegion(), optimizeMeshGRegionGmsh());
Thomas Toulorge
committed
// Ensure that all volume Jacobians are positive
m->setAllVolumesPositive();

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done optimizing 3D mesh (%g s)", t2 - t1);

Christophe Geuzaine
committed
Msg::StatusBar(true, "Adapting 3D mesh...");
for(int i = 0; i < 10; i++)
std::for_each(m->firstRegion(), m->lastRegion(), adaptMeshGRegion());

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done adaptating 3D mesh (%g s)", t2 - t1);
void RecombineMesh(GModel *m)
{

Christophe Geuzaine
committed
Msg::StatusBar(true, "Recombining 2D mesh...");
for(GModel::fiter it = m->firstFace(); it != m->lastFace(); ++it){
GFace *gf = *it;
recombineIntoQuads(gf);
}
CTX::instance()->mesh.changed = ENT_ALL;
double t2 = Cpu();

Christophe Geuzaine
committed
Msg::StatusBar(true, "Done recombining 2D mesh (%g s)", t2 - t1);
//#include <google/profiler.h>
if(CTX::instance()->lock) {
Msg::Info("I'm busy! Ask me that later...");
CTX::instance()->lock = 1;
Msg::ResetErrorCounter();
// Initialize pseudo random mesh generator with the same seed
srand(1);
std::for_each(m->firstRegion(), m->lastRegion(), deMeshGRegion());
std::for_each(m->firstFace(), m->lastFace(), deMeshGFace());
}
std::for_each(m->firstRegion(), m->lastRegion(), deMeshGRegion());
}
}
// Orient the line and surface meshes so that they match the orientation of
// the geometrical entities and/or the user orientation constraints
if(m->getMeshStatus() >= 1)
std::for_each(m->firstEdge(), m->lastEdge(), orientMeshGEdge());
if(m->getMeshStatus() >= 2)
std::for_each(m->firstFace(), m->lastFace(), orientMeshGFace());
// Optimize quality of 3D tet mesh
if(m->getMeshStatus() == 3){
for(int i = 0; i < std::max(CTX::instance()->mesh.optimize,

Christophe Geuzaine
committed
CTX::instance()->mesh.optimizeNetgen); i++){
if(CTX::instance()->mesh.optimize > i) OptimizeMesh(m);
if(CTX::instance()->mesh.optimizeNetgen > i) OptimizeMeshNetgen(m);
if(m->getMeshStatus() == 2 && CTX::instance()->mesh.algoSubdivide == 1)
RefineMesh(m, CTX::instance()->mesh.secondOrderLinear, true);
else if(m->getMeshStatus() == 3 && CTX::instance()->mesh.algoSubdivide == 2)
RefineMesh(m, CTX::instance()->mesh.secondOrderLinear, false, true);
// Compute homology if necessary
if(!Msg::GetErrorCount()) m->computeHomology();
if(m->getMeshStatus() && CTX::instance()->mesh.order > 1)
SetOrderN(m, CTX::instance()->mesh.order, CTX::instance()->mesh.secondOrderLinear,

Christophe Geuzaine
committed
CTX::instance()->mesh.secondOrderIncomplete);
// Optimize high order elements
if(CTX::instance()->mesh.hoOptimize){
#if defined(HAVE_OPTHOM)
if(CTX::instance()->mesh.hoOptimize < 0){
Thomas Toulorge
committed
ElasticAnalogy(GModel::current(), false);
OptHomParameters p;
p.nbLayers = CTX::instance()->mesh.hoNLayers;
p.BARRIER_MIN = CTX::instance()->mesh.hoThresholdMin;
p.BARRIER_MAX = CTX::instance()->mesh.hoThresholdMax;
p.dim = GModel::current()->getDim();
p.optPrimSurfMesh = CTX::instance()->mesh.hoOptPrimSurfMesh;
HighOrderMeshOptimizer(GModel::current(), p);
}
#else
Msg::Error("High-order mesh optimization requires the OPTHOM module");
#endif
}
Msg::Info("%d vertices %d elements",

Christophe Geuzaine
committed
m->getNumMeshVertices(), m->getNumMeshElements());
Msg::PrintErrorCounter("Mesh generation error summary");
CTX::instance()->lock = 0;
CTX::instance()->mesh.changed = ENT_ALL;