1 // Copyright (C) 2004-2010 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.
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"
31 #include <Basics_Utils.hxx>
33 #include "SMESH_Gen.hxx"
34 #include "SMESH_Mesh.hxx"
35 #include "SMESH_Comment.hxx"
36 #include "SMESH_MesherHelper.hxx"
37 #include "SMESH_MeshEditor.hxx"
39 #include "SMDS_MeshElement.hxx"
40 #include "SMDS_MeshNode.hxx"
41 #include "SMDS_FaceOfNodes.hxx"
42 #include "SMDS_VolumeOfNodes.hxx"
44 #include <StdMeshers_QuadToTriaAdaptor.hxx>
46 #include <BRepAdaptor_Surface.hxx>
47 #include <BRepBndLib.hxx>
48 #include <BRepClass3d_SolidClassifier.hxx>
49 #include <BRepTools.hxx>
50 #include <BRep_Tool.hxx>
51 #include <Bnd_Box.hxx>
52 #include <GeomAPI_ProjectPointOnSurf.hxx>
53 #include <OSD_File.hxx>
54 #include <Precision.hxx>
55 #include <Quantity_Parameter.hxx>
56 #include <Standard_ProgramError.hxx>
57 #include <Standard_ErrorHandler.hxx>
58 #include <Standard_Failure.hxx>
60 #include <TopExp_Explorer.hxx>
61 #include <TopTools_IndexedMapOfShape.hxx>
62 #include <TopTools_ListIteratorOfListOfShape.hxx>
64 //#include <BRepClass_FaceClassifier.hxx>
65 #include <TopTools_MapOfShape.hxx>
66 #include <BRepGProp.hxx>
67 #include <GProp_GProps.hxx>
69 #include "utilities.h"
74 #include <sys/sysinfo.h>
79 //#include <Standard_Stream.hxx>
82 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
103 typedef const list<const SMDS_MeshFace*> TTriaList;
105 static void removeFile( const TCollection_AsciiString& fileName )
108 OSD_File( fileName ).Remove();
110 catch ( Standard_ProgramError ) {
111 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
115 //=============================================================================
119 //=============================================================================
121 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
122 : SMESH_3D_Algo(hypId, studyId, gen)
124 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
126 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
127 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
130 _compatibleHypothesis.push_back("GHS3D_Parameters");
131 _requireShape = false; // can work without shape
134 //=============================================================================
138 //=============================================================================
140 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
142 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
145 //=============================================================================
149 //=============================================================================
151 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
152 const TopoDS_Shape& aShape,
153 Hypothesis_Status& aStatus )
155 aStatus = SMESH_Hypothesis::HYP_OK;
157 // there is only one compatible Hypothesis so far
160 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
162 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
164 _keepFiles = _hyp->GetKeepFiles();
169 //=======================================================================
170 //function : findShape
172 //=======================================================================
174 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
176 const TopoDS_Shape shape[],
179 TopAbs_State * state = 0)
182 int j, iShape, nbNode = 4;
184 for ( j=0; j<nbNode; j++ ) {
185 gp_XYZ p ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
186 if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) {
193 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
194 if (state) *state = SC.State();
195 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
196 for (iShape = 0; iShape < nShape; iShape++) {
197 aShape = shape[iShape];
198 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
199 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
200 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
201 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
202 if (state) *state = SC.State();
203 if (SC.State() == TopAbs_IN)
211 //=======================================================================
212 //function : readMapIntLine
214 //=======================================================================
216 static char* readMapIntLine(char* ptr, int tab[]) {
218 std::cout << std::endl;
220 for ( int i=0; i<17; i++ ) {
221 intVal = strtol(ptr, &ptr, 10);
228 //=======================================================================
229 //function : countShape
231 //=======================================================================
233 template < class Mesh, class Shape >
234 static int countShape( Mesh* mesh, Shape shape ) {
235 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
237 for ( ; expShape.More(); expShape.Next() ) {
243 //=======================================================================
244 //function : writeFaces
246 //=======================================================================
248 static bool writeFaces (ofstream & theFile,
249 SMESHDS_Mesh* theMesh,
250 const TopoDS_Shape& theShape,
251 vector<StdMeshers_QuadToTriaAdaptor>& theQuad2Trias,
252 const map <int,int> & theSmdsToGhs3dIdMap)
256 // NB_ELEMS DUMMY_INT
257 // Loop from 1 to NB_ELEMS
258 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
261 SMESHDS_SubMesh* theSubMesh;
262 const SMDS_MeshElement* aFace;
263 const char* space = " ";
264 const int dummyint = 0;
265 map<int,int>::const_iterator itOnMap;
266 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
267 int nbNodes, aSmdsID;
269 // count triangles bound to geometry
272 TopTools_IndexedMapOfShape facesMap;
273 TopExp::MapShapes( theShape, TopAbs_FACE, facesMap );
275 // 2 adaptors for each face in facesMap, as a face can belong to 2 solids
276 typedef vector< StdMeshers_QuadToTriaAdaptor* > TwoAdaptors;
277 vector< TwoAdaptors > qttaByFace;
278 if ( theQuad2Trias.empty() )
280 // case w/o quadrangles
281 for ( int i = 1; i <= facesMap.Extent(); ++i )
282 if (( theSubMesh = theMesh->MeshElements( facesMap(i))))
283 nbTriangles += theSubMesh->NbElements();
287 // case with quadrangles
288 qttaByFace.resize( facesMap.Extent() );
289 for ( unsigned i = 0; i < theQuad2Trias.size(); ++i )
291 TopoDS_Shape solid = theQuad2Trias[i].GetShape();
292 TopExp_Explorer expface( solid, TopAbs_FACE );
293 for ( ; expface.More(); expface.Next() )
294 if (( theSubMesh = theMesh->MeshElements( expface.Current()) ))
296 const int faceIndex = facesMap.Add( expface.Current() );
297 TwoAdaptors& aTwoAdaptors = qttaByFace[ faceIndex-1 ];
298 const bool newFaceEncounters = aTwoAdaptors.empty();
299 aTwoAdaptors.push_back( & theQuad2Trias[i] );
301 // on a shared face encountered for the second time
302 // we count only triangles of pyramids
303 const int countTrias = int( newFaceEncounters );
304 itOnSubMesh = theSubMesh->GetElements();
305 while ( itOnSubMesh->more() )
307 aFace = itOnSubMesh->next();
308 if ( aFace->NbCornerNodes() != 4 )
309 nbTriangles += countTrias;
310 else if ( TTriaList* trias = theQuad2Trias[i].GetTriangles( aFace ))
311 nbTriangles += trias->size();
313 nbTriangles += countTrias;
319 std::cout << " " << facesMap.Extent() << " shapes of 2D dimension" << std::endl;
320 std::cout << std::endl;
322 theFile << space << nbTriangles << space << dummyint << std::endl;
324 vector< const SMDS_MeshElement* > trias;
325 trias.resize( 8 ); // 4 triangles from 2 pyramids basing on one quadranle
327 for ( int i = 1; i <= facesMap.Extent(); i++ )
329 aShape = facesMap(i);
330 theSubMesh = theMesh->MeshElements( aShape );
331 if ( !theSubMesh ) continue;
332 TwoAdaptors& aTwoAdaptors = qttaByFace[ i-1 ];
333 itOnSubMesh = theSubMesh->GetElements();
334 while ( itOnSubMesh->more() )
336 aFace = itOnSubMesh->next();
337 if ( aFace->NbCornerNodes() == 4 )
340 for ( unsigned j = 0; j < aTwoAdaptors.size(); ++j )
341 if ( TTriaList* t = aTwoAdaptors[j]->GetTriangles( aFace ))
342 for ( TTriaList::const_iterator tIt = t->begin(); tIt != t->end(); ++tIt)
343 trias[nbTriangles++] = *tIt;
344 if ( nbTriangles == 0 )
355 for ( int j = 0; j < nbTriangles; ++j )
358 nbNodes = aFace->NbNodes();
360 theFile << space << nbNodes;
362 itOnSubFace = aFace->nodesIterator();
363 while ( itOnSubFace->more() ) {
365 aSmdsID = itOnSubFace->next()->GetID();
366 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
367 // if ( itOnMap == theSmdsToGhs3dIdMap.end() ) {
368 // cout << "not found node: " << aSmdsID << endl;
371 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
373 theFile << space << (*itOnMap).second;
376 // (NB_NODES + 1) times: DUMMY_INT
377 for ( int j=0; j<=nbNodes; j++)
378 theFile << space << dummyint;
380 theFile << std::endl;
388 //=======================================================================
389 //function : writeFaces
390 //purpose : Write Faces in case if generate 3D mesh w/o geometry
391 //=======================================================================
393 static bool writeFaces (ofstream & theFile,
394 SMESHDS_Mesh * theMesh,
395 StdMeshers_QuadToTriaAdaptor& theQuad2Trias,
396 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
400 // NB_ELEMS DUMMY_INT
401 // Loop from 1 to NB_ELEMS
402 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
404 int nbNodes, nbTriangles = 0;
405 map< const SMDS_MeshNode*,int >::iterator it;
406 SMDS_ElemIteratorPtr nodeIt;
407 const SMDS_MeshElement* elem;
409 const char* space = " ";
410 const int dummyint = 0;
414 nbTriangles = theQuad2Trias.TotalNbOfTriangles();
415 if ( nbTriangles == 0 )
416 // theQuad2Trias not computed as there are no quadrangles in mesh
417 nbTriangles = theMesh->NbFaces();
419 if ( nbTriangles == 0 )
422 std::cout << "The initial 2D mesh contains " << nbTriangles << " faces and ";
424 // NB_ELEMS DUMMY_INT
425 theFile << space << nbTriangles << space << dummyint << std::endl;
428 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
430 vector< const SMDS_MeshElement* > trias;
431 trias.resize( 8 ); // 8 == 4 triangles from 2 pyramids basing on one quadranle
433 // Loop from 1 to NB_ELEMS
435 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
436 while ( eIt->more() )
440 if ( elem->NbCornerNodes() == 4 )
443 if ( TTriaList* t = theQuad2Trias.GetTriangles( elem ))
444 for ( TTriaList::const_iterator tIt = t->begin(); tIt != t->end(); ++tIt)
445 trias[nbTriangles++] = *tIt;
446 if ( nbTriangles == 0 )
458 for ( int j = 0; j < nbTriangles; ++j )
462 nbNodes = elem->NbNodes();
463 theFile << space << nbNodes;
465 // NODE_NB_1 NODE_NB_2 ...
466 nodeIt = elem->nodesIterator();
467 while ( nodeIt->more() )
470 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
471 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
472 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
473 theFile << space << it->second;
476 // (NB_NODES + 1) times: DUMMY_INT
477 for ( int i=0; i<=nbNodes; i++)
478 theFile << space << dummyint;
479 theFile << std::endl;
483 // put nodes to theNodeByGhs3dId vector
484 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
485 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
486 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
488 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
494 //=======================================================================
495 //function : writePoints
497 //=======================================================================
499 static bool writePoints (ofstream & theFile,
500 SMESH_MesherHelper& theHelper,
501 map <int,int> & theSmdsToGhs3dIdMap,
502 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap,
503 map<vector<double>,double> & theEnforcedVertices)
508 // Loop from 1 to NB_NODES
511 SMESHDS_Mesh * theMesh = theHelper.GetMeshDS();
512 int nbNodes = theMesh->NbNodes();
515 int nbEnforcedVertices = theEnforcedVertices.size();
517 // Issue 020674: EDF 870 SMESH: Mesh generated by Netgen not usable by GHS3D
518 // The problem is in nodes on degenerated edges, we need to skip them
519 if ( theHelper.HasDegeneratedEdges() )
521 // here we decrease total nb of nodes by nb of nodes on degenerated edges
523 for (TopExp_Explorer e(theMesh->ShapeToMesh(), TopAbs_EDGE ); e.More(); e.Next())
525 SMESH_subMesh* sm = theHelper.GetMesh()->GetSubMesh( e.Current() );
526 if ( checkedSM.insert( sm->GetId() ).second && theHelper.IsDegenShape(sm->GetId() )) {
527 if ( sm->GetSubMeshDS() )
528 nbNodes -= sm->GetSubMeshDS()->NbNodes();
532 const char* space = " ";
533 const int dummyint = 0;
536 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
537 const SMDS_MeshNode* node;
540 std::cout << std::endl;
541 std::cout << "The initial 2D mesh contains :" << std::endl;
542 std::cout << " " << nbNodes << " nodes" << std::endl;
543 if (nbEnforcedVertices > 0) {
544 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
546 std::cout << std::endl;
547 std::cout << "Start writing in 'points' file ..." << std::endl;
548 theFile << space << nbNodes << std::endl;
550 // Loop from 1 to NB_NODES
555 if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE &&
556 theHelper.IsDegenShape( node->GetPosition()->GetShapeId() )) // Issue 020674
559 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
560 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
565 << space << node->X()
566 << space << node->Y()
567 << space << node->Z()
568 << space << dummyint;
570 theFile << std::endl;
574 // Iterate over the enforced vertices
575 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
576 const TopoDS_Shape shapeToMesh = theMesh->ShapeToMesh();
577 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
578 double x = vertexIt->first[0];
579 double y = vertexIt->first[1];
580 double z = vertexIt->first[2];
581 // Test if point is inside shape to mesh
582 gp_Pnt myPoint(x,y,z);
583 BRepClass3d_SolidClassifier scl(shapeToMesh);
584 scl.Perform(myPoint, 1e-7);
585 TopAbs_State result = scl.State();
586 if ( result == TopAbs_IN ) {
587 MESSAGE("Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second);
588 // X Y Z PHY_SIZE DUMMY_INT
593 << space << vertexIt->second
594 << space << dummyint;
596 theFile << std::endl;
599 MESSAGE("Enforced vertex (" << x << "," << y <<"," << z << ") is not inside the geometry: it was not added ");
603 std::cout << std::endl;
604 std::cout << "End writing in 'points' file." << std::endl;
609 //=======================================================================
610 //function : writePoints
612 //=======================================================================
614 static bool writePoints (ofstream & theFile,
615 SMESH_Mesh * theMesh,
616 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId,
617 map<vector<double>,double> & theEnforcedVertices)
622 // Loop from 1 to NB_NODES
625 int nbNodes = theNodeByGhs3dId.size();
629 int nbEnforcedVertices = theEnforcedVertices.size();
631 const char* space = " ";
632 const int dummyint = 0;
634 const SMDS_MeshNode* node;
637 std::cout << std::endl;
638 std::cout << "The initial 2D mesh contains :" << std::endl;
639 std::cout << " " << nbNodes << " nodes" << std::endl;
640 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
641 std::cout << std::endl;
642 std::cout << "Start writing in 'points' file ..." << std::endl;
643 theFile << space << nbNodes << std::endl;
645 // Loop from 1 to NB_NODES
647 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
648 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
649 for ( ; nodeIt != after; ++nodeIt )
655 << space << node->X()
656 << space << node->Y()
657 << space << node->Z()
658 << space << dummyint;
660 theFile << std::endl;
664 // Iterate over the enforced vertices
665 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
666 auto_ptr< SMESH_ElementSearcher > pntCls ( SMESH_MeshEditor( theMesh ).GetElementSearcher());
667 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
668 double x = vertexIt->first[0];
669 double y = vertexIt->first[1];
670 double z = vertexIt->first[2];
671 // Test if point is inside shape to mesh
672 gp_Pnt myPoint(x,y,z);
673 TopAbs_State result = pntCls->GetPointState( myPoint );
674 if ( result == TopAbs_IN ) {
675 std::cout << "Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second << std::endl;
677 // X Y Z PHY_SIZE DUMMY_INT
682 << space << vertexIt->second
683 << space << dummyint;
685 theFile << std::endl;
688 std::cout << std::endl;
689 std::cout << "End writing in 'points' file." << std::endl;
694 //=======================================================================
695 //function : findShapeID
696 //purpose : find the solid corresponding to GHS3D sub-domain following
697 // the technique proposed in GHS3D manual in chapter
698 // "B.4 Subdomain (sub-region) assignment"
699 //=======================================================================
701 static int findShapeID(SMESH_Mesh& mesh,
702 const SMDS_MeshNode* node1,
703 const SMDS_MeshNode* node2,
704 const SMDS_MeshNode* node3,
705 const bool toMeshHoles)
707 const int invalidID = 0;
708 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
710 // face the nodes belong to
711 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
715 std::cout << "bnd face " << face->GetID() << " - ";
717 // geom face the face assigned to
718 SMESH_MeshEditor editor(&mesh);
719 int geomFaceID = editor.FindShape( face );
722 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
723 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
725 TopoDS_Face geomFace = TopoDS::Face( shape );
727 // solids bounded by geom face
728 TopTools_IndexedMapOfShape solids, shells;
729 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
730 for ( ; ansIt.More(); ansIt.Next() ) {
731 switch ( ansIt.Value().ShapeType() ) {
733 solids.Add( ansIt.Value() ); break;
735 shells.Add( ansIt.Value() ); break;
739 // analyse found solids
740 if ( solids.Extent() == 0 || shells.Extent() == 0)
743 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
744 if ( solids.Extent() == 1 )
747 return meshDS->ShapeToIndex( solid1 );
749 // - Are we at a hole boundary face?
750 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
751 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
753 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
754 // check if any edge of shells(1) belongs to another shell
755 for ( ; eExp.More() && !touch; eExp.Next() ) {
756 ansIt = mesh.GetAncestors( eExp.Current() );
757 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
758 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
759 touch = ( !ansIt.Value().IsSame( shells(1) ));
763 return meshDS->ShapeToIndex( solid1 );
766 // find orientation of geom face within the first solid
767 TopExp_Explorer fExp( solid1, TopAbs_FACE );
768 for ( ; fExp.More(); fExp.Next() )
769 if ( geomFace.IsSame( fExp.Current() )) {
770 geomFace = TopoDS::Face( fExp.Current() );
774 return invalidID; // face not found
776 // normale to triangle
777 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
778 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
779 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
780 gp_Vec vec12( node1Pnt, node2Pnt );
781 gp_Vec vec13( node1Pnt, node3Pnt );
782 gp_Vec meshNormal = vec12 ^ vec13;
783 if ( meshNormal.SquareMagnitude() < DBL_MIN )
786 // get normale to geomFace at any node
787 bool geomNormalOK = false;
789 const SMDS_MeshNode* nodes[3] = { node1, node2, node3 };
790 SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
791 for ( int i = 0; !geomNormalOK && i < 3; ++i )
793 // find UV of i-th node on geomFace
794 const SMDS_MeshNode* nNotOnSeamEdge = 0;
795 if ( helper.IsSeamShape( nodes[i]->GetPosition()->GetShapeId() ))
796 if ( helper.IsSeamShape( nodes[(i+1)%3]->GetPosition()->GetShapeId() ))
797 nNotOnSeamEdge = nodes[(i+2)%3];
799 nNotOnSeamEdge = nodes[(i+1)%3];
801 gp_XY uv = helper.GetNodeUV( geomFace, nodes[i], nNotOnSeamEdge, &uvOK );
802 // check that uv is correct
805 TopoDS_Shape nodeShape = helper.GetSubShapeByNode( nodes[i], meshDS );
806 if ( !nodeShape.IsNull() )
807 switch ( nodeShape.ShapeType() )
809 case TopAbs_FACE: tol = BRep_Tool::Tolerance( TopoDS::Face( nodeShape )); break;
810 case TopAbs_EDGE: tol = BRep_Tool::Tolerance( TopoDS::Edge( nodeShape )); break;
811 case TopAbs_VERTEX: tol = BRep_Tool::Tolerance( TopoDS::Vertex( nodeShape )); break;
814 gp_Pnt nodePnt ( nodes[i]->X(), nodes[i]->Y(), nodes[i]->Z() );
815 BRepAdaptor_Surface surface( geomFace );
816 uvOK = ( nodePnt.Distance( surface.Value( uv.X(), uv.Y() )) < 2 * tol );
818 // normale to geomFace at UV
820 surface.D1( uv.X(), uv.Y(), nodePnt, du, dv );
821 geomNormal = du ^ dv;
822 if ( geomFace.Orientation() == TopAbs_REVERSED )
823 geomNormal.Reverse();
824 geomNormalOK = ( geomNormal.SquareMagnitude() > DBL_MIN * 1e3 );
832 bool isReverse = ( meshNormal * geomNormal ) < 0;
834 return meshDS->ShapeToIndex( solid1 );
836 if ( solids.Extent() == 1 )
837 return HOLE_ID; // we are inside a hole
839 return meshDS->ShapeToIndex( solids(2) );
842 //=======================================================================
843 //function : readResultFile
845 //=======================================================================
847 static bool readResultFile(const int fileOpen,
849 const char* fileName,
852 TopoDS_Shape tabShape[],
855 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
857 int nbEnforcedVertices)
859 MESSAGE("GHS3DPlugin_GHS3D::readResultFile()");
860 Kernel_Utils::Localizer loc;
868 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
871 int nbElems, nbNodes, nbInputNodes;
872 int nodeId/*, triangleId*/;
874 int ID, shapeID, ghs3dShapeID;
877 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
879 int *tab, *tabID, *nodeID, *nodeAssigne;
881 const SMDS_MeshNode **node;
884 //tabID = new int[nbShape];
886 coord = new double[3];
887 node = new const SMDS_MeshNode*[4];
890 SMDS_MeshNode * aNewNode;
891 map <int,const SMDS_MeshNode*>::iterator itOnNode;
892 SMDS_MeshElement* aTet;
897 // Read the file state
898 fileStat = fstat(fileOpen, &status);
899 length = status.st_size;
901 // Mapping the result file into memory
903 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
904 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
905 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
906 0, (DWORD)length, NULL);
907 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
909 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
913 ptr = readMapIntLine(ptr, tab);
918 nbInputNodes = tab[2];
920 nodeAssigne = new int[ nbNodes+1 ];
923 aSolid = tabShape[0];
925 // Reading the nodeId
926 for (int i=0; i < 4*nbElems; i++)
927 nodeId = strtol(ptr, &ptr, 10);
929 MESSAGE("nbInputNodes: "<<nbInputNodes);
930 MESSAGE("nbEnforcedVertices: "<<nbEnforcedVertices);
931 // Reading the nodeCoor and update the nodeMap
932 for (int iNode=1; iNode <= nbNodes; iNode++) {
933 for (int iCoor=0; iCoor < 3; iCoor++)
934 coord[ iCoor ] = strtod(ptr, &ptr);
935 nodeAssigne[ iNode ] = 1;
936 if ( iNode > (nbInputNodes-nbEnforcedVertices) ) {
937 // Creating SMESH nodes
938 // - for enforced vertices
939 // - for vertices of forced edges
941 nodeAssigne[ iNode ] = 0;
942 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
943 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
947 // Reading the number of triangles which corresponds to the number of sub-domains
948 nbTriangle = strtol(ptr, &ptr, 10);
950 tabID = new int[nbTriangle];
951 for (int i=0; i < nbTriangle; i++) {
953 // find the solid corresponding to GHS3D sub-domain following
954 // the technique proposed in GHS3D manual in chapter
955 // "B.4 Subdomain (sub-region) assignment"
956 int nodeId1 = strtol(ptr, &ptr, 10);
957 int nodeId2 = strtol(ptr, &ptr, 10);
958 int nodeId3 = strtol(ptr, &ptr, 10);
959 if ( nbTriangle > 1 ) {
960 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
961 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
962 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
965 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
966 // -- 0020330: Pb with ghs3d as a submesh
967 // check that found shape is to be meshed
968 if ( tabID[i] > 0 ) {
969 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] );
970 bool isToBeMeshed = false;
971 for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS )
972 isToBeMeshed = foundShape.IsSame( tabShape[ iS ]);
976 // END -- 0020330: Pb with ghs3d as a submesh
978 std::cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << std::endl;
981 catch ( Standard_Failure & ex)
984 std::cout << i+1 << " subdomain: Exception caugt: " << ex.GetMessageString() << std::endl;
989 std::cout << i+1 << " subdomain: unknown exception caught " << std::endl;
997 if ( nbTriangle <= nbShape ) // no holes
998 toMeshHoles = true; // not avoid creating tetras in holes
1000 // Associating the tetrahedrons to the shapes
1001 shapeID = compoundID;
1002 for (int iElem = 0; iElem < nbElems; iElem++) {
1003 for (int iNode = 0; iNode < 4; iNode++) {
1004 ID = strtol(tetraPtr, &tetraPtr, 10);
1005 itOnNode = theGhs3dIdToNodeMap.find(ID);
1006 node[ iNode ] = itOnNode->second;
1007 nodeID[ iNode ] = ID;
1009 // We always run GHS3D with "to mesh holes"==TRUE but we must not create
1010 // tetras within holes depending on hypo option,
1011 // so we first check if aTet is inside a hole and then create it
1012 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
1013 if ( nbTriangle > 1 ) {
1014 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
1015 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
1016 if ( tabID[ ghs3dShapeID ] == 0 ) {
1018 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
1019 if ( toMeshHoles || state == TopAbs_IN )
1020 shapeID = theMeshDS->ShapeToIndex( aSolid );
1021 tabID[ ghs3dShapeID ] = shapeID;
1024 shapeID = tabID[ ghs3dShapeID ];
1026 else if ( nbShape > 1 ) {
1027 // Case where nbTriangle == 1 while nbShape == 2 encountered
1028 // with compound of 2 boxes and "To mesh holes"==False,
1029 // so there are no subdomains specified for each tetrahedron.
1030 // Try to guess a solid by a node already bound to shape
1032 for ( int i=0; i<4 && shapeID==0; i++ ) {
1033 if ( nodeAssigne[ nodeID[i] ] == 1 &&
1034 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
1035 node[i]->GetPosition()->GetShapeId() > 1 )
1037 shapeID = node[i]->GetPosition()->GetShapeId();
1041 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
1042 shapeID = theMeshDS->ShapeToIndex( aSolid );
1045 // set new nodes and tetrahedron onto the shape
1046 for ( int i=0; i<4; i++ ) {
1047 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
1048 if ( shapeID != HOLE_ID )
1049 theMeshDS->SetNodeInVolume( node[i], shapeID );
1050 nodeAssigne[ nodeID[i] ] = shapeID;
1053 if ( toMeshHoles || shapeID != HOLE_ID ) {
1054 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
1055 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
1058 shapeIDs.insert( shapeID );
1061 // Remove nodes of tetras inside holes if !toMeshHoles
1062 if ( !toMeshHoles ) {
1063 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
1064 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
1065 ID = itOnNode->first;
1066 if ( nodeAssigne[ ID ] == HOLE_ID )
1067 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
1072 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
1074 UnmapViewOfFile(mapPtr);
1075 CloseHandle(hMapObject);
1078 munmap(mapPtr, length);
1087 delete [] nodeAssigne;
1091 if ( shapeIDs.size() != nbShape ) {
1092 std::cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << std::endl;
1093 for (int i=0; i<nbShape; i++) {
1094 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
1095 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
1096 std::cout << " Solid #" << shapeID << " not found" << std::endl;
1104 //=======================================================================
1105 //function : readResultFile
1107 //=======================================================================
1109 static bool readResultFile(const int fileOpen,
1111 const char* fileName,
1113 SMESH_Mesh& theMesh,
1114 TopoDS_Shape aSolid,
1115 vector <const SMDS_MeshNode*>& theNodeByGhs3dId,
1116 int nbEnforcedVertices)
1118 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
1120 Kernel_Utils::Localizer loc;
1129 int nbElems, nbNodes, nbInputNodes;
1130 int nodeId, triangleId;
1136 const SMDS_MeshNode **node;
1139 coord = new double[3];
1140 node = new const SMDS_MeshNode*[4];
1142 SMDS_MeshNode * aNewNode;
1143 map <int,const SMDS_MeshNode*>::iterator IdNode;
1144 SMDS_MeshElement* aTet;
1146 // Read the file state
1147 fileStat = fstat(fileOpen, &status);
1148 length = status.st_size;
1150 // Mapping the result file into memory
1152 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
1153 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
1154 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
1155 0, (DWORD)length, NULL);
1156 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
1158 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
1162 ptr = readMapIntLine(ptr, tab);
1167 nbInputNodes = tab[2];
1169 theNodeByGhs3dId.resize( nbNodes );
1171 // Reading the nodeId
1172 for (int i=0; i < 4*nbElems; i++)
1173 nodeId = strtol(ptr, &ptr, 10);
1175 // Issue 0020682. Avoid creating nodes and tetras at place where
1176 // volumic elements already exist
1177 SMESH_ElementSearcher* elemSearcher = 0;
1178 vector< const SMDS_MeshElement* > foundVolumes;
1179 if ( theMesh.NbVolumes() > 0 )
1180 elemSearcher = SMESH_MeshEditor( &theMesh ).GetElementSearcher();
1182 // Reading the nodeCoor and update the nodeMap
1183 shapeID = theMeshDS->ShapeToIndex( aSolid );
1184 for (int iNode=0; iNode < nbNodes; iNode++) {
1185 for (int iCoor=0; iCoor < 3; iCoor++)
1186 coord[ iCoor ] = strtod(ptr, &ptr);
1187 if ((iNode+1) > (nbInputNodes-nbEnforcedVertices)) {
1188 // Issue 0020682. Avoid creating nodes and tetras at place where
1189 // volumic elements already exist
1190 if ( elemSearcher &&
1191 elemSearcher->FindElementsByPoint( gp_Pnt(coord[0],coord[1],coord[2]),
1192 SMDSAbs_Volume, foundVolumes ))
1194 theNodeByGhs3dId[ iNode ] = 0;
1198 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
1199 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
1200 theNodeByGhs3dId[ iNode ] = aNewNode;
1205 // Reading the triangles
1206 nbTriangle = strtol(ptr, &ptr, 10);
1208 for (int i=0; i < 3*nbTriangle; i++)
1209 triangleId = strtol(ptr, &ptr, 10);
1213 // Associating the tetrahedrons to the shapes
1214 for (int iElem = 0; iElem < nbElems; iElem++) {
1215 for (int iNode = 0; iNode < 4; iNode++) {
1216 ID = strtol(tetraPtr, &tetraPtr, 10);
1217 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
1221 // Issue 0020682. Avoid creating nodes and tetras at place where
1222 // volumic elements already exist
1223 if ( !node[1] || !node[0] || !node[2] || !node[3] )
1225 if ( elemSearcher->FindElementsByPoint(( SMESH_MeshEditor::TNodeXYZ(node[0]) +
1226 SMESH_MeshEditor::TNodeXYZ(node[1]) +
1227 SMESH_MeshEditor::TNodeXYZ(node[2]) +
1228 SMESH_MeshEditor::TNodeXYZ(node[3]) ) / 4.,
1229 SMDSAbs_Volume, foundVolumes ))
1232 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
1233 shapeID = theMeshDS->ShapeToIndex( aSolid );
1234 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
1237 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
1239 UnmapViewOfFile(mapPtr);
1240 CloseHandle(hMapObject);
1243 munmap(mapPtr, length);
1254 //=============================================================================
1256 *Here we are going to use the GHS3D mesher
1258 //=============================================================================
1260 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1261 const TopoDS_Shape& theShape)
1264 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1266 // we count the number of shapes
1267 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
1269 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1270 for ( ; expBox.More(); expBox.Next() )
1273 // create bounding box for every shape inside the compound
1276 TopoDS_Shape* tabShape;
1278 tabShape = new TopoDS_Shape[_nbShape];
1279 tabBox = new double*[_nbShape];
1280 for (int i=0; i<_nbShape; i++)
1281 tabBox[i] = new double[6];
1282 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1284 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
1285 tabShape[iShape] = expBox.Current();
1286 Bnd_Box BoundingBox;
1287 BRepBndLib::Add(expBox.Current(), BoundingBox);
1288 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1289 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1290 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1291 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1295 // a unique working file name
1296 // to avoid access to the same files by eg different users
1297 TCollection_AsciiString aGenericName
1298 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1300 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1301 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1302 aFacesFileName = aGenericName + ".faces"; // in faces
1303 aPointsFileName = aGenericName + ".points"; // in points
1304 aResultFileName = aGenericName + ".noboite";// out points and volumes
1305 aBadResFileName = aGenericName + ".boite"; // out bad result
1306 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1307 aLogFileName = aGenericName + ".log"; // log
1309 // -----------------
1311 // -----------------
1313 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1314 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1317 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1319 INFOS( "Can't write into " << aFacesFileName);
1320 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1322 map <int,int> aSmdsToGhs3dIdMap;
1323 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1324 map<vector<double>,double> enforcedVertices;
1325 int nbEnforcedVertices = 0;
1327 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1328 nbEnforcedVertices = enforcedVertices.size();
1333 SMESH_MesherHelper helper( theMesh );
1334 helper.SetSubShape( theShape );
1336 // make prisms on quadrangles
1337 vector<StdMeshers_QuadToTriaAdaptor> aQuad2Trias;
1338 if ( theMesh.NbQuadrangles() > 0 )
1340 aQuad2Trias.resize( _nbShape );
1341 for (_iShape = 0, expBox.ReInit(); expBox.More(); expBox.Next())
1342 aQuad2Trias[ _iShape++ ].Compute( theMesh, expBox.Current() );
1345 Ok = (writePoints( aPointsFile, helper, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap, enforcedVertices) &&
1346 writeFaces ( aFacesFile, meshDS, theShape, aQuad2Trias, aSmdsToGhs3dIdMap ));
1348 // Write aSmdsToGhs3dIdMap to temp file
1349 TCollection_AsciiString aSmdsToGhs3dIdMapFileName;
1350 aSmdsToGhs3dIdMapFileName = aGenericName + ".ids"; // ids relation
1351 ofstream aIdsFile ( aSmdsToGhs3dIdMapFileName.ToCString() , ios::out);
1353 aIdsFile.rdbuf()->is_open();
1355 INFOS( "Can't write into " << aSmdsToGhs3dIdMapFileName);
1356 return error(SMESH_Comment("Can't write into ") << aSmdsToGhs3dIdMapFileName);
1358 aIdsFile << "Smds Ghs3d" << std::endl;
1359 map <int,int>::const_iterator myit;
1360 for (myit=aSmdsToGhs3dIdMap.begin() ; myit != aSmdsToGhs3dIdMap.end() ; ++myit) {
1361 aIdsFile << myit->first << " " << myit->second << std::endl;
1365 aPointsFile.close();
1369 if ( !_keepFiles ) {
1370 removeFile( aFacesFileName );
1371 removeFile( aPointsFileName );
1372 removeFile( aSmdsToGhs3dIdMapFileName );
1374 return error(COMPERR_BAD_INPUT_MESH);
1376 removeFile( aResultFileName ); // needed for boundary recovery module usage
1378 // -----------------
1380 // -----------------
1382 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1383 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1384 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1386 std::cout << std::endl;
1387 std::cout << "Ghs3d execution..." << std::endl;
1388 std::cout << cmd << std::endl;
1390 system( cmd.ToCString() ); // run
1392 std::cout << std::endl;
1393 std::cout << "End of Ghs3d execution !" << std::endl;
1399 // Mapping the result file
1402 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1403 if ( fileOpen < 0 ) {
1404 std::cout << std::endl;
1405 std::cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << std::endl;
1406 std::cout << "Log: " << aLogFileName << std::endl;
1411 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1412 Ok = readResultFile( fileOpen,
1414 aResultFileName.ToCString(),
1416 theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1417 toMeshHoles, nbEnforcedVertices );
1420 // ---------------------
1421 // remove working files
1422 // ---------------------
1427 removeFile( aLogFileName );
1429 else if ( OSD_File( aLogFileName ).Size() > 0 )
1431 // get problem description from the log file
1432 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1433 storeErrorDescription( aLogFileName, conv );
1437 // the log file is empty
1438 removeFile( aLogFileName );
1439 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1440 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1443 if ( !_keepFiles ) {
1444 removeFile( aFacesFileName );
1445 removeFile( aPointsFileName );
1446 removeFile( aResultFileName );
1447 removeFile( aBadResFileName );
1448 removeFile( aBbResFileName );
1449 removeFile( aSmdsToGhs3dIdMapFileName );
1451 std::cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1453 std::cout << "not ";
1454 std::cout << "treated !" << std::endl;
1455 std::cout << std::endl;
1457 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1464 //=============================================================================
1466 *Here we are going to use the GHS3D mesher w/o geometry
1468 //=============================================================================
1469 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1470 SMESH_MesherHelper* aHelper)
1472 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1474 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1475 TopoDS_Shape theShape = aHelper->GetSubShape();
1477 // a unique working file name
1478 // to avoid access to the same files by eg different users
1479 TCollection_AsciiString aGenericName
1480 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1482 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1483 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1484 aFacesFileName = aGenericName + ".faces"; // in faces
1485 aPointsFileName = aGenericName + ".points"; // in points
1486 aResultFileName = aGenericName + ".noboite";// out points and volumes
1487 aBadResFileName = aGenericName + ".boite"; // out bad result
1488 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1489 aLogFileName = aGenericName + ".log"; // log
1491 // -----------------
1493 // -----------------
1495 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1496 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1498 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1501 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1503 GHS3DPlugin_Hypothesis::TEnforcedVertexValues enforcedVertices;
1504 int nbEnforcedVertices = 0;
1506 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1507 nbEnforcedVertices = enforcedVertices.size();
1512 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1514 StdMeshers_QuadToTriaAdaptor aQuad2Trias;
1515 if ( theMesh.NbQuadrangles() > 0 )
1516 aQuad2Trias.Compute( theMesh );
1518 Ok = (writeFaces ( aFacesFile, meshDS, aQuad2Trias, aNodeByGhs3dId ) &&
1519 writePoints( aPointsFile, &theMesh, aNodeByGhs3dId,enforcedVertices));
1522 aPointsFile.close();
1525 if ( !_keepFiles ) {
1526 removeFile( aFacesFileName );
1527 removeFile( aPointsFileName );
1529 return error(COMPERR_BAD_INPUT_MESH);
1531 removeFile( aResultFileName ); // needed for boundary recovery module usage
1533 // -----------------
1535 // -----------------
1537 TCollection_AsciiString cmd =
1538 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1539 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1540 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1542 system( cmd.ToCString() ); // run
1548 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1549 if ( fileOpen < 0 ) {
1550 std::cout << std::endl;
1551 std::cout << "Error when opening the " << aResultFileName.ToCString() << " file" << std::endl;
1552 std::cout << "Log: " << aLogFileName << std::endl;
1553 std::cout << std::endl;
1557 Ok = readResultFile( fileOpen,
1559 aResultFileName.ToCString(),
1561 theMesh, theShape ,aNodeByGhs3dId, nbEnforcedVertices );
1564 // ---------------------
1565 // remove working files
1566 // ---------------------
1571 removeFile( aLogFileName );
1573 else if ( OSD_File( aLogFileName ).Size() > 0 )
1575 // get problem description from the log file
1576 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1577 storeErrorDescription( aLogFileName, conv );
1580 // the log file is empty
1581 removeFile( aLogFileName );
1582 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1583 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1588 removeFile( aFacesFileName );
1589 removeFile( aPointsFileName );
1590 removeFile( aResultFileName );
1591 removeFile( aBadResFileName );
1592 removeFile( aBbResFileName );
1598 //================================================================================
1600 * \brief Provide human readable text by error code reported by ghs3d
1602 //================================================================================
1604 static string translateError(const int errNum)
1608 return "The surface mesh includes a face of type other than edge, "
1609 "triangle or quadrilateral. This face type is not supported.";
1611 return "Not enough memory for the face table.";
1613 return "Not enough memory.";
1615 return "Not enough memory.";
1617 return "Face is ignored.";
1619 return "End of file. Some data are missing in the file.";
1621 return "Read error on the file. There are wrong data in the file.";
1623 return "the metric file is inadequate (dimension other than 3).";
1625 return "the metric file is inadequate (values not per vertices).";
1627 return "the metric file contains more than one field.";
1629 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1630 "value of number of mesh vertices in the \".noboite\" file.";
1632 return "Too many sub-domains.";
1634 return "the number of vertices is negative or null.";
1636 return "the number of faces is negative or null.";
1638 return "A face has a null vertex.";
1640 return "incompatible data.";
1642 return "the number of vertices is negative or null.";
1644 return "the number of vertices is negative or null (in the \".mesh\" file).";
1646 return "the number of faces is negative or null.";
1648 return "A face appears more than once in the input surface mesh.";
1650 return "An edge appears more than once in the input surface mesh.";
1652 return "A face has a vertex negative or null.";
1654 return "NOT ENOUGH MEMORY.";
1656 return "Not enough available memory.";
1658 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1659 "in terms of quality or the input list of points is wrong.";
1661 return "Some vertices are too close to one another or coincident.";
1663 return "Some vertices are too close to one another or coincident.";
1665 return "A vertex cannot be inserted.";
1667 return "There are at least two points considered as coincident.";
1669 return "Some vertices are too close to one another or coincident.";
1671 return "The surface mesh regeneration step has failed.";
1673 return "Constrained edge cannot be enforced.";
1675 return "Constrained face cannot be enforced.";
1677 return "Missing faces.";
1679 return "No guess to start the definition of the connected component(s).";
1681 return "The surface mesh includes at least one hole. The domain is not well defined.";
1683 return "Impossible to define a component.";
1685 return "The surface edge intersects another surface edge.";
1687 return "The surface edge intersects the surface face.";
1689 return "One boundary point lies within a surface face.";
1691 return "One surface edge intersects a surface face.";
1693 return "One boundary point lies within a surface edge.";
1695 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1696 "to too many swaps.";
1698 return "Edge is unique (i.e., bounds a hole in the surface).";
1700 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1702 return "Too many components, too many sub-domain.";
1704 return "The surface mesh includes at least one hole. "
1705 "Therefore there is no domain properly defined.";
1707 return "Statistics.";
1709 return "Statistics.";
1711 return "Warning, it is dramatically tedious to enforce the boundary items.";
1713 return "Not enough memory at this time, nevertheless, the program continues. "
1714 "The expected mesh will be correct but not really as large as required.";
1716 return "see above error code, resulting quality may be poor.";
1718 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1720 return "Unknown face type.";
1723 return "End of file. Some data are missing in the file.";
1725 return "A too small volume element is detected.";
1727 return "There exists at least a null or negative volume element.";
1729 return "There exist null or negative volume elements.";
1731 return "A too small volume element is detected. A face is considered being degenerated.";
1733 return "Some element is suspected to be very bad shaped or wrong.";
1735 return "A too bad quality face is detected. This face is considered degenerated.";
1737 return "A too bad quality face is detected. This face is degenerated.";
1739 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1741 return "Abnormal error occured, contact hotline.";
1743 return "Not enough memory for the face table.";
1745 return "The algorithm cannot run further. "
1746 "The surface mesh is probably very bad in terms of quality.";
1748 return "Bad vertex number.";
1753 //================================================================================
1755 * \brief Retrieve from a string given number of integers
1757 //================================================================================
1759 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1762 ids.reserve( nbIds );
1765 while ( !isdigit( *ptr )) ++ptr;
1766 if ( ptr[-1] == '-' ) --ptr;
1767 ids.push_back( strtol( ptr, &ptr, 10 ));
1773 //================================================================================
1775 * \brief Retrieve problem description form a log file
1776 * \retval bool - always false
1778 //================================================================================
1780 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1781 const _Ghs2smdsConvertor & toSmdsConvertor )
1785 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1787 int file = ::open (logFile.ToCString(), O_RDONLY);
1790 return error( SMESH_Comment("See ") << logFile << " for problem description");
1793 // struct stat status;
1794 // fstat(file, &status);
1795 // size_t length = status.st_size;
1796 off_t length = lseek( file, 0, SEEK_END);
1797 lseek( file, 0, SEEK_SET);
1800 vector< char > buf( length );
1801 int nBytesRead = ::read (file, & buf[0], length);
1803 char* ptr = & buf[0];
1804 char* bufEnd = ptr + nBytesRead;
1806 SMESH_Comment errDescription;
1808 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1810 // look for errors "ERR #"
1812 set<string> foundErrorStr; // to avoid reporting same error several times
1813 set<int> elemErrorNums; // not to report different types of errors with bad elements
1814 while ( ++ptr < bufEnd )
1816 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1819 list<const SMDS_MeshElement*> badElems;
1820 vector<int> nodeIds;
1824 int errNum = strtol(ptr, &ptr, 10);
1825 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1827 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1828 ptr = getIds(ptr, NODE, nodeIds);
1829 ptr = getIds(ptr, TRIA, nodeIds);
1830 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1832 case 1000: // ERR 1000 : 1 3 2
1833 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1834 ptr = getIds(ptr, TRIA, nodeIds);
1835 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1838 // Edge (e1, e2) appears more than once in the input surface mesh
1839 ptr = getIds(ptr, EDGE, nodeIds);
1840 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1843 // Face (f 1, f 2, f 3) has a vertex negative or null
1844 ptr = getIds(ptr, TRIA, nodeIds);
1845 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1848 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1849 ptr = getIds(ptr, NODE, nodeIds);
1850 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1851 ptr = getIds(ptr, NODE, nodeIds);
1852 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1855 // Vertex v1 cannot be inserted (warning).
1856 ptr = getIds(ptr, NODE, nodeIds);
1857 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1860 // There are at least two points whose distance is dist, i.e., considered as coincident
1861 case 2103: // ERR 2103 : 16 WITH 3
1862 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1863 ptr = getIds(ptr, NODE, nodeIds);
1864 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1865 ptr = getIds(ptr, NODE, nodeIds);
1866 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1869 // Constrained edge (e1, e2) cannot be enforced (warning).
1870 ptr = getIds(ptr, EDGE, nodeIds);
1871 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1874 // Constrained face (f 1, f 2, f 3) cannot be enforced
1875 ptr = getIds(ptr, TRIA, nodeIds);
1876 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1878 case 3103: // ERR 3103 : 1 2 WITH 7 3
1879 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1880 ptr = getIds(ptr, EDGE, nodeIds);
1881 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1882 ptr = getIds(ptr, EDGE, nodeIds);
1883 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1885 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1886 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1887 ptr = getIds(ptr, EDGE, nodeIds);
1888 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1889 ptr = getIds(ptr, TRIA, nodeIds);
1890 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1892 case 3105: // ERR 3105 : 8 IN 2 3 5
1893 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1894 ptr = getIds(ptr, NODE, nodeIds);
1895 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1896 ptr = getIds(ptr, TRIA, nodeIds);
1897 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1900 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1901 ptr = getIds(ptr, EDGE, nodeIds);
1902 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1903 ptr = getIds(ptr, TRIA, nodeIds);
1904 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1906 case 3107: // ERR 3107 : 2 IN 4 1
1907 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1908 ptr = getIds(ptr, NODE, nodeIds);
1909 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1910 ptr = getIds(ptr, EDGE, nodeIds);
1911 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1913 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1914 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1915 ptr = getIds(ptr, EDGE, nodeIds);
1916 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1918 case 9000: // ERR 9000
1919 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1920 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1921 // A too small volume element is detected. Are reported the index of the element,
1922 // its four vertex indices, its volume and the tolerance threshold value
1923 ptr = getIds(ptr, ID, nodeIds);
1924 ptr = getIds(ptr, VOL, nodeIds);
1925 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1926 // even if all nodes found, volume it most probably invisible,
1927 // add its faces to demenstrate it anyhow
1929 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1930 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1931 faceNodes[2] = nodeIds[3]; // 013
1932 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1933 faceNodes[1] = nodeIds[2]; // 023
1934 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1935 faceNodes[0] = nodeIds[1]; // 123
1936 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1939 case 9001: // ERR 9001
1940 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1941 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1942 // %% NUMBER OF NULL VOLUME TETS : 0
1943 // There exists at least a null or negative volume element
1946 // There exist n null or negative volume elements
1949 // A too small volume element is detected
1952 // A too bad quality face is detected. This face is considered degenerated,
1953 // its index, its three vertex indices together with its quality value are reported
1954 break; // same as next
1955 case 9112: // ERR 9112
1956 // FACE 2 WITH VERTICES : 4 2 5
1957 // SMALL INRADIUS : 0.
1958 // A too bad quality face is detected. This face is degenerated,
1959 // its index, its three vertex indices together with its inradius are reported
1960 ptr = getIds(ptr, ID, nodeIds);
1961 ptr = getIds(ptr, TRIA, nodeIds);
1962 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1963 // add triangle edges as it most probably has zero area and hence invisible
1965 vector<int> edgeNodes(2);
1966 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1967 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1968 edgeNodes[1] = nodeIds[2]; // 0-2
1969 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1970 edgeNodes[0] = nodeIds[1]; // 1-2
1971 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1976 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1978 continue; // not to report same error several times
1980 // const SMDS_MeshElement* nullElem = 0;
1981 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1983 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1984 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1985 // if ( oneMoreErrorType )
1986 // continue; // not to report different types of errors with bad elements
1989 // store bad elements
1990 //if ( allElemsOk ) {
1991 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1992 for ( ; elem != badElems.end(); ++elem )
1993 addBadInputElement( *elem );
1997 string text = translateError( errNum );
1998 if ( errDescription.find( text ) == text.npos ) {
1999 if ( !errDescription.empty() )
2000 errDescription << "\n";
2001 errDescription << text;
2006 if ( errDescription.empty() ) { // no errors found
2007 char msg[] = "connection to server failed";
2008 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
2009 errDescription << "Licence problems.";
2012 if ( errDescription.empty() )
2013 errDescription << "See " << logFile << " for problem description";
2015 errDescription << "\nSee " << logFile << " for more information";
2017 return error( errDescription );
2020 //================================================================================
2022 * \brief Creates _Ghs2smdsConvertor
2024 //================================================================================
2026 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
2027 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
2031 //================================================================================
2033 * \brief Creates _Ghs2smdsConvertor
2035 //================================================================================
2037 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
2038 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
2042 //================================================================================
2044 * \brief Return SMDS element by ids of GHS3D nodes
2046 //================================================================================
2048 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
2050 size_t nbNodes = ghsNodes.size();
2051 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
2052 for ( size_t i = 0; i < nbNodes; ++i ) {
2053 int ghsNode = ghsNodes[ i ];
2054 if ( _ghs2NodeMap ) {
2055 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
2056 if ( in == _ghs2NodeMap->end() )
2058 nodes[ i ] = in->second;
2061 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
2063 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
2069 if ( nbNodes == 2 ) {
2070 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
2072 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
2075 if ( nbNodes == 3 ) {
2076 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
2078 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
2082 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
2088 //=============================================================================
2092 //=============================================================================
2093 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
2094 const TopoDS_Shape& aShape,
2095 MapShapeNbElems& aResMap)
2097 int nbtri = 0, nbqua = 0;
2098 double fullArea = 0.0;
2099 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2100 TopoDS_Face F = TopoDS::Face( exp.Current() );
2101 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2102 MapShapeNbElemsItr anIt = aResMap.find(sm);
2103 if( anIt==aResMap.end() ) {
2104 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2105 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
2106 "Submesh can not be evaluated",this));
2109 std::vector<int> aVec = (*anIt).second;
2110 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
2111 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
2113 BRepGProp::SurfaceProperties(F,G);
2114 double anArea = G.Mass();
2118 // collect info from edges
2119 int nb0d_e = 0, nb1d_e = 0;
2120 bool IsQuadratic = false;
2121 bool IsFirst = true;
2122 TopTools_MapOfShape tmpMap;
2123 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2124 TopoDS_Edge E = TopoDS::Edge(exp.Current());
2125 if( tmpMap.Contains(E) )
2128 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
2129 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
2130 std::vector<int> aVec = (*anIt).second;
2131 nb0d_e += aVec[SMDSEntity_Node];
2132 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
2134 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
2140 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
2143 BRepGProp::VolumeProperties(aShape,G);
2144 double aVolume = G.Mass();
2145 double tetrVol = 0.1179*ELen*ELen*ELen;
2146 double CoeffQuality = 0.9;
2147 int nbVols = int(aVolume/tetrVol/CoeffQuality);
2148 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2149 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2150 std::vector<int> aVec(SMDSEntity_Last);
2151 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2153 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2154 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2155 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2158 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2159 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2160 aVec[SMDSEntity_Pyramid] = nbqua;
2162 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2163 aResMap.insert(std::make_pair(sm,aVec));