1 // Copyright (C) 2004-2022 CEA/DEN, EDF R&D
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 //=============================================================================
21 // File : GHS3DPlugin_GHS3D.cxx
23 // Author : Edward AGAPOV, modified by Lioka RAZAFINDRAZAKA (CEA) 09/02/2007
25 //=============================================================================
27 #include "GHS3DPlugin_GHS3D.hxx"
28 #include "GHS3DPlugin_Hypothesis.hxx"
29 #include "MG_Tetra_API.hxx"
31 #include <SMDS_FaceOfNodes.hxx>
32 #include <SMDS_LinearEdge.hxx>
33 #include <SMDS_MeshElement.hxx>
34 #include <SMDS_MeshNode.hxx>
35 #include <SMDS_VolumeOfNodes.hxx>
36 #include <SMESHDS_Group.hxx>
37 #include <SMESHDS_Mesh.hxx>
38 #include <SMESH_Comment.hxx>
39 #include <SMESH_File.hxx>
40 #include <SMESH_Group.hxx>
41 #include <SMESH_HypoFilter.hxx>
42 #include <SMESH_Mesh.hxx>
43 #include <SMESH_MeshAlgos.hxx>
44 #include <SMESH_MeshEditor.hxx>
45 #include <SMESH_MesherHelper.hxx>
46 #include <SMESH_OctreeNode.hxx>
47 #include <SMESH_subMeshEventListener.hxx>
48 #include <StdMeshers_QuadToTriaAdaptor.hxx>
49 #include <StdMeshers_ViscousLayers.hxx>
51 #include <BRepAdaptor_Surface.hxx>
52 #include <BRepBndLib.hxx>
53 #include <BRepBuilderAPI_MakeVertex.hxx>
54 #include <BRepClass3d.hxx>
55 #include <BRepClass3d_SolidClassifier.hxx>
56 #include <BRepExtrema_DistShapeShape.hxx>
57 #include <BRepGProp.hxx>
58 #include <BRepTools.hxx>
59 #include <BRep_Tool.hxx>
60 #include <Bnd_Box.hxx>
61 #include <GProp_GProps.hxx>
62 #include <GeomAPI_ProjectPointOnSurf.hxx>
63 #include <Precision.hxx>
64 #include <Standard_ErrorHandler.hxx>
65 #include <Standard_Failure.hxx>
66 #include <Standard_ProgramError.hxx>
68 #include <TopExp_Explorer.hxx>
69 #include <TopTools_IndexedMapOfShape.hxx>
70 #include <TopTools_ListIteratorOfListOfShape.hxx>
71 #include <TopTools_MapOfShape.hxx>
73 #include <TopoDS_Shell.hxx>
74 #include <TopoDS_Solid.hxx>
76 #include <Basics_Utils.hxx>
77 #include <utilities.h>
82 #include <boost/filesystem.hpp>
84 namespace boofs = boost::filesystem;
90 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
96 // flags returning state of enforced entities, returned from writeGMFFile
97 enum InvalidEnforcedFlags { FLAG_BAD_ENF_VERT = 1,
98 FLAG_BAD_ENF_NODE = 2,
99 FLAG_BAD_ENF_EDGE = 4,
100 FLAG_BAD_ENF_TRIA = 8
102 static std::string flagsToErrorStr( int anInvalidEnforcedFlags )
105 if ( anInvalidEnforcedFlags != 0 )
107 if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_VERT )
108 str = "There are enforced vertices incorrectly defined.\n";
109 if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_NODE )
110 str += "There are enforced nodes incorrectly defined.\n";
111 if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_EDGE )
112 str += "There are enforced edge incorrectly defined.\n";
113 if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_TRIA )
114 str += "There are enforced triangles incorrectly defined.\n";
119 // change results files permissions to user only (using boost to be used without C++17)
120 static void chmodUserOnly(const char* filename)
122 boofs::permissions(filename, boofs::remove_perms | boofs::group_all | boofs::others_all );
125 typedef const list<const SMDS_MeshFace*> TTriaList;
127 static const char theDomainGroupNamePrefix[] = "Domain_";
129 static void removeFile( const TCollection_AsciiString& fileName )
132 SMESH_File( fileName.ToCString() ).remove();
135 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
139 //=============================================================================
143 //=============================================================================
145 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, SMESH_Gen* gen)
146 : SMESH_3D_Algo(hypId, gen), _isLibUsed( false )
149 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
150 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
153 _compatibleHypothesis.push_back( GHS3DPlugin_Hypothesis::GetHypType());
154 _compatibleHypothesis.push_back( StdMeshers_ViscousLayers::GetHypType() );
155 _requireShape = false; // can work without shape
157 _computeCanceled = false;
158 _progressAdvance = 1e-4;
161 //=============================================================================
165 //=============================================================================
167 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
171 //=============================================================================
175 //=============================================================================
177 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
178 const TopoDS_Shape& aShape,
179 Hypothesis_Status& aStatus )
181 aStatus = SMESH_Hypothesis::HYP_OK;
184 _viscousLayersHyp = 0;
186 _removeLogOnSuccess = true;
187 _logInStandardOutput = false;
189 const list <const SMESHDS_Hypothesis * >& hyps =
190 GetUsedHypothesis(aMesh, aShape, /*ignoreAuxiliary=*/false);
191 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
192 for ( ; h != hyps.end(); ++h )
195 _hyp = dynamic_cast< const GHS3DPlugin_Hypothesis*> ( *h );
196 if ( !_viscousLayersHyp )
197 _viscousLayersHyp = dynamic_cast< const StdMeshers_ViscousLayers*> ( *h );
201 _keepFiles = _hyp->GetKeepFiles();
202 _removeLogOnSuccess = _hyp->GetRemoveLogOnSuccess();
203 _logInStandardOutput = _hyp->GetStandardOutputLog();
206 if ( _viscousLayersHyp )
207 error( _viscousLayersHyp->CheckHypothesis( aMesh, aShape, aStatus ));
209 return aStatus == HYP_OK;
213 //=======================================================================
214 //function : entryToShape
216 //=======================================================================
218 TopoDS_Shape GHS3DPlugin_GHS3D::entryToShape(std::string entry)
220 if ( SMESH_Gen_i::GetSMESHGen()->getStudyServant()->_is_nil() )
221 throw SALOME_Exception("MG-Tetra plugin can't work w/o publishing in the study");
223 GEOM::GEOM_Object_var aGeomObj;
224 TopoDS_Shape S = TopoDS_Shape();
225 SALOMEDS::SObject_var aSObj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->FindObjectID( entry.c_str() );
226 if (!aSObj->_is_nil() ) {
227 CORBA::Object_var obj = aSObj->GetObject();
228 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
231 if ( !aGeomObj->_is_nil() )
232 S = SMESH_Gen_i::GetSMESHGen()->GeomObjectToShape( aGeomObj.in() );
236 //================================================================================
238 * \brief returns id of a solid if a triangle defined by the nodes is a temporary face
239 * either on a side facet of pyramid or a top of pentahedron and defines sub-domian
240 * outside the volume; else returns HOLE_ID
242 //================================================================================
244 static int checkTmpFace(const SMDS_MeshNode* node1,
245 const SMDS_MeshNode* node2,
246 const SMDS_MeshNode* node3)
248 // find a pyramid sharing the 3 nodes
249 SMDS_ElemIteratorPtr vIt1 = node1->GetInverseElementIterator(SMDSAbs_Volume);
250 while ( vIt1->more() )
252 const SMDS_MeshElement* vol = vIt1->next();
253 const int nbNodes = vol->NbCornerNodes();
254 if ( nbNodes != 5 && nbNodes != 6 ) continue;
256 if ( (i2 = vol->GetNodeIndex( node2 )) >= 0 &&
257 (i3 = vol->GetNodeIndex( node3 )) >= 0 )
261 // Triangle defines sub-domian inside the pyramid if it's
262 // normal points out of the vol
264 // make i2 and i3 hold indices of base nodes of the vol while
265 // keeping the nodes order in the triangle
268 i2 = i3, i3 = vol->GetNodeIndex( node1 );
269 else if ( i3 == iApex )
270 i3 = i2, i2 = vol->GetNodeIndex( node1 );
272 int i3base = (i2+1) % 4; // next index after i2 within the pyramid base
273 bool isDomainInPyramid = ( i3base != i3 );
274 return isDomainInPyramid ? HOLE_ID : vol->getshapeId();
278 int i1 = vol->GetNodeIndex( node1 );
279 if (( i1 == 5 && i2 == 4 && i3 == 3 ) ||
280 ( i1 == 4 && i2 == 3 && i3 == 5 ) ||
281 ( i1 == 3 && i2 == 5 && i3 == 4 ))
284 return vol->getshapeId(); // triangle is a prism top
291 //=======================================================================
292 //function : findShapeID
293 //purpose : find the solid corresponding to MG-Tetra sub-domain following
294 // the technique proposed in MG-Tetra manual (available within
295 // MG-Tetra installation) in chapter "B.4 Subdomain (sub-region) assignment".
296 // In brief: normal of the triangle defined by the given nodes
297 // points out of the domain it is associated to
298 //=======================================================================
300 static int findShapeID(SMESH_Mesh& mesh,
301 const SMDS_MeshNode* node1,
302 const SMDS_MeshNode* node2,
303 const SMDS_MeshNode* node3,
304 const bool toMeshHoles)
306 const int invalidID = 0;
307 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
309 // face the nodes belong to
310 vector<const SMDS_MeshNode *> nodes(3);
314 const SMDS_MeshElement * face = meshDS->FindElement( nodes, SMDSAbs_Face, /*noMedium=*/true);
316 return checkTmpFace(node1, node2, node3);
318 std::cout << "bnd face " << face->GetID() << " - ";
320 // geom face the face assigned to
321 SMESH_MeshEditor editor(&mesh);
322 int geomFaceID = editor.FindShape( face );
324 return checkTmpFace(node1, node2, node3);
325 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
326 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
328 TopoDS_Face geomFace = TopoDS::Face( shape );
330 // solids bounded by geom face
331 TopTools_IndexedMapOfShape solids, shells;
332 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
333 for ( ; ansIt.More(); ansIt.Next() ) {
334 switch ( ansIt.Value().ShapeType() ) {
336 solids.Add( ansIt.Value() ); break;
338 shells.Add( ansIt.Value() ); break;
342 // analyse found solids
343 if ( solids.Extent() == 0 || shells.Extent() == 0)
346 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
347 if ( solids.Extent() == 1 )
350 return meshDS->ShapeToIndex( solid1 );
352 // - Are we at a hole boundary face?
353 if ( shells(1).IsSame( BRepClass3d::OuterShell( solid1 )) )
354 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touch another shell?
356 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
357 // check if any edge of shells(1) belongs to another shell
358 for ( ; eExp.More() && !touch; eExp.Next() ) {
359 ansIt = mesh.GetAncestors( eExp.Current() );
360 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
361 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
362 touch = ( !ansIt.Value().IsSame( shells(1) ));
366 return meshDS->ShapeToIndex( solid1 );
369 // find orientation of geom face within the first solid
370 TopExp_Explorer fExp( solid1, TopAbs_FACE );
371 for ( ; fExp.More(); fExp.Next() )
372 if ( geomFace.IsSame( fExp.Current() )) {
373 geomFace = TopoDS::Face( fExp.Current() );
377 return invalidID; // face not found
379 // normale to triangle
380 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
381 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
382 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
383 gp_Vec vec12( node1Pnt, node2Pnt );
384 gp_Vec vec13( node1Pnt, node3Pnt );
385 gp_Vec meshNormal = vec12 ^ vec13;
386 if ( meshNormal.SquareMagnitude() < DBL_MIN )
389 // get normale to geomFace at any node
390 bool geomNormalOK = false;
392 SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
393 for ( int i = 0; !geomNormalOK && i < 3; ++i )
395 // find UV of i-th node on geomFace
396 const SMDS_MeshNode* nNotOnSeamEdge = 0;
397 if ( helper.IsSeamShape( nodes[i]->getshapeId() )) {
398 if ( helper.IsSeamShape( nodes[(i+1)%3]->getshapeId() ))
399 nNotOnSeamEdge = nodes[(i+2)%3];
401 nNotOnSeamEdge = nodes[(i+1)%3];
404 gp_XY uv = helper.GetNodeUV( geomFace, nodes[i], nNotOnSeamEdge, &uvOK );
405 // check that uv is correct
408 TopoDS_Shape nodeShape = helper.GetSubShapeByNode( nodes[i], meshDS );
409 if ( !nodeShape.IsNull() )
410 switch ( nodeShape.ShapeType() )
412 case TopAbs_FACE: tol = BRep_Tool::Tolerance( TopoDS::Face( nodeShape )); break;
413 case TopAbs_EDGE: tol = BRep_Tool::Tolerance( TopoDS::Edge( nodeShape )); break;
414 case TopAbs_VERTEX: tol = BRep_Tool::Tolerance( TopoDS::Vertex( nodeShape )); break;
417 gp_Pnt nodePnt ( nodes[i]->X(), nodes[i]->Y(), nodes[i]->Z() );
418 BRepAdaptor_Surface surface( geomFace );
419 uvOK = ( nodePnt.Distance( surface.Value( uv.X(), uv.Y() )) < 2 * tol );
421 // normale to geomFace at UV
423 surface.D1( uv.X(), uv.Y(), nodePnt, du, dv );
424 geomNormal = du ^ dv;
425 if ( geomFace.Orientation() == TopAbs_REVERSED )
426 geomNormal.Reverse();
427 geomNormalOK = ( geomNormal.SquareMagnitude() > DBL_MIN * 1e3 );
435 bool isReverse = ( meshNormal * geomNormal ) < 0;
437 return meshDS->ShapeToIndex( solid1 );
439 if ( solids.Extent() == 1 )
440 return HOLE_ID; // we are inside a hole
442 return meshDS->ShapeToIndex( solids(2) );
445 //=======================================================================
446 //function : addElemInMeshGroup
447 //purpose : Update or create groups in mesh
448 //=======================================================================
450 static void addElemInMeshGroup(SMESH_Mesh* theMesh,
451 const SMDS_MeshElement* anElem,
452 std::string& groupName,
453 std::set<std::string>& /*groupsToRemove*/)
455 if ( !anElem ) return; // issue 0021776
457 bool groupDone = false;
458 SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
459 while (grIt->more()) {
460 SMESH_Group * group = grIt->next();
461 if ( !group ) continue;
462 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
463 if ( !groupDS ) continue;
464 if ( groupDS->GetType()==anElem->GetType() &&groupName.compare(group->GetName())==0) {
465 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
466 aGroupDS->SMDSGroup().Add(anElem);
474 SMESH_Group* aGroup = theMesh->AddGroup(anElem->GetType(), groupName.c_str());
475 aGroup->SetName( groupName.c_str() );
476 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
477 aGroupDS->SMDSGroup().Add(anElem);
481 throw SALOME_Exception(LOCALIZED("A given element was not added to a group"));
485 //=======================================================================
486 //function : updateMeshGroups
487 //purpose : Update or create groups in mesh
488 //=======================================================================
490 static void updateMeshGroups(SMESH_Mesh* theMesh, std::set<std::string> groupsToRemove)
492 SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
493 while (grIt->more()) {
494 SMESH_Group * group = grIt->next();
495 if ( !group ) continue;
496 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
497 if ( !groupDS ) continue;
498 std::string currentGroupName = (string)group->GetName();
499 if (groupDS->IsEmpty() && groupsToRemove.find(currentGroupName) != groupsToRemove.end()) {
500 // Previous group created by enforced elements
501 theMesh->RemoveGroup(groupDS->GetID());
506 //=======================================================================
507 //function : removeEmptyGroupsOfDomains
508 //purpose : remove empty groups named "Domain_nb" created due to
509 // "To make groups of domains" option.
510 //=======================================================================
512 static void removeEmptyGroupsOfDomains(SMESH_Mesh* mesh,
513 bool notEmptyAsWell = false)
515 const char* refName = theDomainGroupNamePrefix;
516 const size_t refLen = strlen( theDomainGroupNamePrefix );
518 std::list<int> groupIDs = mesh->GetGroupIds();
519 std::list<int>::const_iterator id = groupIDs.begin();
520 for ( ; id != groupIDs.end(); ++id )
522 SMESH_Group* group = mesh->GetGroup( *id );
523 if ( !group || ( !group->GetGroupDS()->IsEmpty() && !notEmptyAsWell ))
525 const char* name = group->GetName();
528 if ( strncmp( name, refName, refLen ) == 0 && // starts from refName;
529 isdigit( *( name + refLen )) && // refName is followed by a digit;
530 strtol( name + refLen, &end, 10) >= 0 && // there are only digits ...
531 *end == '\0') // ... till a string end.
533 mesh->RemoveGroup( *id );
538 //================================================================================
540 * \brief Create the groups corresponding to domains
542 //================================================================================
544 static void makeDomainGroups( std::vector< std::vector< const SMDS_MeshElement* > >& elemsOfDomain,
545 SMESH_MesherHelper* theHelper)
547 // int nbDomains = 0;
548 // for ( size_t i = 0; i < elemsOfDomain.size(); ++i )
549 // nbDomains += ( elemsOfDomain[i].size() > 0 );
551 // if ( nbDomains > 1 )
552 for ( size_t iDomain = 0; iDomain < elemsOfDomain.size(); ++iDomain )
554 std::vector< const SMDS_MeshElement* > & elems = elemsOfDomain[ iDomain ];
555 if ( elems.empty() ) continue;
557 // find existing groups
558 std::vector< SMESH_Group* > groupOfType( SMDSAbs_NbElementTypes, (SMESH_Group*)NULL );
559 const std::string domainName = ( SMESH_Comment( theDomainGroupNamePrefix ) << iDomain );
560 SMESH_Mesh::GroupIteratorPtr groupIt = theHelper->GetMesh()->GetGroups();
561 while ( groupIt->more() )
563 SMESH_Group* group = groupIt->next();
564 if ( domainName == group->GetName() &&
565 dynamic_cast< SMESHDS_Group* >( group->GetGroupDS()) )
566 groupOfType[ group->GetGroupDS()->GetType() ] = group;
568 // create and fill the groups
572 SMESH_Group* group = groupOfType[ elems[ iElem ]->GetType() ];
574 group = theHelper->GetMesh()->AddGroup( elems[ iElem ]->GetType(),
575 domainName.c_str() );
576 SMDS_MeshGroup& groupDS =
577 static_cast< SMESHDS_Group* >( group->GetGroupDS() )->SMDSGroup();
579 while ( iElem < elems.size() && groupDS.Add( elems[iElem] ))
582 } while ( iElem < elems.size() );
586 //=======================================================================
587 //function : readGMFFile
588 //purpose : read GMF file w/o geometry associated to mesh
589 //=======================================================================
591 static bool readGMFFile(MG_Tetra_API* MGOutput,
593 GHS3DPlugin_GHS3D* theAlgo,
594 SMESH_MesherHelper* theHelper,
595 std::vector <const SMDS_MeshNode*> & theNodeByGhs3dId,
596 std::vector <const SMDS_MeshElement*> & theFaceByGhs3dId,
597 map<const SMDS_MeshNode*,int> & /*theNodeToGhs3dIdMap*/,
598 std::vector<std::string> & aNodeGroupByGhs3dId,
599 std::vector<std::string> & anEdgeGroupByGhs3dId,
600 std::vector<std::string> & aFaceGroupByGhs3dId,
601 std::set<std::string> & groupsToRemove,
602 bool toMakeGroupsOfDomains=false,
603 bool toMeshHoles=true)
606 SMESHDS_Mesh* theMeshDS = theHelper->GetMeshDS();
607 const bool hasGeom = ( theHelper->GetMesh()->HasShapeToMesh() );
609 int nbInitialNodes = (int) theNodeByGhs3dId.size();
612 const bool isQuadMesh =
613 theHelper->GetMesh()->NbEdges( ORDER_QUADRATIC ) ||
614 theHelper->GetMesh()->NbFaces( ORDER_QUADRATIC ) ||
615 theHelper->GetMesh()->NbVolumes( ORDER_QUADRATIC );
616 std::cout << "theNodeByGhs3dId.size(): " << nbInitialNodes << std::endl;
617 std::cout << "theHelper->GetMesh()->NbNodes(): " << theMeshDS->NbNodes() << std::endl;
618 std::cout << "isQuadMesh: " << isQuadMesh << std::endl;
621 // ---------------------------------
622 // Read generated elements and nodes
623 // ---------------------------------
625 int nbElem = 0, nbRef = 0;
627 std::vector< const SMDS_MeshNode*> GMFNode;
629 std::map<int, std::set<int> > subdomainId2tetraId;
631 std::map <GmfKwdCod,int> tabRef;
632 const bool force3d = !hasGeom;
635 tabRef[GmfVertices] = 3; // for new nodes and enforced nodes
636 tabRef[GmfCorners] = 1;
637 tabRef[GmfEdges] = 2; // for enforced edges
638 tabRef[GmfRidges] = 1;
639 tabRef[GmfTriangles] = 3; // for enforced faces
640 tabRef[GmfQuadrilaterals] = 4;
641 tabRef[GmfTetrahedra] = 4; // for new tetras
642 tabRef[GmfHexahedra] = 8;
645 int InpMsh = MGOutput->GmfOpenMesh( theFile, GmfRead, &ver, &dim);
649 // Read ids of domains
650 vector< int > solidIDByDomain;
653 int solid1; // id used in case of 1 domain or some reading failure
654 if ( theHelper->GetSubShape().ShapeType() == TopAbs_SOLID )
655 solid1 = theHelper->GetSubShapeID();
657 solid1 = theMeshDS->ShapeToIndex
658 ( TopExp_Explorer( theHelper->GetSubShape(), TopAbs_SOLID ).Current() );
660 int nbDomains = MGOutput->GmfStatKwd( InpMsh, GmfSubDomainFromGeom );
663 solidIDByDomain.resize( nbDomains+1, theHelper->GetSubShapeID() );
664 int faceNbNodes, faceIndex, orientation, domainNb;
665 MGOutput->GmfGotoKwd( InpMsh, GmfSubDomainFromGeom );
666 for ( int i = 0; i < nbDomains; ++i )
669 MGOutput->GmfGetLin( InpMsh, GmfSubDomainFromGeom,
670 &faceNbNodes, &faceIndex, &orientation, &domainNb, i);
671 solidIDByDomain[ domainNb ] = 1;
672 if ( 0 < faceIndex && faceIndex-1 < (int)theFaceByGhs3dId.size() )
674 const SMDS_MeshElement* face = theFaceByGhs3dId[ faceIndex-1 ];
675 const SMDS_MeshNode* nn[3] = { face->GetNode(0),
678 if ( orientation < 0 )
679 std::swap( nn[1], nn[2] );
680 solidIDByDomain[ domainNb ] =
681 findShapeID( *theHelper->GetMesh(), nn[0], nn[1], nn[2], toMeshHoles );
682 if ( solidIDByDomain[ domainNb ] > 0 )
685 std::cout << "solid " << solidIDByDomain[ domainNb ] << std::endl;
687 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( solidIDByDomain[ domainNb ] );
688 if ( ! theHelper->IsSubShape( foundShape, theHelper->GetSubShape() ))
689 solidIDByDomain[ domainNb ] = HOLE_ID;
694 if ( solidIDByDomain.size() < 2 )
695 solidIDByDomain.resize( 2, solid1 );
698 // Issue 0020682. Avoid creating nodes and tetras at place where
699 // volumic elements already exist
700 SMESH_ElementSearcher* elemSearcher = 0;
701 std::vector< const SMDS_MeshElement* > foundVolumes;
702 if ( !hasGeom && theHelper->GetMesh()->NbVolumes() > 0 )
703 elemSearcher = SMESH_MeshAlgos::GetElementSearcher( *theMeshDS );
704 unique_ptr< SMESH_ElementSearcher > elemSearcherDeleter( elemSearcher );
706 // IMP 0022172: [CEA 790] create the groups corresponding to domains
707 std::vector< std::vector< const SMDS_MeshElement* > > elemsOfDomain;
709 int nbVertices = MGOutput->GmfStatKwd( InpMsh, GmfVertices ) - nbInitialNodes;
710 if ( nbVertices < 0 )
712 GMFNode.resize( nbVertices + 1 );
714 std::map <GmfKwdCod,int>::const_iterator it = tabRef.begin();
715 for ( ; it != tabRef.end() ; ++it)
717 if(theAlgo->computeCanceled()) {
721 GmfKwdCod token = it->first;
724 nbElem = MGOutput->GmfStatKwd( InpMsh, token);
726 MGOutput->GmfGotoKwd( InpMsh, token);
727 std::cout << "Read " << nbElem;
732 std::vector<int> id (nbElem*tabRef[token]); // node ids
733 std::vector<int> domainID( nbElem ); // domain
735 if (token == GmfVertices) {
736 (nbElem <= 1) ? tmpStr = " vertex" : tmpStr = " vertices";
737 // std::cout << nbInitialNodes << " from input mesh " << std::endl;
739 // Remove orphan nodes from previous enforced mesh which was cleared
740 // if ( nbElem < nbMeshNodes ) {
741 // const SMDS_MeshNode* node;
742 // SMDS_NodeIteratorPtr nodeIt = theMeshDS->nodesIterator();
743 // while ( nodeIt->more() )
745 // node = nodeIt->next();
746 // if (theNodeToGhs3dIdMap.find(node) != theNodeToGhs3dIdMap.end())
747 // theMeshDS->RemoveNode(node);
756 const SMDS_MeshNode * aGMFNode;
758 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
759 if(theAlgo->computeCanceled()) {
762 if (ver == GmfFloat) {
763 MGOutput->GmfGetLin( InpMsh, token, &VerTab_f[0], &VerTab_f[1], &VerTab_f[2], &dummy);
769 MGOutput->GmfGetLin( InpMsh, token, &x, &y, &z, &dummy);
771 if (iElem >= nbInitialNodes) {
773 elemSearcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_Volume, foundVolumes))
776 aGMFNode = theHelper->AddNode(x, y, z);
778 aGMFID = iElem -nbInitialNodes +1;
779 GMFNode[ aGMFID ] = aGMFNode;
780 if (aGMFID-1 < (int)aNodeGroupByGhs3dId.size() && !aNodeGroupByGhs3dId.at(aGMFID-1).empty())
781 addElemInMeshGroup(theHelper->GetMesh(), aGMFNode, aNodeGroupByGhs3dId.at(aGMFID-1), groupsToRemove);
785 else if (token == GmfCorners && nbElem > 0) {
786 (nbElem <= 1) ? tmpStr = " corner" : tmpStr = " corners";
787 for ( int iElem = 0; iElem < nbElem; iElem++ )
788 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]]);
790 else if (token == GmfRidges && nbElem > 0) {
791 (nbElem <= 1) ? tmpStr = " ridge" : tmpStr = " ridges";
792 for ( int iElem = 0; iElem < nbElem; iElem++ )
793 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]]);
795 else if (token == GmfEdges && nbElem > 0) {
796 (nbElem <= 1) ? tmpStr = " edge" : tmpStr = " edges";
797 for ( int iElem = 0; iElem < nbElem; iElem++ )
798 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &domainID[iElem]);
800 else if (token == GmfTriangles && nbElem > 0) {
801 (nbElem <= 1) ? tmpStr = " triangle" : tmpStr = " triangles";
802 for ( int iElem = 0; iElem < nbElem; iElem++ )
803 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &domainID[iElem]);
805 else if (token == GmfQuadrilaterals && nbElem > 0) {
806 (nbElem <= 1) ? tmpStr = " Quadrilateral" : tmpStr = " Quadrilaterals";
807 for ( int iElem = 0; iElem < nbElem; iElem++ )
808 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3], &domainID[iElem]);
810 else if (token == GmfTetrahedra && nbElem > 0) {
811 (nbElem <= 1) ? tmpStr = " Tetrahedron" : tmpStr = " Tetrahedra";
812 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
813 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3], &domainID[iElem]);
815 subdomainId2tetraId[dummy].insert(iElem+1);
819 else if (token == GmfHexahedra && nbElem > 0) {
820 (nbElem <= 1) ? tmpStr = " Hexahedron" : tmpStr = " Hexahedra";
821 for ( int iElem = 0; iElem < nbElem; iElem++ )
822 MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
823 &id[iElem*tabRef[token]+4], &id[iElem*tabRef[token]+5], &id[iElem*tabRef[token]+6], &id[iElem*tabRef[token]+7], &domainID[iElem]);
825 std::cout << tmpStr << std::endl;
826 std::cout << std::endl;
833 case GmfQuadrilaterals:
837 std::vector< const SMDS_MeshNode* > node( nbRef );
838 std::vector< int > nodeID( nbRef );
839 std::vector< SMDS_MeshNode* > enfNode( nbRef );
840 const SMDS_MeshElement* aCreatedElem;
842 for ( int iElem = 0; iElem < nbElem; iElem++ )
844 if(theAlgo->computeCanceled()) {
847 // Check if elem is already in input mesh. If yes => skip
848 bool fullyCreatedElement = false; // if at least one of the nodes was created
849 for ( int iRef = 0; iRef < nbRef; iRef++ )
851 aGMFNodeID = id[iElem*tabRef[token]+iRef]; // read nbRef aGMFNodeID
852 if (aGMFNodeID <= nbInitialNodes) // input nodes
855 node[ iRef ] = theNodeByGhs3dId[aGMFNodeID];
859 fullyCreatedElement = true;
860 aGMFNodeID -= nbInitialNodes;
861 nodeID[ iRef ] = aGMFNodeID ;
862 node [ iRef ] = GMFNode[ aGMFNodeID ];
869 if (fullyCreatedElement) {
870 aCreatedElem = theHelper->AddEdge( node[0], node[1], noID, force3d );
871 if (anEdgeGroupByGhs3dId.size() && !anEdgeGroupByGhs3dId[iElem].empty())
872 addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, anEdgeGroupByGhs3dId[iElem], groupsToRemove);
876 if (fullyCreatedElement) {
877 aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], noID, force3d );
878 if (aFaceGroupByGhs3dId.size() && !aFaceGroupByGhs3dId[iElem].empty())
879 addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, aFaceGroupByGhs3dId[iElem], groupsToRemove);
882 case GmfQuadrilaterals:
883 if (fullyCreatedElement) {
884 aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], node[3], noID, force3d );
890 solidID = solidIDByDomain[ domainID[iElem]];
891 if ( solidID != HOLE_ID )
893 aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
895 theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
896 for ( int iN = 0; iN < 4; ++iN )
897 if ( node[iN]->getshapeId() < 1 )
898 theMeshDS->SetNodeInVolume( node[iN], solidID );
903 if ( elemSearcher ) {
904 // Issue 0020682. Avoid creating nodes and tetras at place where
905 // volumic elements already exist
906 if ( !node[1] || !node[0] || !node[2] || !node[3] )
908 if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
909 SMESH_TNodeXYZ(node[1]) +
910 SMESH_TNodeXYZ(node[2]) +
911 SMESH_TNodeXYZ(node[3]) ) / 4.,
912 SMDSAbs_Volume, foundVolumes ))
915 aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
922 solidID = solidIDByDomain[ domainID[iElem]];
923 if ( solidID != HOLE_ID )
925 aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
926 node[4], node[7], node[6], node[5],
928 theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
929 for ( int iN = 0; iN < 8; ++iN )
930 if ( node[iN]->getshapeId() < 1 )
931 theMeshDS->SetNodeInVolume( node[iN], solidID );
936 if ( elemSearcher ) {
937 // Issue 0020682. Avoid creating nodes and tetras at place where
938 // volumic elements already exist
939 if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] || !node[6] || !node[7])
941 if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
942 SMESH_TNodeXYZ(node[1]) +
943 SMESH_TNodeXYZ(node[2]) +
944 SMESH_TNodeXYZ(node[3]) +
945 SMESH_TNodeXYZ(node[4]) +
946 SMESH_TNodeXYZ(node[5]) +
947 SMESH_TNodeXYZ(node[6]) +
948 SMESH_TNodeXYZ(node[7])) / 8.,
949 SMDSAbs_Volume, foundVolumes ))
952 aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
953 node[4], node[7], node[6], node[5],
960 // care about medium nodes
962 aCreatedElem->IsQuadratic() &&
963 ( solidID = aCreatedElem->getshapeId() ) > 0 )
965 int iN = aCreatedElem->NbCornerNodes(), nbN = aCreatedElem->NbNodes();
966 for ( ; iN < nbN; ++iN )
968 const SMDS_MeshNode* n = aCreatedElem->GetNode(iN);
969 if ( n->getshapeId() < 1 )
970 theMeshDS->SetNodeInVolume( n, solidID );
974 if ( aCreatedElem && toMakeGroupsOfDomains )
976 if ( domainID[iElem] >= (int) elemsOfDomain.size() )
977 elemsOfDomain.resize( domainID[iElem] + 1 );
978 elemsOfDomain[ domainID[iElem] ].push_back( aCreatedElem );
980 } // loop on elements of one type
987 // remove nodes in holes
990 for ( int i = 1; i <= nbVertices; ++i )
991 if ( GMFNode[i]->NbInverseElements() == 0 )
992 theMeshDS->RemoveFreeNode( GMFNode[i], /*sm=*/0, /*fromGroups=*/false );
995 MGOutput->GmfCloseMesh( InpMsh);
997 // 0022172: [CEA 790] create the groups corresponding to domains
998 if ( toMakeGroupsOfDomains )
999 makeDomainGroups( elemsOfDomain, theHelper );
1002 std::map<int, std::set<int> >::const_iterator subdomainIt = subdomainId2tetraId.begin();
1003 TCollection_AsciiString aSubdomainFileName = theFile;
1004 aSubdomainFileName = aSubdomainFileName + ".subdomain";
1005 ofstream aSubdomainFile ( aSubdomainFileName.ToCString() , ios::out);
1007 aSubdomainFile << "Nb subdomains " << subdomainId2tetraId.size() << std::endl;
1008 for(;subdomainIt != subdomainId2tetraId.end() ; ++subdomainIt) {
1009 int subdomainId = subdomainIt->first;
1010 std::set<int> tetraIds = subdomainIt->second;
1011 std::set<int>::const_iterator tetraIdsIt = tetraIds.begin();
1012 aSubdomainFile << subdomainId << std::endl;
1013 for(;tetraIdsIt != tetraIds.end() ; ++tetraIdsIt) {
1014 aSubdomainFile << (*tetraIdsIt) << " ";
1016 aSubdomainFile << std::endl;
1018 aSubdomainFile.close();
1025 static bool writeGMFFile(MG_Tetra_API* MGInput,
1026 const char* theMeshFileName,
1027 const char* theRequiredFileName,
1028 const char* theSolFileName,
1029 const SMESH_ProxyMesh& theProxyMesh,
1030 SMESH_MesherHelper& theHelper,
1031 std::vector <const SMDS_MeshNode*> & theNodeByGhs3dId,
1032 std::vector <const SMDS_MeshElement*> & theFaceByGhs3dId,
1033 std::map<const SMDS_MeshNode*,int> & aNodeToGhs3dIdMap,
1034 std::vector<std::string> & aNodeGroupByGhs3dId,
1035 std::vector<std::string> & anEdgeGroupByGhs3dId,
1036 std::vector<std::string> & aFaceGroupByGhs3dId,
1037 GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap & theEnforcedNodes,
1038 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedEdges,
1039 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedTriangles,
1040 std::map<std::vector<double>, std::string> & enfVerticesWithGroup,
1041 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues & theEnforcedVertices,
1042 int & theInvalidEnforcedFlags)
1045 int idx, idxRequired = 0, idxSol = 0;
1046 const int dummyint = 0;
1047 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues::const_iterator vertexIt;
1048 std::vector<double> enfVertexSizes;
1049 const SMDS_MeshElement* elem;
1050 TIDSortedElemSet anElemSet, theKeptEnforcedEdges, theKeptEnforcedTriangles;
1051 SMDS_ElemIteratorPtr nodeIt;
1052 std::vector <const SMDS_MeshNode*> theEnforcedNodeByGhs3dId;
1053 map<const SMDS_MeshNode*,int> anEnforcedNodeToGhs3dIdMap, anExistingEnforcedNodeToGhs3dIdMap;
1054 std::vector< const SMDS_MeshElement* > foundElems;
1055 map<const SMDS_MeshNode*,TopAbs_State> aNodeToTopAbs_StateMap;
1057 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap::iterator elemIt;
1058 TIDSortedElemSet::iterator elemSetIt;
1060 SMESH_Mesh* theMesh = theHelper.GetMesh();
1061 const bool hasGeom = theMesh->HasShapeToMesh();
1062 SMESHUtils::Deleter< SMESH_ElementSearcher > pntCls
1063 ( SMESH_MeshAlgos::GetElementSearcher(*theMesh->GetMeshDS()));
1065 int nbEnforcedVertices = (int) theEnforcedVertices.size();
1066 theInvalidEnforcedFlags = 0;
1069 smIdType nbFaces = theProxyMesh.NbFaces();
1071 theFaceByGhs3dId.reserve( nbFaces );
1073 // groups management
1074 int usedEnforcedNodes = 0;
1075 std::string gn = "";
1080 idx = MGInput->GmfOpenMesh( theMeshFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1084 /* ========================== FACES ========================== */
1085 /* TRIANGLES ========================== */
1086 SMDS_ElemIteratorPtr eIt =
1087 hasGeom ? theProxyMesh.GetFaces( theHelper.GetSubShape()) : theProxyMesh.GetFaces();
1088 while ( eIt->more() )
1091 anElemSet.insert(elem);
1092 nodeIt = elem->nodesIterator();
1093 nbNodes = elem->NbCornerNodes();
1094 while ( nodeIt->more() && nbNodes--)
1097 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1098 int newId = (int) aNodeToGhs3dIdMap.size() + 1; // MG-Tetra ids count from 1
1099 aNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1102 if ( !anElemSet.empty() &&
1103 (*anElemSet.begin())->IsQuadratic() &&
1104 theProxyMesh.NbProxySubMeshes() > 0 )
1106 // add medium nodes of proxy triangles to theHelper (#16843)
1107 for ( elemSetIt = anElemSet.begin(); elemSetIt != anElemSet.end(); ++elemSetIt )
1108 theHelper.AddTLinks( static_cast< const SMDS_MeshFace* >( *elemSetIt ));
1111 /* EDGES ========================== */
1113 // Iterate over the enforced edges
1114 for(elemIt = theEnforcedEdges.begin() ; elemIt != theEnforcedEdges.end() ; ++elemIt) {
1115 elem = elemIt->first;
1117 nodeIt = elem->nodesIterator();
1119 while ( nodeIt->more() && nbNodes-- ) {
1121 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1122 // Test if point is inside shape to mesh
1123 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1124 TopAbs_State result = pntCls->GetPointState( myPoint );
1125 if ( result == TopAbs_OUT ) {
1127 theInvalidEnforcedFlags |= FLAG_BAD_ENF_EDGE;
1130 aNodeToTopAbs_StateMap.insert( make_pair( node, result ));
1133 nodeIt = elem->nodesIterator();
1136 while ( nodeIt->more() && nbNodes-- ) {
1138 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1139 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1140 nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1142 std::cout << "Node at "<<node->X()<<", "<<node->Y()<<", "<<node->Z()<<std::endl;
1143 std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
1145 if (nbFoundElems ==0) {
1146 if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1147 newId = int( aNodeToGhs3dIdMap.size() + anEnforcedNodeToGhs3dIdMap.size() + 1 ); // MG-Tetra ids count from 1
1148 anEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1151 else if (nbFoundElems ==1) {
1152 const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1153 newId = (*aNodeToGhs3dIdMap.find(existingNode)).second;
1154 anExistingEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1159 std::cout << "MG-Tetra node ID: "<<newId<<std::endl;
1163 theKeptEnforcedEdges.insert(elem);
1165 theInvalidEnforcedFlags |= FLAG_BAD_ENF_EDGE;
1169 /* ENFORCED TRIANGLES ========================== */
1171 // Iterate over the enforced triangles
1172 for(elemIt = theEnforcedTriangles.begin() ; elemIt != theEnforcedTriangles.end() ; ++elemIt) {
1173 elem = elemIt->first;
1175 nodeIt = elem->nodesIterator();
1177 while ( nodeIt->more() && nbNodes--) {
1179 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1180 // Test if point is inside shape to mesh
1181 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1182 TopAbs_State result = pntCls->GetPointState( myPoint );
1183 if ( result == TopAbs_OUT ) {
1185 theInvalidEnforcedFlags |= FLAG_BAD_ENF_TRIA;
1188 aNodeToTopAbs_StateMap.insert( make_pair( node, result ));
1191 nodeIt = elem->nodesIterator();
1194 while ( nodeIt->more() && nbNodes--) {
1196 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1197 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1198 nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1200 std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
1202 if (nbFoundElems ==0) {
1203 if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1204 newId = int( aNodeToGhs3dIdMap.size() + anEnforcedNodeToGhs3dIdMap.size() + 1 ); // MG-Tetra ids count from 1
1205 anEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1208 else if (nbFoundElems ==1) {
1209 const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1210 newId = (*aNodeToGhs3dIdMap.find(existingNode)).second;
1211 anExistingEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1216 std::cout << "MG-Tetra node ID: "<<newId<<std::endl;
1220 theKeptEnforcedTriangles.insert(elem);
1222 theInvalidEnforcedFlags |= FLAG_BAD_ENF_TRIA;
1226 // put nodes to theNodeByGhs3dId vector
1228 std::cout << "aNodeToGhs3dIdMap.size(): "<<aNodeToGhs3dIdMap.size()<<std::endl;
1230 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
1231 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
1232 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
1234 // std::cout << "n2id->first: "<<n2id->first<<std::endl;
1235 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // MG-Tetra ids count from 1
1238 // put nodes to anEnforcedNodeToGhs3dIdMap vector
1240 std::cout << "anEnforcedNodeToGhs3dIdMap.size(): "<<anEnforcedNodeToGhs3dIdMap.size()<<std::endl;
1242 theEnforcedNodeByGhs3dId.resize( anEnforcedNodeToGhs3dIdMap.size());
1243 n2id = anEnforcedNodeToGhs3dIdMap.begin();
1244 for ( ; n2id != anEnforcedNodeToGhs3dIdMap.end(); ++ n2id)
1246 if (n2id->second > (int)aNodeToGhs3dIdMap.size()) {
1247 theEnforcedNodeByGhs3dId[ n2id->second - aNodeToGhs3dIdMap.size() - 1 ] = n2id->first; // MG-Tetra ids count from 1
1252 /* ========================== NODES ========================== */
1253 vector<const SMDS_MeshNode*> theOrderedNodes, theRequiredNodes;
1254 std::set< std::vector<double> > nodesCoords;
1255 vector<const SMDS_MeshNode*>::const_iterator ghs3dNodeIt = theNodeByGhs3dId.begin();
1256 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
1258 (theNodeByGhs3dId.size() <= 1) ? tmpStr = " node" : " nodes";
1259 std::cout << theNodeByGhs3dId.size() << tmpStr << " from mesh ..." << std::endl;
1260 for ( ; ghs3dNodeIt != after; ++ghs3dNodeIt )
1262 const SMDS_MeshNode* node = *ghs3dNodeIt;
1263 std::vector<double> coords;
1264 coords.push_back(node->X());
1265 coords.push_back(node->Y());
1266 coords.push_back(node->Z());
1267 nodesCoords.insert(coords);
1268 theOrderedNodes.push_back(node);
1271 // Iterate over the enforced nodes given by enforced elements
1272 ghs3dNodeIt = theEnforcedNodeByGhs3dId.begin();
1273 after = theEnforcedNodeByGhs3dId.end();
1274 (theEnforcedNodeByGhs3dId.size() <= 1) ? tmpStr = " node" : " nodes";
1275 std::cout << theEnforcedNodeByGhs3dId.size() << tmpStr << " from enforced elements ..." << std::endl;
1276 for ( ; ghs3dNodeIt != after; ++ghs3dNodeIt )
1278 const SMDS_MeshNode* node = *ghs3dNodeIt;
1279 std::vector<double> coords;
1280 coords.push_back(node->X());
1281 coords.push_back(node->Y());
1282 coords.push_back(node->Z());
1284 std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1287 if (nodesCoords.find(coords) != nodesCoords.end()) {
1288 // node already exists in original mesh
1290 std::cout << " found" << std::endl;
1295 if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1296 // node already exists in enforced vertices
1298 std::cout << " found" << std::endl;
1303 // gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1304 // nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1305 // if (nbFoundElems ==0) {
1306 // std::cout << " not found" << std::endl;
1307 // if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1308 // nodesCoords.insert(coords);
1309 // theOrderedNodes.push_back(node);
1313 // std::cout << " found in initial mesh" << std::endl;
1314 // const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1315 // nodesCoords.insert(coords);
1316 // theOrderedNodes.push_back(existingNode);
1320 std::cout << " not found" << std::endl;
1323 nodesCoords.insert(coords);
1324 theOrderedNodes.push_back(node);
1325 // theRequiredNodes.push_back(node);
1329 // Iterate over the enforced nodes
1330 GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap::const_iterator enfNodeIt;
1331 (theEnforcedNodes.size() <= 1) ? tmpStr = " node" : " nodes";
1332 std::cout << theEnforcedNodes.size() << tmpStr << " from enforced nodes ..." << std::endl;
1333 for(enfNodeIt = theEnforcedNodes.begin() ; enfNodeIt != theEnforcedNodes.end() ; ++enfNodeIt)
1335 const SMDS_MeshNode* node = enfNodeIt->first;
1336 std::vector<double> coords;
1337 coords.push_back(node->X());
1338 coords.push_back(node->Y());
1339 coords.push_back(node->Z());
1341 std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1344 // Test if point is inside shape to mesh
1345 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1346 TopAbs_State result = pntCls->GetPointState( myPoint );
1347 if ( result == TopAbs_OUT ) {
1349 std::cout << " out of volume" << std::endl;
1351 theInvalidEnforcedFlags |= FLAG_BAD_ENF_NODE;
1355 if (nodesCoords.find(coords) != nodesCoords.end()) {
1357 std::cout << " found in nodesCoords" << std::endl;
1359 // theRequiredNodes.push_back(node);
1363 if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1365 std::cout << " found in theEnforcedVertices" << std::endl;
1370 // nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1371 // if (nbFoundElems ==0) {
1372 // std::cout << " not found" << std::endl;
1373 // if (result == TopAbs_IN) {
1374 // nodesCoords.insert(coords);
1375 // theRequiredNodes.push_back(node);
1379 // std::cout << " found in initial mesh" << std::endl;
1380 // const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1381 // // nodesCoords.insert(coords);
1382 // theRequiredNodes.push_back(existingNode);
1387 // if (pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems) == 0)
1390 // if ( result != TopAbs_IN )
1394 std::cout << " not found" << std::endl;
1396 nodesCoords.insert(coords);
1397 // theOrderedNodes.push_back(node);
1398 theRequiredNodes.push_back(node);
1400 int requiredNodes = (int) theRequiredNodes.size();
1403 std::vector<std::vector<double> > ReqVerTab;
1404 if (nbEnforcedVertices) {
1405 // ReqVerTab.clear();
1406 (nbEnforcedVertices <= 1) ? tmpStr = " node" : " nodes";
1407 std::cout << nbEnforcedVertices << tmpStr << " from enforced vertices ..." << std::endl;
1408 // Iterate over the enforced vertices
1409 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
1410 double x = vertexIt->first[0];
1411 double y = vertexIt->first[1];
1412 double z = vertexIt->first[2];
1413 // Test if point is inside shape to mesh
1414 gp_Pnt myPoint(x,y,z);
1415 TopAbs_State result = pntCls->GetPointState( myPoint );
1416 if ( result == TopAbs_OUT )
1418 std::cout << "Warning: enforced vertex at ( " << x << "," << y << "," << z << " ) is out of the meshed domain!!!" << std::endl;
1419 theInvalidEnforcedFlags |= FLAG_BAD_ENF_VERT;
1422 std::vector<double> coords;
1423 coords.push_back(x);
1424 coords.push_back(y);
1425 coords.push_back(z);
1426 ReqVerTab.push_back(coords);
1427 enfVertexSizes.push_back(vertexIt->second);
1434 std::cout << "Begin writting required nodes in GmfVertices" << std::endl;
1435 std::cout << "Nb vertices: " << theOrderedNodes.size() << std::endl;
1436 MGInput->GmfSetKwd( idx, GmfVertices, int( theOrderedNodes.size()/*+solSize*/));
1437 for (ghs3dNodeIt = theOrderedNodes.begin();ghs3dNodeIt != theOrderedNodes.end();++ghs3dNodeIt) {
1438 MGInput->GmfSetLin( idx, GmfVertices, (*ghs3dNodeIt)->X(), (*ghs3dNodeIt)->Y(), (*ghs3dNodeIt)->Z(), dummyint);
1441 std::cout << "End writting required nodes in GmfVertices" << std::endl;
1443 if (requiredNodes + solSize) {
1444 std::cout << "Begin writting in req and sol file" << std::endl;
1445 aNodeGroupByGhs3dId.resize( requiredNodes + solSize );
1446 idxRequired = MGInput->GmfOpenMesh( theRequiredFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1450 idxSol = MGInput->GmfOpenMesh( theSolFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1454 int TypTab[] = {GmfSca};
1455 double ValTab[] = {0.0};
1456 MGInput->GmfSetKwd( idxRequired, GmfVertices, requiredNodes + solSize);
1457 MGInput->GmfSetKwd( idxSol, GmfSolAtVertices, requiredNodes + solSize, 1, TypTab);
1458 // int usedEnforcedNodes = 0;
1459 // std::string gn = "";
1460 for (ghs3dNodeIt = theRequiredNodes.begin();ghs3dNodeIt != theRequiredNodes.end();++ghs3dNodeIt) {
1461 MGInput->GmfSetLin( idxRequired, GmfVertices, (*ghs3dNodeIt)->X(), (*ghs3dNodeIt)->Y(), (*ghs3dNodeIt)->Z(), dummyint);
1462 MGInput->GmfSetLin( idxSol, GmfSolAtVertices, ValTab);
1463 if (theEnforcedNodes.find((*ghs3dNodeIt)) != theEnforcedNodes.end())
1464 gn = theEnforcedNodes.find((*ghs3dNodeIt))->second;
1465 aNodeGroupByGhs3dId[usedEnforcedNodes] = gn;
1466 usedEnforcedNodes++;
1469 for (int i=0;i<solSize;i++) {
1470 std::cout << ReqVerTab[i][0] <<" "<< ReqVerTab[i][1] << " "<< ReqVerTab[i][2] << std::endl;
1472 std::cout << "enfVertexSizes.at("<<i<<"): " << enfVertexSizes.at(i) << std::endl;
1474 double solTab[] = {enfVertexSizes.at(i)};
1475 MGInput->GmfSetLin( idxRequired, GmfVertices, ReqVerTab[i][0], ReqVerTab[i][1], ReqVerTab[i][2], dummyint);
1476 MGInput->GmfSetLin( idxSol, GmfSolAtVertices, solTab);
1477 aNodeGroupByGhs3dId[usedEnforcedNodes] = enfVerticesWithGroup.find(ReqVerTab[i])->second;
1479 std::cout << "aNodeGroupByGhs3dId["<<usedEnforcedNodes<<"] = \""<<aNodeGroupByGhs3dId[usedEnforcedNodes]<<"\""<<std::endl;
1481 usedEnforcedNodes++;
1483 std::cout << "End writting in req and sol file" << std::endl;
1486 int nedge[2], ntri[3];
1489 int usedEnforcedEdges = 0;
1490 if (theKeptEnforcedEdges.size()) {
1491 anEdgeGroupByGhs3dId.resize( theKeptEnforcedEdges.size() );
1492 // idxRequired = MGInput->GmfOpenMesh( theRequiredFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1493 // if (!idxRequired)
1495 MGInput->GmfSetKwd( idx, GmfEdges, (int) theKeptEnforcedEdges.size());
1496 // MGInput->GmfSetKwd( idxRequired, GmfEdges, theKeptEnforcedEdges.size());
1497 for(elemSetIt = theKeptEnforcedEdges.begin() ; elemSetIt != theKeptEnforcedEdges.end() ; ++elemSetIt) {
1498 elem = (*elemSetIt);
1499 nodeIt = elem->nodesIterator();
1501 while ( nodeIt->more() ) {
1503 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1504 map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToGhs3dIdMap.find(node);
1505 if (it == anEnforcedNodeToGhs3dIdMap.end()) {
1506 it = anExistingEnforcedNodeToGhs3dIdMap.find(node);
1507 if (it == anEnforcedNodeToGhs3dIdMap.end())
1508 throw "Node not found";
1510 nedge[index] = it->second;
1513 MGInput->GmfSetLin( idx, GmfEdges, nedge[0], nedge[1], dummyint);
1514 anEdgeGroupByGhs3dId[usedEnforcedEdges] = theEnforcedEdges.find(elem)->second;
1515 // MGInput->GmfSetLin( idxRequired, GmfEdges, nedge[0], nedge[1], dummyint);
1516 usedEnforcedEdges++;
1521 if (usedEnforcedEdges) {
1522 MGInput->GmfSetKwd( idx, GmfRequiredEdges, usedEnforcedEdges);
1523 for (int enfID=1;enfID<=usedEnforcedEdges;enfID++) {
1524 MGInput->GmfSetLin( idx, GmfRequiredEdges, enfID);
1529 int usedEnforcedTriangles = 0;
1530 if (anElemSet.size()+theKeptEnforcedTriangles.size()) {
1531 aFaceGroupByGhs3dId.resize( anElemSet.size()+theKeptEnforcedTriangles.size() );
1532 MGInput->GmfSetKwd( idx, GmfTriangles, int( anElemSet.size()+theKeptEnforcedTriangles.size() ));
1534 for(elemSetIt = anElemSet.begin() ; elemSetIt != anElemSet.end() ; ++elemSetIt,++k) {
1535 elem = (*elemSetIt);
1536 theFaceByGhs3dId.push_back( elem );
1537 nodeIt = elem->nodesIterator();
1539 for ( int j = 0; j < 3; ++j ) {
1541 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1542 map< const SMDS_MeshNode*,int >::iterator it = aNodeToGhs3dIdMap.find(node);
1543 if (it == aNodeToGhs3dIdMap.end())
1544 throw "Node not found";
1545 ntri[index] = it->second;
1548 MGInput->GmfSetLin( idx, GmfTriangles, ntri[0], ntri[1], ntri[2], dummyint);
1549 aFaceGroupByGhs3dId[k] = "";
1551 if ( !theHelper.GetMesh()->HasShapeToMesh() )
1552 SMESHUtils::FreeVector( theFaceByGhs3dId );
1553 if (theKeptEnforcedTriangles.size()) {
1554 for(elemSetIt = theKeptEnforcedTriangles.begin() ; elemSetIt != theKeptEnforcedTriangles.end() ; ++elemSetIt,++k) {
1555 elem = (*elemSetIt);
1556 nodeIt = elem->nodesIterator();
1558 for ( int j = 0; j < 3; ++j ) {
1560 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1561 map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToGhs3dIdMap.find(node);
1562 if (it == anEnforcedNodeToGhs3dIdMap.end()) {
1563 it = anExistingEnforcedNodeToGhs3dIdMap.find(node);
1564 if (it == anEnforcedNodeToGhs3dIdMap.end())
1565 throw "Node not found";
1567 ntri[index] = it->second;
1570 MGInput->GmfSetLin( idx, GmfTriangles, ntri[0], ntri[1], ntri[2], dummyint);
1571 aFaceGroupByGhs3dId[k] = theEnforcedTriangles.find(elem)->second;
1572 usedEnforcedTriangles++;
1578 if (usedEnforcedTriangles) {
1579 MGInput->GmfSetKwd( idx, GmfRequiredTriangles, usedEnforcedTriangles);
1580 for (int enfID=1;enfID<=usedEnforcedTriangles;enfID++)
1581 MGInput->GmfSetLin( idx, GmfRequiredTriangles, int( anElemSet.size()+enfID ));
1584 // close input files and change results files permissions to user only
1585 MGInput->GmfCloseMesh(idx);
1586 chmodUserOnly(theMeshFileName);
1589 MGInput->GmfCloseMesh(idxRequired);
1590 chmodUserOnly(theRequiredFileName);
1594 MGInput->GmfCloseMesh(idxSol);
1595 chmodUserOnly(theRequiredFileName);
1602 //=============================================================================
1604 *Here we are going to use the MG-Tetra mesher with geometry
1606 //=============================================================================
1607 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1608 const TopoDS_Shape& theShape)
1611 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1613 // a unique working file name
1614 // to avoid access to the same files by eg different users
1615 _genericName = GHS3DPlugin_Hypothesis::GetFileName(_hyp);
1616 TCollection_AsciiString aGenericName((char*) _genericName.c_str() );
1617 TCollection_AsciiString aGenericNameRequired = aGenericName + "_required";
1619 TCollection_AsciiString aLogFileName = aGenericName + ".log"; // log
1620 TCollection_AsciiString aResultFileName;
1622 TCollection_AsciiString aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1623 aGMFFileName = aGenericName + ".mesh"; // GMF mesh file
1624 aResultFileName = aGenericName + "Vol.mesh"; // GMF mesh file
1625 aResSolFileName = aGenericName + "Vol.sol"; // GMF mesh file
1626 aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1627 aSolFileName = aGenericNameRequired + ".sol"; // GMF solution file
1629 std::map <int,int> aNodeId2NodeIndexMap, aSmdsToGhs3dIdMap, anEnforcedNodeIdToGhs3dIdMap;
1630 std::map <int, int> nodeID2nodeIndexMap;
1631 std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1632 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues coordsSizeMap = GHS3DPlugin_Hypothesis::GetEnforcedVerticesCoordsSize(_hyp);
1633 GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = GHS3DPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1634 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = GHS3DPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1635 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = GHS3DPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1636 GHS3DPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = GHS3DPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1638 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList enfVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1639 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1640 std::vector<double> coords;
1642 for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1644 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertex* enfVertex = (*enfVerIt);
1645 if (enfVertex->coords.size()) {
1646 coordsSizeMap.insert(make_pair(enfVertex->coords,enfVertex->size));
1647 enfVerticesWithGroup.insert(make_pair(enfVertex->coords,enfVertex->groupName));
1650 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
1651 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1653 if (it.Value().ShapeType() == TopAbs_VERTEX){
1654 gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1655 coords.push_back(aPnt.X());
1656 coords.push_back(aPnt.Y());
1657 coords.push_back(aPnt.Z());
1658 if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1659 coordsSizeMap.insert(make_pair(coords,enfVertex->size));
1660 enfVerticesWithGroup.insert(make_pair(coords,enfVertex->groupName));
1666 size_t nbEnforcedVertices = coordsSizeMap.size();
1667 size_t nbEnforcedNodes = enforcedNodes.size();
1670 (nbEnforcedNodes <= 1) ? tmpStr = "node" : "nodes";
1671 std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1672 (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : "vertices";
1673 std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1675 SMESH_MesherHelper helper( theMesh );
1676 helper.SetSubShape( theShape );
1678 std::vector <const SMDS_MeshNode*> aNodeByGhs3dId, anEnforcedNodeByGhs3dId;
1679 std::vector <const SMDS_MeshElement*> aFaceByGhs3dId;
1680 std::map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
1681 std::vector<std::string> aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId;
1683 MG_Tetra_API mgTetra( _computeCanceled, _progress );
1685 _isLibUsed = mgTetra.IsLibrary();
1686 if ( theMesh.NbQuadrangles() > 0 )
1687 _progressAdvance /= 10;
1688 if ( _viscousLayersHyp )
1689 _progressAdvance /= 10;
1691 // proxyMesh must live till readGMFFile() as a proxy face can be used by
1692 // MG-Tetra for domain indication
1693 SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1695 // make prisms on quadrangles and viscous layers
1696 if ( theMesh.NbQuadrangles() > 0 || _viscousLayersHyp )
1698 vector<SMESH_ProxyMesh::Ptr> components;
1699 for (expBox.ReInit(); expBox.More(); expBox.Next())
1701 if ( _viscousLayersHyp )
1703 proxyMesh = _viscousLayersHyp->Compute( theMesh, expBox.Current() );
1706 if ( theMesh.NbQuadrangles() == 0 )
1707 components.push_back( proxyMesh );
1709 if ( theMesh.NbQuadrangles() > 0 )
1711 StdMeshers_QuadToTriaAdaptor* q2t = new StdMeshers_QuadToTriaAdaptor;
1712 Ok = q2t->Compute( theMesh, expBox.Current(), proxyMesh.get() );
1713 components.push_back( SMESH_ProxyMesh::Ptr( q2t ));
1718 proxyMesh.reset( new SMESH_ProxyMesh( components ));
1720 // build viscous layers
1721 // else if ( _viscousLayersHyp )
1723 // proxyMesh = _viscousLayersHyp->Compute( theMesh, theShape );
1724 // if ( !proxyMesh )
1728 int anInvalidEnforcedFlags = 0;
1729 Ok = writeGMFFile(&mgTetra,
1730 aGMFFileName.ToCString(),
1731 aRequiredVerticesFileName.ToCString(),
1732 aSolFileName.ToCString(),
1734 aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
1735 aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
1736 enforcedNodes, enforcedEdges, enforcedTriangles,
1737 enfVerticesWithGroup, coordsSizeMap, anInvalidEnforcedFlags);
1739 // Write aSmdsToGhs3dIdMap to temp file
1740 TCollection_AsciiString aSmdsToGhs3dIdMapFileName;
1741 aSmdsToGhs3dIdMapFileName = aGenericName + ".ids"; // ids relation
1742 ofstream aIdsFile ( aSmdsToGhs3dIdMapFileName.ToCString() , ios::out);
1743 if ( !aIdsFile.rdbuf()->is_open() ) {
1744 INFOS( "Can't write into " << aSmdsToGhs3dIdMapFileName);
1745 //return error(SMESH_Comment("Can't write into ") << aSmdsToGhs3dIdMapFileName);
1749 INFOS( "Writing ids relation into " << aSmdsToGhs3dIdMapFileName);
1750 aIdsFile << "Smds MG-Tetra" << std::endl;
1751 map <int,int>::const_iterator myit;
1752 for (myit=aSmdsToGhs3dIdMap.begin() ; myit != aSmdsToGhs3dIdMap.end() ; ++myit) {
1753 aIdsFile << myit->first << " " << myit->second << std::endl;
1757 chmodUserOnly(aSmdsToGhs3dIdMapFileName.ToCString());
1759 if ( !_keepFiles ) {
1760 removeFile( aGMFFileName );
1761 removeFile( aRequiredVerticesFileName );
1762 removeFile( aSolFileName );
1763 removeFile( aSmdsToGhs3dIdMapFileName );
1765 return error(COMPERR_BAD_INPUT_MESH);
1767 removeFile( aResultFileName ); // needed for boundary recovery module usage
1769 // -----------------
1770 // run MG-Tetra mesher
1771 // -----------------
1773 TCollection_AsciiString cmd = GHS3DPlugin_Hypothesis::CommandToRun( _hyp, true, mgTetra.IsExecutable() ).c_str();
1775 if ( mgTetra.IsExecutable() )
1777 cmd += TCollection_AsciiString(" --in ") + aGMFFileName;
1778 if ( nbEnforcedVertices + nbEnforcedNodes)
1779 cmd += TCollection_AsciiString(" --required_vertices ") + aGenericNameRequired;
1780 cmd += TCollection_AsciiString(" --out ") + aResultFileName;
1782 if ( !_logInStandardOutput )
1784 mgTetra.SetLogFile( aLogFileName.ToCString() );
1785 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1789 BRIEF_INFOS("MG-Tetra execution...")
1792 _computeCanceled = false;
1795 Ok = mgTetra.Compute( cmd.ToCString(), errStr ); // run
1797 if ( _logInStandardOutput && mgTetra.IsLibrary() ) {
1799 BRIEF_INFOS(mgTetra.GetLog());
1804 BRIEF_INFOS("End of MG-Tetra execution !");
1808 // change results files permissions to user only
1809 chmodUserOnly(aLogFileName.ToCString());
1812 chmodUserOnly(aResultFileName.ToCString());
1813 chmodUserOnly(aResSolFileName.ToCString());
1820 GHS3DPlugin_Hypothesis::TSetStrings groupsToRemove = GHS3DPlugin_Hypothesis::GetGroupsToRemove(_hyp);
1822 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1823 const bool toMakeGroupsOfDomains = GHS3DPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
1825 helper.IsQuadraticSubMesh( theShape );
1826 helper.SetElementsOnShape( false );
1828 Ok = readGMFFile(&mgTetra,
1829 aResultFileName.ToCString(),
1831 &helper, aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
1832 aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
1833 groupsToRemove, toMakeGroupsOfDomains, toMeshHoles);
1835 removeEmptyGroupsOfDomains( helper.GetMesh(), /*notEmptyAsWell =*/ !toMakeGroupsOfDomains );
1839 // ---------------------
1840 // remove working files
1841 // ---------------------
1845 if ( anInvalidEnforcedFlags )
1846 error( COMPERR_WARNING, flagsToErrorStr( anInvalidEnforcedFlags ));
1847 if ( _removeLogOnSuccess )
1848 removeFile( aLogFileName );
1849 // if ( _hyp && _hyp->GetToMakeGroupsOfDomains() )
1850 // error( COMPERR_WARNING, "'toMakeGroupsOfDomains' is ignored since the mesh is on shape" );
1852 else if ( mgTetra.HasLog() )
1854 if( _computeCanceled )
1855 error( "interruption initiated by user" );
1858 // get problem description from the log file
1859 _Ghs2smdsConvertor conv( aNodeByGhs3dId, proxyMesh );
1860 error( getErrorDescription( _logInStandardOutput ? 0 : aLogFileName.ToCString(),
1861 mgTetra.GetLog(), conv ));
1864 else if ( !errStr.empty() )
1866 // the log file is empty
1867 removeFile( aLogFileName );
1868 INFOS( "MG-Tetra Error, " << errStr);
1869 error(COMPERR_ALGO_FAILED, errStr);
1872 if ( !_keepFiles ) {
1873 if (! Ok && _computeCanceled )
1874 removeFile( aLogFileName );
1875 removeFile( aGMFFileName );
1876 removeFile( aRequiredVerticesFileName );
1877 removeFile( aSolFileName );
1878 removeFile( aResSolFileName );
1879 removeFile( aResultFileName );
1880 removeFile( aSmdsToGhs3dIdMapFileName );
1882 if ( mgTetra.IsExecutable() )
1884 std::cout << "<" << aResultFileName.ToCString() << "> MG-Tetra output file ";
1886 std::cout << "not ";
1887 std::cout << "treated !" << std::endl;
1888 std::cout << std::endl;
1892 std::cout << "MG-Tetra " << ( Ok ? "succeeded" : "failed") << std::endl;
1897 //=============================================================================
1899 *Here we are going to use the MG-Tetra mesher w/o geometry
1901 //=============================================================================
1903 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1904 SMESH_MesherHelper* theHelper)
1906 theHelper->IsQuadraticSubMesh( theHelper->GetSubShape() );
1908 // a unique working file name
1909 // to avoid access to the same files by eg different users
1910 _genericName = GHS3DPlugin_Hypothesis::GetFileName(_hyp);
1911 TCollection_AsciiString aGenericName((char*) _genericName.c_str() );
1912 TCollection_AsciiString aGenericNameRequired = aGenericName + "_required";
1914 TCollection_AsciiString aLogFileName = aGenericName + ".log"; // log
1915 TCollection_AsciiString aResultFileName;
1918 TCollection_AsciiString aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1919 aGMFFileName = aGenericName + ".mesh"; // GMF mesh file
1920 aResultFileName = aGenericName + "Vol.mesh"; // GMF mesh file
1921 aResSolFileName = aGenericName + "Vol.sol"; // GMF mesh file
1922 aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1923 aSolFileName = aGenericNameRequired + ".sol"; // GMF solution file
1925 std::map <int, int> nodeID2nodeIndexMap;
1926 std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1927 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues coordsSizeMap;
1928 TopoDS_Shape GeomShape;
1929 std::vector<double> coords;
1931 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertex* enfVertex;
1933 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList enfVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1934 GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1936 for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1938 enfVertex = (*enfVerIt);
1939 if (enfVertex->coords.size()) {
1940 coordsSizeMap.insert(make_pair(enfVertex->coords,enfVertex->size));
1941 enfVerticesWithGroup.insert(make_pair(enfVertex->coords,enfVertex->groupName));
1944 GeomShape = entryToShape(enfVertex->geomEntry);
1945 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1947 if (it.Value().ShapeType() == TopAbs_VERTEX){
1948 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1949 coords.push_back(aPnt.X());
1950 coords.push_back(aPnt.Y());
1951 coords.push_back(aPnt.Z());
1952 if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1953 coordsSizeMap.insert(make_pair(coords,enfVertex->size));
1954 enfVerticesWithGroup.insert(make_pair(coords,enfVertex->groupName));
1961 GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = GHS3DPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1962 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = GHS3DPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1963 GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = GHS3DPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1964 GHS3DPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = GHS3DPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1968 size_t nbEnforcedVertices = coordsSizeMap.size();
1969 size_t nbEnforcedNodes = enforcedNodes.size();
1970 (nbEnforcedNodes <= 1) ? tmpStr = "node" : tmpStr = "nodes";
1971 std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1972 (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : tmpStr = "vertices";
1973 std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1975 std::vector <const SMDS_MeshNode*> aNodeByGhs3dId, anEnforcedNodeByGhs3dId;
1976 std::vector <const SMDS_MeshElement*> aFaceByGhs3dId;
1977 std::map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
1978 std::vector<std::string> aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId;
1981 MG_Tetra_API mgTetra( _computeCanceled, _progress );
1983 _isLibUsed = mgTetra.IsLibrary();
1984 if ( theMesh.NbQuadrangles() > 0 )
1985 _progressAdvance /= 10;
1987 // proxyMesh must live till readGMFFile() as a proxy face can be used by
1988 // MG-Tetra for domain indication
1989 SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1990 if ( theMesh.NbQuadrangles() > 0 )
1992 StdMeshers_QuadToTriaAdaptor* aQuad2Trias = new StdMeshers_QuadToTriaAdaptor;
1993 Ok = aQuad2Trias->Compute( theMesh );
1994 proxyMesh.reset( aQuad2Trias );
1999 int anInvalidEnforcedFlags = 0;
2000 Ok = writeGMFFile(&mgTetra,
2001 aGMFFileName.ToCString(), aRequiredVerticesFileName.ToCString(), aSolFileName.ToCString(),
2002 *proxyMesh, *theHelper,
2003 aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
2004 aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
2005 enforcedNodes, enforcedEdges, enforcedTriangles,
2006 enfVerticesWithGroup, coordsSizeMap, anInvalidEnforcedFlags);
2008 // -----------------
2009 // run MG-Tetra mesher
2010 // -----------------
2012 TCollection_AsciiString cmd = GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false, mgTetra.IsExecutable() ).c_str();
2014 if ( mgTetra.IsExecutable() )
2016 cmd += TCollection_AsciiString(" --in ") + aGMFFileName;
2017 if ( nbEnforcedVertices + nbEnforcedNodes)
2018 cmd += TCollection_AsciiString(" --required_vertices ") + aGenericNameRequired;
2019 cmd += TCollection_AsciiString(" --out ") + aResultFileName;
2021 if ( !_logInStandardOutput )
2023 mgTetra.SetLogFile( aLogFileName.ToCString() );
2024 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
2028 BRIEF_INFOS("MG-Tetra execution...")
2031 _computeCanceled = false;
2034 Ok = mgTetra.Compute( cmd.ToCString(), errStr ); // run
2036 if ( _logInStandardOutput && mgTetra.IsLibrary() ) {
2038 BRIEF_INFOS(mgTetra.GetLog());
2043 BRIEF_INFOS("End of MG-Tetra execution !");
2050 GHS3DPlugin_Hypothesis::TSetStrings groupsToRemove = GHS3DPlugin_Hypothesis::GetGroupsToRemove(_hyp);
2051 const bool toMakeGroupsOfDomains = GHS3DPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
2053 Ok = Ok && readGMFFile(&mgTetra,
2054 aResultFileName.ToCString(),
2056 theHelper, aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
2057 aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
2058 groupsToRemove, toMakeGroupsOfDomains);
2060 updateMeshGroups(theHelper->GetMesh(), groupsToRemove);
2061 removeEmptyGroupsOfDomains( theHelper->GetMesh(), /*notEmptyAsWell =*/ !toMakeGroupsOfDomains );
2064 GHS3DPlugin_Hypothesis* that = (GHS3DPlugin_Hypothesis*)this->_hyp;
2066 that->ClearGroupsToRemove();
2068 // ---------------------
2069 // remove working files
2070 // ---------------------
2074 if ( anInvalidEnforcedFlags )
2075 error( COMPERR_WARNING, flagsToErrorStr( anInvalidEnforcedFlags ));
2076 if ( _removeLogOnSuccess )
2077 removeFile( aLogFileName );
2079 //if ( !toMakeGroupsOfDomains && _hyp && _hyp->GetToMakeGroupsOfDomains() )
2080 //error( COMPERR_WARNING, "'toMakeGroupsOfDomains' is ignored since 'toMeshHoles' is OFF." );
2082 else if ( mgTetra.HasLog() )
2084 if( _computeCanceled )
2085 error( "interruption initiated by user" );
2088 // get problem description from the log file
2089 _Ghs2smdsConvertor conv( aNodeByGhs3dId, proxyMesh );
2090 error( getErrorDescription( _logInStandardOutput ? 0 : aLogFileName.ToCString(),
2091 mgTetra.GetLog(), conv ));
2095 // the log file is empty
2096 removeFile( aLogFileName );
2097 INFOS( "MG-Tetra Error, " << errStr);
2098 error(COMPERR_ALGO_FAILED, errStr);
2103 if (! Ok && _computeCanceled)
2104 removeFile( aLogFileName );
2105 removeFile( aGMFFileName );
2106 removeFile( aResultFileName );
2107 removeFile( aRequiredVerticesFileName );
2108 removeFile( aSolFileName );
2109 removeFile( aResSolFileName );
2114 void GHS3DPlugin_GHS3D::CancelCompute()
2116 _computeCanceled = true;
2117 #if !defined WIN32 && !defined __APPLE__
2118 std::string cmd = "ps xo pid,args | grep " + _genericName;
2119 //cmd += " | grep -e \"^ *[0-9]\\+ \\+" + GHS3DPlugin_Hypothesis::GetExeName() + "\"";
2120 cmd += " | awk '{print $1}' | xargs kill -9 > /dev/null 2>&1";
2121 system( cmd.c_str() );
2125 //================================================================================
2127 * \brief Provide human readable text by error code reported by MG-Tetra
2129 //================================================================================
2131 static const char* translateError(const int errNum)
2135 return "The surface mesh includes a face of type other than edge, "
2136 "triangle or quadrilateral. This face type is not supported.";
2138 return "Not enough memory for the face table.";
2140 return "Not enough memory.";
2142 return "Not enough memory.";
2144 return "Face is ignored.";
2146 return "End of file. Some data are missing in the file.";
2148 return "Read error on the file. There are wrong data in the file.";
2150 return "the metric file is inadequate (dimension other than 3).";
2152 return "the metric file is inadequate (values not per vertices).";
2154 return "the metric file contains more than one field.";
2156 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
2157 "value of number of mesh vertices in the \".noboite\" file.";
2159 return "Too many sub-domains.";
2161 return "the number of vertices is negative or null.";
2163 return "the number of faces is negative or null.";
2165 return "A face has a null vertex.";
2167 return "incompatible data.";
2169 return "the number of vertices is negative or null.";
2171 return "the number of vertices is negative or null (in the \".mesh\" file).";
2173 return "the number of faces is negative or null.";
2175 return "A face appears more than once in the input surface mesh.";
2177 return "An edge appears more than once in the input surface mesh.";
2179 return "A face has a vertex negative or null.";
2181 return "NOT ENOUGH MEMORY.";
2183 return "Not enough available memory.";
2185 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
2186 "in terms of quality or the input list of points is wrong.";
2188 return "Some vertices are too close to one another or coincident.";
2190 return "Some vertices are too close to one another or coincident.";
2192 return "A vertex cannot be inserted.";
2194 return "There are at least two points considered as coincident.";
2196 return "Some vertices are too close to one another or coincident.";
2198 return "The surface mesh regeneration step has failed.";
2200 return "Constrained edge cannot be enforced.";
2202 return "Constrained face cannot be enforced.";
2204 return "Missing faces.";
2206 return "No guess to start the definition of the connected component(s).";
2208 return "The surface mesh includes at least one hole. The domain is not well defined.";
2210 return "Impossible to define a component.";
2212 return "The surface edge intersects another surface edge.";
2214 return "The surface edge intersects the surface face.";
2216 return "One boundary point lies within a surface face.";
2218 return "One surface edge intersects a surface face.";
2220 return "One boundary point lies within a surface edge.";
2222 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
2223 "to too many swaps.";
2225 return "Edge is unique (i.e., bounds a hole in the surface).";
2227 return "Presumably, the surface mesh is not compatible with the domain being processed.";
2229 return "Too many components, too many sub-domain.";
2231 return "The surface mesh includes at least one hole. "
2232 "Therefore there is no domain properly defined.";
2234 return "Statistics.";
2236 return "Statistics.";
2238 return "Warning, it is dramatically tedious to enforce the boundary items.";
2240 return "Not enough memory at this time, nevertheless, the program continues. "
2241 "The expected mesh will be correct but not really as large as required.";
2243 return "see above error code, resulting quality may be poor.";
2245 return "Not enough memory at this time, nevertheless, the program continues (warning).";
2247 return "Unknown face type.";
2250 return "End of file. Some data are missing in the file.";
2252 return "A too small volume element is detected.";
2254 return "There exists at least a null or negative volume element.";
2256 return "There exist null or negative volume elements.";
2258 return "A too small volume element is detected. A face is considered being degenerated.";
2260 return "Some element is suspected to be very bad shaped or wrong.";
2262 return "A too bad quality face is detected. This face is considered degenerated.";
2264 return "A too bad quality face is detected. This face is degenerated.";
2266 return "Presumably, the surface mesh is not compatible with the domain being processed.";
2268 return "Abnormal error occured, contact hotline.";
2270 return "Not enough memory for the face table.";
2272 return "The algorithm cannot run further. "
2273 "The surface mesh is probably very bad in terms of quality.";
2275 return "Bad vertex number.";
2277 return "Cannot close mesh file NomFil.";
2279 return "There are wrong data.";
2281 return "The number of faces is negative or null.";
2283 return "The number of vertices is negative or null in the '.sol' file.";
2285 return "The number of tetrahedra is negative or null.";
2287 return "The number of vertices is negative or null.";
2289 return "A face has a vertex negative or null.";
2291 return "The field is not a size in file NomFil.";
2293 return "A count is wrong in the enclosing box in the .boite.mesh input "
2294 "file (option '--read_boite').";
2296 return "A tetrahedron has a vertex with a negative number.";
2298 return "the 'MeshVersionFormatted' is not 1 or 2 in the '.mesh' file or the '.sol'.";
2300 return "The number of values in the '.sol' (metric file) is incompatible with "
2301 "the expected value of number of mesh vertices in the '.mesh' file.";
2303 return "Not enough memory.";
2305 return "Not enough memory for the face table.";
2307 return "Insufficient memory ressources detected due to a bad quality "
2308 "surface mesh leading to too many swaps.";
2310 return "The surface coordinates of a vertex are differing from the "
2311 "volume coordinates, probably due to a precision problem.";
2313 return "Invalid dimension. Dimension 3 expected.";
2315 return "A point has a tag 0. This point is probably outside the domain which has been meshed.";
2317 return "The vertices of an element are too close to one another or coincident.";
2319 return "There are at least two points whose distance is very small, and considered as coincident.";
2321 return "Two vertices are too close to one another or coincident.";
2323 return "A vertex cannot be inserted.";
2325 return "Two vertices are too close to one another or coincident. Note : When "
2326 "this error occurs during the overconstrained processing phase, this is only "
2327 "a warning which means that it is difficult to break some overconstrained facets.";
2329 return "Two surface edges are intersecting.";
2331 return "A surface edge intersects a surface face.";
2333 return "A boundary point lies within a surface face.";
2335 return "A boundary point lies within a surface edge.";
2337 return "A surface mesh appears more than once in the input surface mesh.";
2339 return "An edge appears more than once in the input surface mesh.";
2341 return "Surface with unvalid triangles.";
2343 return "The metric in the '.sol' file contains more than one field.";
2345 return "The surface mesh includes at least one hole. The domain is not well defined.";
2347 return "Presumably, the surface mesh is not compatible with the domain being processed (warning).";
2349 return "Probable faces overlapping somewher.";
2351 return "The quadratic version does not work with prescribed free edges.";
2353 return "The quadratic version does not work with a volume mesh.";
2355 return "The metric in the '.sol' file is inadequate (values not per vertices).";
2357 return "The number of vertices in the '.sol' is different from the one in the "
2358 "'.mesh' file for the required vertices (option '--required_vertices').";
2360 return "More than one type in file NomFil. The type must be equal to 1 in the '.sol'"
2361 "for the required vertices (option '--required_vertices').";
2363 return "Bad vertex number.";
2365 return "No guess to start the definition of the connected component(s).";
2367 return "Some initial points cannot be inserted.";
2369 return "A too bad quality face is detected. This face is considered degenerated.";
2371 return "A too bad quality face is detected. This face is degenerated.";
2373 return "The algorithm cannot run further.";
2375 return "A too small volume element is detected.";
2377 return "A tetrahedra is suspected to be very bad shaped or wrong.";
2379 return "There is at least a null or negative volume element. The resulting mesh"
2380 "may be inappropriate.";
2382 return "There are some null or negative volume element. The resulting mesh may"
2383 "be inappropriate.";
2385 return "An edge is unique (i.e., bounds a hole in the surface).";
2387 return "Abnormal or internal error.";
2389 return "Too many components with respect to too many sub-domain.";
2391 return "An internal error has been encountered or a signal has been received. "
2392 "Current mesh will not be saved.";
2394 return "Impossible to define a component.";
2396 return "There are some overconstrained edges.";
2398 return "There are some overconstrained facets.";
2400 return "Give the number of missing faces (information given when regeneration phase failed).";
2402 return "A constrained face cannot be enforced (information given when regeneration phase failed).";
2404 return "A constrained edge cannot be enforced.";
2406 return "It is dramatically tedious to enforce the boundary items.";
2408 return "The surface mesh regeneration step has failed. A .boite.mesh and .boite.map files are created.";
2410 return "Invalid resulting mesh.";
2412 return "P2 correction not successful.";
2414 return "Program has received an interruption or a termination signal sent by the "
2415 "user or the system administrator. Current mesh will not be saved.";
2420 //================================================================================
2422 * \brief Retrieve from a string given number of integers
2424 //================================================================================
2426 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
2429 ids.reserve( nbIds );
2432 while ( !isdigit( *ptr )) ++ptr;
2433 if ( ptr[-1] == '-' ) --ptr;
2434 ids.push_back((int) strtol( ptr, &ptr, 10 ));
2440 //================================================================================
2442 * \brief Retrieve problem description form a log file
2443 * \retval bool - always false
2445 //================================================================================
2447 SMESH_ComputeErrorPtr
2448 GHS3DPlugin_GHS3D::getErrorDescription(const char* logFile,
2449 const std::string& log,
2450 const _Ghs2smdsConvertor & toSmdsConvertor,
2451 const bool isOk/* = false*/ )
2453 SMESH_BadInputElements* badElemsErr =
2454 new SMESH_BadInputElements( toSmdsConvertor.getMesh(), COMPERR_ALGO_FAILED );
2455 SMESH_ComputeErrorPtr err( badElemsErr );
2457 char* ptr = const_cast<char*>( log.c_str() );
2458 char* buf = ptr, * bufEnd = ptr + log.size();
2461 SMESH_Comment errDescription;
2463 enum { NODE = 1, EDGE, TRIA, VOL, SKIP_ID = 1 };
2465 // look for MeshGems version
2466 // Since "MG-TETRA -- MeshGems 1.1-3 (January, 2013)" error codes change.
2467 // To discriminate old codes from new ones we add 1000000 to the new codes.
2468 // This way value of the new codes is same as absolute value of codes printed
2469 // in the log after "MGMESSAGE" string.
2470 int versionAddition = 0;
2473 while ( ++verPtr < bufEnd )
2475 if ( strncmp( verPtr, "MG-TETRA -- MeshGems ", 21 ) != 0 )
2477 if ( strcmp( verPtr, "MG-TETRA -- MeshGems 1.1-3 " ) >= 0 )
2478 versionAddition = 1000000;
2484 // look for errors "ERR #"
2486 set<string> foundErrorStr; // to avoid reporting same error several times
2487 set<int> elemErrorNums; // not to report different types of errors with bad elements
2488 while ( ++ptr < bufEnd )
2490 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
2493 list<const SMDS_MeshElement*>& badElems = badElemsErr->myBadElements;
2494 vector<int> nodeIds;
2498 int errNum = int( strtol(ptr, &ptr, 10) + versionAddition );
2499 // we treat errors enumerated in [SALOME platform 0019316] issue
2500 // and all errors from a new (Release 1.1) MeshGems User Manual
2502 case 0015: // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
2503 case 1005620 : // a too bad quality face is detected. This face is considered degenerated.
2504 ptr = getIds(ptr, SKIP_ID, nodeIds);
2505 ptr = getIds(ptr, TRIA, nodeIds);
2506 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2508 case 1005621 : // a too bad quality face is detected. This face is degenerated.
2509 // hence the is degenerated it is invisible, add its edges in addition
2510 ptr = getIds(ptr, SKIP_ID, nodeIds);
2511 ptr = getIds(ptr, TRIA, nodeIds);
2512 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2514 vector<int> edgeNodes( nodeIds.begin(), --nodeIds.end() ); // 01
2515 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2516 edgeNodes[1] = nodeIds[2]; // 02
2517 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2518 edgeNodes[0] = nodeIds[1]; // 12
2521 case 1000: // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
2523 case 1002: // Face (f 1, f 2, f 3) has a vertex negative or null
2524 case 3019: // Constrained face (f 1, f 2, f 3) cannot be enforced
2525 case 1002211: // a face has a vertex negative or null.
2526 case 1005200 : // a surface mesh appears more than once in the input surface mesh.
2527 case 1008423 : // a constrained face cannot be enforced (regeneration phase failed).
2528 ptr = getIds(ptr, TRIA, nodeIds);
2529 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2531 case 1001: // Edge (e1, e2) appears more than once in the input surface mesh
2532 case 3009: // Constrained edge (e1, e2) cannot be enforced (warning).
2533 // ERR 3109 : EDGE 5 6 UNIQUE
2534 case 3109: // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
2535 case 1005210 : // an edge appears more than once in the input surface mesh.
2536 case 1005820 : // an edge is unique (i.e., bounds a hole in the surface).
2537 case 1008441 : // a constrained edge cannot be enforced.
2538 ptr = getIds(ptr, EDGE, nodeIds);
2539 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2541 case 2004: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
2542 case 2014: // at least two points whose distance is dist, i.e., considered as coincident
2543 case 2103: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
2544 // ERR 2103 : 16 WITH 3
2545 case 1005105 : // two vertices are too close to one another or coincident.
2546 case 1005107: // Two vertices are too close to one another or coincident.
2547 ptr = getIds(ptr, NODE, nodeIds);
2548 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2549 ptr = getIds(ptr, NODE, nodeIds);
2550 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2552 case 2012: // Vertex v1 cannot be inserted (warning).
2553 case 1005106 : // a vertex cannot be inserted.
2554 ptr = getIds(ptr, NODE, nodeIds);
2555 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2557 case 3103: // The surface edge (e1, e2) intersects another surface edge (e3, e4)
2558 case 1005110 : // two surface edges are intersecting.
2559 // ERR 3103 : 1 2 WITH 7 3
2560 ptr = getIds(ptr, EDGE, nodeIds);
2561 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2562 ptr = getIds(ptr, EDGE, nodeIds);
2563 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2565 case 3104: // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
2566 // ERR 3104 : 9 10 WITH 1 2 3
2567 case 3106: // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
2568 case 1005120 : // a surface edge intersects a surface face.
2569 ptr = getIds(ptr, EDGE, nodeIds);
2570 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2571 ptr = getIds(ptr, TRIA, nodeIds);
2572 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2574 case 3105: // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
2575 // ERR 3105 : 8 IN 2 3 5
2576 case 1005150 : // a boundary point lies within a surface face.
2577 ptr = getIds(ptr, NODE, nodeIds);
2578 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2579 ptr = getIds(ptr, TRIA, nodeIds);
2580 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2582 case 3107: // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
2583 // ERR 3107 : 2 IN 4 1
2584 case 1005160 : // a boundary point lies within a surface edge.
2585 ptr = getIds(ptr, NODE, nodeIds);
2586 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2587 ptr = getIds(ptr, EDGE, nodeIds);
2588 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2590 case 9000: // ERR 9000
2591 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
2592 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
2593 // A too small volume element is detected. Are reported the index of the element,
2594 // its four vertex indices, its volume and the tolerance threshold value
2595 ptr = getIds(ptr, SKIP_ID, nodeIds);
2596 ptr = getIds(ptr, VOL, nodeIds);
2597 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2598 // even if all nodes found, volume it most probably invisible,
2599 // add its faces to demonstrate it anyhow
2601 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
2602 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2603 faceNodes[2] = nodeIds[3]; // 013
2604 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2605 faceNodes[1] = nodeIds[2]; // 023
2606 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2607 faceNodes[0] = nodeIds[1]; // 123
2608 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2611 case 9001: // ERR 9001
2612 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
2613 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
2614 // %% NUMBER OF NULL VOLUME TETS : 0
2615 // There exists at least a null or negative volume element
2618 // There exist n null or negative volume elements
2621 // A too small volume element is detected
2624 // A too bad quality face is detected. This face is considered degenerated,
2625 // its index, its three vertex indices together with its quality value are reported
2626 break; // same as next
2627 case 9112: // ERR 9112
2628 // FACE 2 WITH VERTICES : 4 2 5
2629 // SMALL INRADIUS : 0.
2630 // A too bad quality face is detected. This face is degenerated,
2631 // its index, its three vertex indices together with its inradius are reported
2632 ptr = getIds(ptr, SKIP_ID, nodeIds);
2633 ptr = getIds(ptr, TRIA, nodeIds);
2634 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2635 // add triangle edges as it most probably has zero area and hence invisible
2637 vector<int> edgeNodes(2);
2638 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
2639 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2640 edgeNodes[1] = nodeIds[2]; // 0-2
2641 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2642 edgeNodes[0] = nodeIds[1]; // 1-2
2643 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2646 case 1005103 : // the vertices of an element are too close to one another or coincident.
2647 ptr = getIds(ptr, TRIA, nodeIds);
2648 if ( nodeIds.back() == 0 ) // index of the third vertex of the element (0 for an edge)
2649 nodeIds.resize( EDGE );
2650 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2654 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
2656 continue; // not to report same error several times
2658 // const SMDS_MeshElement* nullElem = 0;
2659 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
2661 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
2662 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
2663 // if ( oneMoreErrorType )
2664 // continue; // not to report different types of errors with bad elements
2668 string text = translateError( errNum );
2669 if ( errDescription.find( text ) == text.npos ) {
2670 if ( !errDescription.empty() )
2671 errDescription << "\n";
2672 errDescription << text;
2677 if ( errDescription.empty() ) { // no errors found
2678 char msgLic1[] = "connection to server failed";
2679 char msgLic2[] = " Dlim ";
2680 char msgLic3[] = "license is not valid";
2681 if ( search( &buf[0], bufEnd, msgLic1, msgLic1 + strlen(msgLic1)) != bufEnd ||
2682 search( &buf[0], bufEnd, msgLic2, msgLic2 + strlen(msgLic2)) != bufEnd )
2683 errDescription << "Network license problem.";
2684 else if ( search( &buf[0], bufEnd, msgLic3, msgLic3 + strlen(msgLic3)) != bufEnd )
2685 errDescription << "License is not valid.";
2688 char msg2[] = "SEGMENTATION FAULT";
2689 if ( search( &buf[0], bufEnd, msg2, msg2 + strlen(msg2)) != bufEnd )
2690 errDescription << "MG-Tetra: SEGMENTATION FAULT. ";
2694 if ( !isOk && logFile && logFile[0] )
2696 if ( errDescription.empty() )
2697 errDescription << "See " << logFile << " for problem description";
2699 errDescription << "\nSee " << logFile << " for more information";
2702 err->myComment = errDescription;
2704 if ( err->myComment.empty() && !err->HasBadElems() )
2705 err = SMESH_ComputeError::New(); // OK
2710 //================================================================================
2712 * \brief Creates _Ghs2smdsConvertor
2714 //================================================================================
2716 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap,
2717 SMESH_ProxyMesh::Ptr mesh)
2718 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 ), _mesh( mesh )
2722 //================================================================================
2724 * \brief Creates _Ghs2smdsConvertor
2726 //================================================================================
2728 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId,
2729 SMESH_ProxyMesh::Ptr mesh)
2730 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId ), _mesh( mesh )
2734 //================================================================================
2736 * \brief Return SMDS element by ids of MG-Tetra nodes
2738 //================================================================================
2740 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
2742 size_t nbNodes = ghsNodes.size();
2743 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
2744 for ( size_t i = 0; i < nbNodes; ++i ) {
2745 int ghsNode = ghsNodes[ i ];
2746 if ( _ghs2NodeMap ) {
2747 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
2748 if ( in == _ghs2NodeMap->end() )
2750 nodes[ i ] = in->second;
2753 if ( ghsNode < 1 || ghsNode > (int)_nodeByGhsId->size() )
2755 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
2761 if ( nbNodes == 2 ) {
2762 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
2763 if ( !edge || edge->GetID() < 1 || _mesh->IsTemporary( edge ))
2764 edge = new SMDS_LinearEdge( nodes[0], nodes[1] );
2767 if ( nbNodes == 3 ) {
2768 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
2769 if ( !face || face->GetID() < 1 || _mesh->IsTemporary( face ))
2770 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
2774 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
2779 //================================================================================
2781 * \brief Return a mesh
2783 //================================================================================
2785 const SMDS_Mesh* _Ghs2smdsConvertor::getMesh() const
2787 return _mesh->GetMeshDS();
2790 //=============================================================================
2794 //=============================================================================
2795 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
2796 const TopoDS_Shape& aShape,
2797 MapShapeNbElems& aResMap)
2799 smIdType nbtri = 0, nbqua = 0;
2800 double fullArea = 0.0;
2801 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2802 TopoDS_Face F = TopoDS::Face( exp.Current() );
2803 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2804 MapShapeNbElemsItr anIt = aResMap.find(sm);
2805 if( anIt==aResMap.end() ) {
2806 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2807 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
2808 "Submesh can not be evaluated",this));
2811 std::vector<smIdType> aVec = (*anIt).second;
2812 nbtri += std::max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
2813 nbqua += std::max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
2815 BRepGProp::SurfaceProperties(F,G);
2816 double anArea = G.Mass();
2820 // collect info from edges
2821 smIdType nb0d_e = 0, nb1d_e = 0;
2822 bool IsQuadratic = false;
2823 bool IsFirst = true;
2824 TopTools_MapOfShape tmpMap;
2825 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2826 TopoDS_Edge E = TopoDS::Edge(exp.Current());
2827 if( tmpMap.Contains(E) )
2830 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
2831 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
2832 std::vector<smIdType> aVec = (*anIt).second;
2833 nb0d_e += aVec[SMDSEntity_Node];
2834 nb1d_e += std::max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
2836 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
2842 double ELen = sqrt(2.* ( fullArea/double(nbtri+nbqua*2) ) / sqrt(3.0) );
2845 BRepGProp::VolumeProperties(aShape,G);
2846 double aVolume = G.Mass();
2847 double tetrVol = 0.1179*ELen*ELen*ELen;
2848 double CoeffQuality = 0.9;
2849 smIdType nbVols = smIdType(aVolume/tetrVol/CoeffQuality);
2850 smIdType nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2851 smIdType nb1d_in = (smIdType) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2852 std::vector<smIdType> aVec(SMDSEntity_Last);
2853 for(smIdType i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2855 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2856 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2857 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2860 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2861 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2862 aVec[SMDSEntity_Pyramid] = nbqua;
2864 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2865 aResMap.insert(std::make_pair(sm,aVec));
2870 bool GHS3DPlugin_GHS3D::importGMFMesh(const char* theGMFFileName, SMESH_Mesh& theMesh)
2872 SMESH_ComputeErrorPtr err = theMesh.GMFToMesh( theGMFFileName, /*makeRequiredGroups =*/ true );
2874 theMesh.GetMeshDS()->Modified();
2876 return ( !err || err->IsOK());
2881 //================================================================================
2883 * \brief Sub-mesh event listener setting enforced elements as soon as an enforced
2886 struct _EnforcedMeshRestorer : public SMESH_subMeshEventListener
2888 _EnforcedMeshRestorer():
2889 SMESH_subMeshEventListener( /*isDeletable = */true, Name() )
2892 //================================================================================
2894 * \brief Returns an ID of listener
2896 static const char* Name() { return "GHS3DPlugin_GHS3D::_EnforcedMeshRestorer"; }
2898 //================================================================================
2900 * \brief Treat events of the subMesh
2902 void ProcessEvent(const int event,
2903 const int eventType,
2904 SMESH_subMesh* /*subMesh*/,
2905 SMESH_subMeshEventListenerData* data,
2906 const SMESH_Hypothesis* /*hyp*/)
2908 if ( SMESH_subMesh::SUBMESH_LOADED == event &&
2909 SMESH_subMesh::COMPUTE_EVENT == eventType &&
2911 !data->mySubMeshes.empty() )
2913 // An enforced mesh (subMesh->_father) has been loaded from hdf file
2914 if ( GHS3DPlugin_Hypothesis* hyp = GetGHSHypothesis( data->mySubMeshes.front() ))
2915 hyp->RestoreEnfElemsByMeshes();
2918 //================================================================================
2920 * \brief Returns GHS3DPlugin_Hypothesis used to compute a subMesh
2922 static GHS3DPlugin_Hypothesis* GetGHSHypothesis( SMESH_subMesh* subMesh )
2924 SMESH_HypoFilter ghsHypFilter
2925 ( SMESH_HypoFilter::HasName( GHS3DPlugin_Hypothesis::GetHypType() ));
2926 return (GHS3DPlugin_Hypothesis* )
2927 subMesh->GetFather()->GetHypothesis( subMesh->GetSubShape(),
2929 /*visitAncestors=*/true);
2933 //================================================================================
2935 * \brief Sub-mesh event listener removing empty groups created due to "To make
2936 * groups of domains".
2938 struct _GroupsOfDomainsRemover : public SMESH_subMeshEventListener
2940 _GroupsOfDomainsRemover():
2941 SMESH_subMeshEventListener( /*isDeletable = */true,
2942 "GHS3DPlugin_GHS3D::_GroupsOfDomainsRemover" ) {}
2944 * \brief Treat events of the subMesh
2946 void ProcessEvent(const int /*event*/,
2947 const int eventType,
2948 SMESH_subMesh* subMesh,
2949 SMESH_subMeshEventListenerData* /*data*/,
2950 const SMESH_Hypothesis* /*hyp*/)
2952 if (SMESH_subMesh::ALGO_EVENT == eventType &&
2953 !subMesh->GetAlgo() )
2955 removeEmptyGroupsOfDomains( subMesh->GetFather(), /*notEmptyAsWell=*/true );
2961 //================================================================================
2963 * \brief Set an event listener to set enforced elements as soon as an enforced
2966 //================================================================================
2968 void GHS3DPlugin_GHS3D::SubmeshRestored(SMESH_subMesh* subMesh)
2970 if ( GHS3DPlugin_Hypothesis* hyp = _EnforcedMeshRestorer::GetGHSHypothesis( subMesh ))
2972 GHS3DPlugin_Hypothesis::TGHS3DEnforcedMeshList enfMeshes = hyp->_GetEnforcedMeshes();
2973 GHS3DPlugin_Hypothesis::TGHS3DEnforcedMeshList::iterator it = enfMeshes.begin();
2974 for(;it != enfMeshes.end();++it) {
2975 GHS3DPlugin_Hypothesis::TGHS3DEnforcedMesh* enfMesh = *it;
2976 if ( SMESH_Mesh* mesh = GetMeshByPersistentID( enfMesh->persistID ))
2978 SMESH_subMesh* smToListen = mesh->GetSubMesh( mesh->GetShapeToMesh() );
2979 // a listener set to smToListen will care of hypothesis stored in SMESH_EventListenerData
2980 subMesh->SetEventListener( new _EnforcedMeshRestorer(),
2981 SMESH_subMeshEventListenerData::MakeData( subMesh ),
2988 //================================================================================
2990 * \brief Sets an event listener removing empty groups created due to "To make
2991 * groups of domains".
2992 * \param subMesh - submesh where algo is set
2994 * This method is called when a submesh gets HYP_OK algo_state.
2995 * After being set, event listener is notified on each event of a submesh.
2997 //================================================================================
2999 void GHS3DPlugin_GHS3D::SetEventListener(SMESH_subMesh* subMesh)
3001 subMesh->SetEventListener( new _GroupsOfDomainsRemover(), 0, subMesh );
3004 //================================================================================
3006 * \brief If possible, returns progress of computation [0.,1.]
3008 //================================================================================
3010 double GHS3DPlugin_GHS3D::GetProgress() const
3014 // this->_progress is advanced by MG_Tetra_API according to messages from MG library
3015 // but sharply. Advance it a bit to get smoother advancement.
3016 GHS3DPlugin_GHS3D* me = const_cast<GHS3DPlugin_GHS3D*>( this );
3017 if ( _progress < 0.1 ) // the first message is at 10%
3018 me->_progress = GetProgressByTic();
3019 else if ( _progress < 0.98 )
3020 me->_progress += _progressAdvance;