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
26 //=============================================================================
28 #include "GHS3DPlugin_GHS3D.hxx"
29 #include "GHS3DPlugin_Hypothesis.hxx"
32 #include <Basics_Utils.hxx>
34 #include "SMESH_Gen.hxx"
35 #include "SMESH_Mesh.hxx"
36 #include "SMESH_Comment.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_MeshEditor.hxx"
40 #include "SMDS_MeshElement.hxx"
41 #include "SMDS_MeshNode.hxx"
42 #include "SMDS_FaceOfNodes.hxx"
43 #include "SMDS_VolumeOfNodes.hxx"
45 #include <BRepAdaptor_Surface.hxx>
46 #include <BRepBndLib.hxx>
47 #include <BRepClass3d_SolidClassifier.hxx>
48 #include <BRepTools.hxx>
49 #include <BRep_Tool.hxx>
50 #include <Bnd_Box.hxx>
51 #include <GeomAPI_ProjectPointOnSurf.hxx>
52 #include <OSD_File.hxx>
53 #include <Precision.hxx>
54 #include <Quantity_Parameter.hxx>
55 #include <Standard_ProgramError.hxx>
56 #include <Standard_ErrorHandler.hxx>
57 #include <Standard_Failure.hxx>
59 #include <TopExp_Explorer.hxx>
60 #include <TopTools_IndexedMapOfShape.hxx>
61 #include <TopTools_ListIteratorOfListOfShape.hxx>
63 //#include <BRepClass_FaceClassifier.hxx>
64 #include <TopTools_MapOfShape.hxx>
65 #include <BRepGProp.hxx>
66 #include <GProp_GProps.hxx>
68 #include "utilities.h"
73 #include <sys/sysinfo.h>
78 //#include <Standard_Stream.hxx>
81 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
102 static void removeFile( const TCollection_AsciiString& fileName )
105 OSD_File( fileName ).Remove();
107 catch ( Standard_ProgramError ) {
108 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
112 //=============================================================================
116 //=============================================================================
118 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
119 : SMESH_3D_Algo(hypId, studyId, gen)
121 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
123 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
124 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
127 _compatibleHypothesis.push_back("GHS3D_Parameters");
128 _requireShape = false; // can work without shape
131 //=============================================================================
135 //=============================================================================
137 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
139 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
142 //=============================================================================
146 //=============================================================================
148 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
149 const TopoDS_Shape& aShape,
150 Hypothesis_Status& aStatus )
152 aStatus = SMESH_Hypothesis::HYP_OK;
154 // there is only one compatible Hypothesis so far
157 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
159 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
161 _keepFiles = _hyp->GetKeepFiles();
166 //=======================================================================
167 //function : findShape
169 //=======================================================================
171 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
173 const TopoDS_Shape shape[],
176 TopAbs_State * state = 0)
179 int j, iShape, nbNode = 4;
181 for ( j=0; j<nbNode; j++ ) {
182 gp_XYZ p ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
183 if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) {
190 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
191 if (state) *state = SC.State();
192 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
193 for (iShape = 0; iShape < nShape; iShape++) {
194 aShape = shape[iShape];
195 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
196 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
197 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
198 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
199 if (state) *state = SC.State();
200 if (SC.State() == TopAbs_IN)
208 //=======================================================================
209 //function : readMapIntLine
211 //=======================================================================
213 static char* readMapIntLine(char* ptr, int tab[]) {
215 std::cout << std::endl;
217 for ( int i=0; i<17; i++ ) {
218 intVal = strtol(ptr, &ptr, 10);
225 //=======================================================================
226 //function : countShape
228 //=======================================================================
230 template < class Mesh, class Shape >
231 static int countShape( Mesh* mesh, Shape shape ) {
232 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
234 for ( ; expShape.More(); expShape.Next() ) {
240 //=======================================================================
241 //function : writeFaces
243 //=======================================================================
245 static bool writeFaces (ofstream & theFile,
246 SMESHDS_Mesh * theMesh,
247 const map <int,int> & theSmdsToGhs3dIdMap)
251 // NB_ELEMS DUMMY_INT
252 // Loop from 1 to NB_ELEMS
253 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
255 int nbShape = countShape( theMesh, TopAbs_FACE );
257 int *tabID; tabID = new int[nbShape];
258 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
260 SMESHDS_SubMesh* theSubMesh;
261 const SMDS_MeshElement* aFace;
262 const char* space = " ";
263 const int dummyint = 0;
264 map<int,int>::const_iterator itOnMap;
265 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
266 int shapeID, nbNodes, aSmdsID;
269 std::cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << std::endl;
270 std::cout << std::endl;
272 theFile << space << theMesh->NbFaces() << space << dummyint << std::endl;
274 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
275 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
277 aShape = expface.Current();
278 shapeID = theMesh->ShapeToIndex( aShape );
280 for ( int j=0; j<=i; j++) {
281 if ( shapeID == tabID[j] ) {
288 tabShape[i] = aShape;
291 for ( int i =0; i < nbShape; i++ ) {
292 if ( tabID[i] != 0 ) {
293 aShape = tabShape[i];
295 theSubMesh = theMesh->MeshElements( aShape );
296 if ( !theSubMesh ) continue;
297 itOnSubMesh = theSubMesh->GetElements();
298 while ( itOnSubMesh->more() ) {
299 aFace = itOnSubMesh->next();
300 nbNodes = aFace->NbNodes();
302 theFile << space << nbNodes;
304 itOnSubFace = aFace->nodesIterator();
305 while ( itOnSubFace->more() ) {
307 aSmdsID = itOnSubFace->next()->GetID();
308 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
309 // if ( itOnMap == theSmdsToGhs3dIdMap.end() ) {
310 // cout << "not found node: " << aSmdsID << endl;
313 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
315 theFile << space << (*itOnMap).second;
318 // (NB_NODES + 1) times: DUMMY_INT
319 for ( int j=0; j<=nbNodes; j++)
320 theFile << space << dummyint;
322 theFile << std::endl;
334 //=======================================================================
335 //function : writeFaces
336 //purpose : Write Faces in case if generate 3D mesh w/o geometry
337 //=======================================================================
339 static bool writeFaces (ofstream & theFile,
340 SMESHDS_Mesh * theMesh,
341 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
345 // NB_ELEMS DUMMY_INT
346 // Loop from 1 to NB_ELEMS
347 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
351 list< const SMDS_MeshElement* > faces;
352 list< const SMDS_MeshElement* >::iterator f;
353 map< const SMDS_MeshNode*,int >::iterator it;
354 SMDS_ElemIteratorPtr nodeIt;
355 const SMDS_MeshElement* elem;
358 const char* space = " ";
359 const int dummyint = 0;
361 //get all faces from mesh
362 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
363 while ( eIt->more() ) {
364 const SMDS_MeshElement* elem = eIt->next();
367 faces.push_back( elem );
374 std::cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
376 // NB_ELEMS DUMMY_INT
377 theFile << space << nbFaces << space << dummyint << std::endl;
379 // Loop from 1 to NB_ELEMS
381 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
383 for ( ; f != faces.end(); ++f )
387 nbNodes = elem->NbNodes();
388 theFile << space << nbNodes;
390 // NODE_NB_1 NODE_NB_2 ...
391 nodeIt = elem->nodesIterator();
392 while ( nodeIt->more() )
395 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
396 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
397 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
398 theFile << space << it->second;
401 // (NB_NODES + 1) times: DUMMY_INT
402 for ( int i=0; i<=nbNodes; i++)
403 theFile << space << dummyint;
405 theFile << std::endl;
408 // put nodes to theNodeByGhs3dId vector
409 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
410 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
411 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
413 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
419 //=======================================================================
420 //function : writePoints
422 //=======================================================================
424 static bool writePoints (ofstream & theFile,
425 SMESH_MesherHelper& theHelper,
426 map <int,int> & theSmdsToGhs3dIdMap,
427 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap,
428 map<vector<double>,double> & theEnforcedVertices)
433 // Loop from 1 to NB_NODES
436 SMESHDS_Mesh * theMesh = theHelper.GetMeshDS();
437 int nbNodes = theMesh->NbNodes();
440 int nbEnforcedVertices = theEnforcedVertices.size();
442 // Issue 020674: EDF 870 SMESH: Mesh generated by Netgen not usable by GHS3D
443 // The problem is in nodes on degenerated edges, we need to skip them
444 if ( theHelper.HasDegeneratedEdges() )
446 // here we decrease total nb of nodes by nb of nodes on degenerated edges
448 for (TopExp_Explorer e(theMesh->ShapeToMesh(), TopAbs_EDGE ); e.More(); e.Next())
450 SMESH_subMesh* sm = theHelper.GetMesh()->GetSubMesh( e.Current() );
451 if ( checkedSM.insert( sm->GetId() ).second && theHelper.IsDegenShape(sm->GetId() )) {
452 if ( sm->GetSubMeshDS() )
453 nbNodes -= sm->GetSubMeshDS()->NbNodes();
457 const char* space = " ";
458 const int dummyint = 0;
461 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
462 const SMDS_MeshNode* node;
465 std::cout << std::endl;
466 std::cout << "The initial 2D mesh contains :" << std::endl;
467 std::cout << " " << nbNodes << " nodes" << std::endl;
468 if (nbEnforcedVertices > 0) {
469 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
471 std::cout << std::endl;
472 std::cout << "Start writing in 'points' file ..." << std::endl;
473 theFile << space << nbNodes << std::endl;
475 // Loop from 1 to NB_NODES
480 if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE &&
481 theHelper.IsDegenShape( node->GetPosition()->GetShapeId() )) // Issue 020674
484 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
485 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
490 << space << node->X()
491 << space << node->Y()
492 << space << node->Z()
493 << space << dummyint;
495 theFile << std::endl;
499 // Iterate over the enforced vertices
500 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
501 const TopoDS_Shape shapeToMesh = theMesh->ShapeToMesh();
502 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
503 double x = vertexIt->first[0];
504 double y = vertexIt->first[1];
505 double z = vertexIt->first[2];
506 // Test if point is inside shape to mesh
507 gp_Pnt myPoint(x,y,z);
508 BRepClass3d_SolidClassifier scl(shapeToMesh);
509 scl.Perform(myPoint, 1e-7);
510 TopAbs_State result = scl.State();
511 if ( result == TopAbs_IN ) {
512 MESSAGE("Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second);
513 // X Y Z PHY_SIZE DUMMY_INT
518 << space << vertexIt->second
519 << space << dummyint;
521 theFile << std::endl;
524 MESSAGE("Enforced vertex (" << x << "," << y <<"," << z << ") is not inside the geometry: it was not added ");
528 std::cout << std::endl;
529 std::cout << "End writing in 'points' file." << std::endl;
534 //=======================================================================
535 //function : writePoints
537 //=======================================================================
539 static bool writePoints (ofstream & theFile,
540 SMESH_Mesh * theMesh,
541 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId,
542 map<vector<double>,double> & theEnforcedVertices)
547 // Loop from 1 to NB_NODES
550 int nbNodes = theNodeByGhs3dId.size();
554 int nbEnforcedVertices = theEnforcedVertices.size();
556 const char* space = " ";
557 const int dummyint = 0;
559 const SMDS_MeshNode* node;
562 std::cout << std::endl;
563 std::cout << "The initial 2D mesh contains :" << std::endl;
564 std::cout << " " << nbNodes << " nodes" << std::endl;
565 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
566 std::cout << std::endl;
567 std::cout << "Start writing in 'points' file ..." << std::endl;
568 theFile << space << nbNodes << std::endl;
570 // Loop from 1 to NB_NODES
572 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
573 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
574 for ( ; nodeIt != after; ++nodeIt )
580 << space << node->X()
581 << space << node->Y()
582 << space << node->Z()
583 << space << dummyint;
585 theFile << std::endl;
589 // Iterate over the enforced vertices
590 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
591 auto_ptr< SMESH_ElementSearcher > pntCls ( SMESH_MeshEditor( theMesh ).GetElementSearcher());
592 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
593 double x = vertexIt->first[0];
594 double y = vertexIt->first[1];
595 double z = vertexIt->first[2];
596 // Test if point is inside shape to mesh
597 gp_Pnt myPoint(x,y,z);
598 TopAbs_State result = pntCls->GetPointState( myPoint );
599 if ( result == TopAbs_IN ) {
600 std::cout << "Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second << std::endl;
602 // X Y Z PHY_SIZE DUMMY_INT
607 << space << vertexIt->second
608 << space << dummyint;
610 theFile << std::endl;
613 std::cout << std::endl;
614 std::cout << "End writing in 'points' file." << std::endl;
619 //=======================================================================
620 //function : findShapeID
621 //purpose : find the solid corresponding to GHS3D sub-domain following
622 // the technique proposed in GHS3D manual in chapter
623 // "B.4 Subdomain (sub-region) assignment"
624 //=======================================================================
626 static int findShapeID(SMESH_Mesh& mesh,
627 const SMDS_MeshNode* node1,
628 const SMDS_MeshNode* node2,
629 const SMDS_MeshNode* node3,
630 const bool toMeshHoles)
632 const int invalidID = 0;
633 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
635 // face the nodes belong to
636 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
640 // geom face the face assigned to
641 SMESH_MeshEditor editor(&mesh);
642 int geomFaceID = editor.FindShape( face );
645 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
646 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
648 TopoDS_Face geomFace = TopoDS::Face( shape );
650 // solids bounded by geom face
651 TopTools_IndexedMapOfShape solids, shells;
652 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
653 for ( ; ansIt.More(); ansIt.Next() ) {
654 switch ( ansIt.Value().ShapeType() ) {
656 solids.Add( ansIt.Value() ); break;
658 shells.Add( ansIt.Value() ); break;
662 // analyse found solids
663 if ( solids.Extent() == 0 || shells.Extent() == 0)
666 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
667 if ( solids.Extent() == 1 )
670 return meshDS->ShapeToIndex( solid1 );
672 // - Are we at a hole boundary face?
673 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
674 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
676 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
677 // check if any edge of shells(1) belongs to another shell
678 for ( ; eExp.More() && !touch; eExp.Next() ) {
679 ansIt = mesh.GetAncestors( eExp.Current() );
680 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
681 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
682 touch = ( !ansIt.Value().IsSame( shells(1) ));
686 return meshDS->ShapeToIndex( solid1 );
689 // find orientation of geom face within the first solid
690 TopExp_Explorer fExp( solid1, TopAbs_FACE );
691 for ( ; fExp.More(); fExp.Next() )
692 if ( geomFace.IsSame( fExp.Current() )) {
693 geomFace = TopoDS::Face( fExp.Current() );
697 return invalidID; // face not found
699 // normale to triangle
700 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
701 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
702 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
703 gp_Vec vec12( node1Pnt, node2Pnt );
704 gp_Vec vec13( node1Pnt, node3Pnt );
705 gp_Vec meshNormal = vec12 ^ vec13;
706 if ( meshNormal.SquareMagnitude() < DBL_MIN )
709 // find UV of node1 on geomFace
710 SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
711 const SMDS_MeshNode* nNotOnSeamEdge = 0;
712 if ( helper.IsSeamShape( node1->GetPosition()->GetShapeId() ))
713 if ( helper.IsSeamShape( node2->GetPosition()->GetShapeId() ))
714 nNotOnSeamEdge = node3;
716 nNotOnSeamEdge = node2;
718 gp_XY uv = helper.GetNodeUV( geomFace, node1, nNotOnSeamEdge, &uvOK );
719 // check that uv is correct
721 TopoDS_Shape nodeShape = helper.GetSubShapeByNode( node1, meshDS );
722 if ( !nodeShape.IsNull() )
723 switch ( nodeShape.ShapeType() )
725 case TopAbs_FACE: tol = BRep_Tool::Tolerance( TopoDS::Face( nodeShape )); break;
726 case TopAbs_EDGE: tol = BRep_Tool::Tolerance( TopoDS::Edge( nodeShape )); break;
727 case TopAbs_VERTEX: tol = BRep_Tool::Tolerance( TopoDS::Vertex( nodeShape )); break;
730 BRepAdaptor_Surface surface( geomFace );
731 if ( !uvOK || node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol )
734 // normale to geomFace at UV
736 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
737 gp_Vec geomNormal = du ^ dv;
738 if ( geomNormal.SquareMagnitude() < DBL_MIN )
739 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
740 if ( geomFace.Orientation() == TopAbs_REVERSED )
741 geomNormal.Reverse();
744 bool isReverse = ( meshNormal * geomNormal ) < 0;
746 return meshDS->ShapeToIndex( solid1 );
748 if ( solids.Extent() == 1 )
749 return HOLE_ID; // we are inside a hole
751 return meshDS->ShapeToIndex( solids(2) );
754 //=======================================================================
755 //function : readResultFile
757 //=======================================================================
759 static bool readResultFile(const int fileOpen,
761 const char* fileName,
764 TopoDS_Shape tabShape[],
767 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
769 int nbEnforcedVertices)
771 MESSAGE("GHS3DPlugin_GHS3D::readResultFile()");
772 Kernel_Utils::Localizer loc;
780 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
783 int nbElems, nbNodes, nbInputNodes;
784 int nodeId/*, triangleId*/;
786 int ID, shapeID, ghs3dShapeID;
789 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
791 int *tab, *tabID, *nodeID, *nodeAssigne;
793 const SMDS_MeshNode **node;
796 //tabID = new int[nbShape];
798 coord = new double[3];
799 node = new const SMDS_MeshNode*[4];
802 SMDS_MeshNode * aNewNode;
803 map <int,const SMDS_MeshNode*>::iterator itOnNode;
804 SMDS_MeshElement* aTet;
809 // Read the file state
810 fileStat = fstat(fileOpen, &status);
811 length = status.st_size;
813 // Mapping the result file into memory
815 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
816 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
817 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
818 0, (DWORD)length, NULL);
819 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
821 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
825 ptr = readMapIntLine(ptr, tab);
830 nbInputNodes = tab[2];
832 nodeAssigne = new int[ nbNodes+1 ];
835 aSolid = tabShape[0];
837 // Reading the nodeId
838 for (int i=0; i < 4*nbElems; i++)
839 nodeId = strtol(ptr, &ptr, 10);
841 MESSAGE("nbInputNodes: "<<nbInputNodes);
842 MESSAGE("nbEnforcedVertices: "<<nbEnforcedVertices);
843 // Reading the nodeCoor and update the nodeMap
844 for (int iNode=1; iNode <= nbNodes; iNode++) {
845 for (int iCoor=0; iCoor < 3; iCoor++)
846 coord[ iCoor ] = strtod(ptr, &ptr);
847 nodeAssigne[ iNode ] = 1;
848 if ( iNode > (nbInputNodes-nbEnforcedVertices) ) {
849 // Creating SMESH nodes
850 // - for enforced vertices
851 // - for vertices of forced edges
853 nodeAssigne[ iNode ] = 0;
854 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
855 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
859 // Reading the number of triangles which corresponds to the number of sub-domains
860 nbTriangle = strtol(ptr, &ptr, 10);
862 tabID = new int[nbTriangle];
863 for (int i=0; i < nbTriangle; i++) {
865 // find the solid corresponding to GHS3D sub-domain following
866 // the technique proposed in GHS3D manual in chapter
867 // "B.4 Subdomain (sub-region) assignment"
868 int nodeId1 = strtol(ptr, &ptr, 10);
869 int nodeId2 = strtol(ptr, &ptr, 10);
870 int nodeId3 = strtol(ptr, &ptr, 10);
871 if ( nbTriangle > 1 ) {
872 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
873 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
874 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
877 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
878 // -- 0020330: Pb with ghs3d as a submesh
879 // check that found shape is to be meshed
880 if ( tabID[i] > 0 ) {
881 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] );
882 bool isToBeMeshed = false;
883 for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS )
884 isToBeMeshed = foundShape.IsSame( tabShape[ iS ]);
888 // END -- 0020330: Pb with ghs3d as a submesh
890 std::cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << std::endl;
892 } catch ( Standard_Failure ) {
899 if ( nbTriangle <= nbShape ) // no holes
900 toMeshHoles = true; // not avoid creating tetras in holes
902 // Associating the tetrahedrons to the shapes
903 shapeID = compoundID;
904 for (int iElem = 0; iElem < nbElems; iElem++) {
905 for (int iNode = 0; iNode < 4; iNode++) {
906 ID = strtol(tetraPtr, &tetraPtr, 10);
907 itOnNode = theGhs3dIdToNodeMap.find(ID);
908 node[ iNode ] = itOnNode->second;
909 nodeID[ iNode ] = ID;
911 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
912 // tetras within holes depending on hypo option,
913 // so we first check if aTet is inside a hole and then create it
914 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
915 if ( nbTriangle > 1 ) {
916 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
917 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
918 if ( tabID[ ghs3dShapeID ] == 0 ) {
920 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
921 if ( toMeshHoles || state == TopAbs_IN )
922 shapeID = theMeshDS->ShapeToIndex( aSolid );
923 tabID[ ghs3dShapeID ] = shapeID;
926 shapeID = tabID[ ghs3dShapeID ];
928 else if ( nbShape > 1 ) {
929 // Case where nbTriangle == 1 while nbShape == 2 encountered
930 // with compound of 2 boxes and "To mesh holes"==False,
931 // so there are no subdomains specified for each tetrahedron.
932 // Try to guess a solid by a node already bound to shape
934 for ( int i=0; i<4 && shapeID==0; i++ ) {
935 if ( nodeAssigne[ nodeID[i] ] == 1 &&
936 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
937 node[i]->GetPosition()->GetShapeId() > 1 )
939 shapeID = node[i]->GetPosition()->GetShapeId();
943 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
944 shapeID = theMeshDS->ShapeToIndex( aSolid );
947 // set new nodes and tetrahedron onto the shape
948 for ( int i=0; i<4; i++ ) {
949 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
950 if ( shapeID != HOLE_ID )
951 theMeshDS->SetNodeInVolume( node[i], shapeID );
952 nodeAssigne[ nodeID[i] ] = shapeID;
955 if ( toMeshHoles || shapeID != HOLE_ID ) {
956 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
957 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
960 shapeIDs.insert( shapeID );
963 // Remove nodes of tetras inside holes if !toMeshHoles
964 if ( !toMeshHoles ) {
965 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
966 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
967 ID = itOnNode->first;
968 if ( nodeAssigne[ ID ] == HOLE_ID )
969 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
974 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
976 UnmapViewOfFile(mapPtr);
977 CloseHandle(hMapObject);
980 munmap(mapPtr, length);
989 delete [] nodeAssigne;
992 if ( shapeIDs.size() != nbShape ) {
993 std::cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << std::endl;
994 for (int i=0; i<nbShape; i++) {
995 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
996 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
997 std::cout << " Solid #" << shapeID << " not found" << std::endl;
1005 //=======================================================================
1006 //function : readResultFile
1008 //=======================================================================
1010 static bool readResultFile(const int fileOpen,
1012 const char* fileName,
1014 SMESH_Mesh& theMesh,
1015 TopoDS_Shape aSolid,
1016 vector <const SMDS_MeshNode*>& theNodeByGhs3dId,
1017 int nbEnforcedVertices)
1019 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
1021 Kernel_Utils::Localizer loc;
1030 int nbElems, nbNodes, nbInputNodes;
1031 int nodeId, triangleId;
1037 const SMDS_MeshNode **node;
1040 coord = new double[3];
1041 node = new const SMDS_MeshNode*[4];
1043 SMDS_MeshNode * aNewNode;
1044 map <int,const SMDS_MeshNode*>::iterator IdNode;
1045 SMDS_MeshElement* aTet;
1047 // Read the file state
1048 fileStat = fstat(fileOpen, &status);
1049 length = status.st_size;
1051 // Mapping the result file into memory
1053 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
1054 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
1055 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
1056 0, (DWORD)length, NULL);
1057 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
1059 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
1063 ptr = readMapIntLine(ptr, tab);
1068 nbInputNodes = tab[2];
1070 theNodeByGhs3dId.resize( nbNodes );
1072 // Reading the nodeId
1073 for (int i=0; i < 4*nbElems; i++)
1074 nodeId = strtol(ptr, &ptr, 10);
1076 // Issue 0020682. Avoid creating nodes and tetras at place where
1077 // volumic elements already exist
1078 SMESH_ElementSearcher* elemSearcher = 0;
1079 vector< const SMDS_MeshElement* > foundVolumes;
1080 if ( theMesh.NbVolumes() > 0 )
1081 elemSearcher = SMESH_MeshEditor( &theMesh ).GetElementSearcher();
1083 // Reading the nodeCoor and update the nodeMap
1084 shapeID = theMeshDS->ShapeToIndex( aSolid );
1085 for (int iNode=0; iNode < nbNodes; iNode++) {
1086 for (int iCoor=0; iCoor < 3; iCoor++)
1087 coord[ iCoor ] = strtod(ptr, &ptr);
1088 if ((iNode+1) > (nbInputNodes-nbEnforcedVertices)) {
1089 // Issue 0020682. Avoid creating nodes and tetras at place where
1090 // volumic elements already exist
1091 if ( elemSearcher &&
1092 elemSearcher->FindElementsByPoint( gp_Pnt(coord[0],coord[1],coord[2]),
1093 SMDSAbs_Volume, foundVolumes ))
1095 theNodeByGhs3dId[ iNode ] = 0;
1099 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
1100 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
1101 theNodeByGhs3dId[ iNode ] = aNewNode;
1106 // Reading the triangles
1107 nbTriangle = strtol(ptr, &ptr, 10);
1109 for (int i=0; i < 3*nbTriangle; i++)
1110 triangleId = strtol(ptr, &ptr, 10);
1114 // Associating the tetrahedrons to the shapes
1115 for (int iElem = 0; iElem < nbElems; iElem++) {
1116 for (int iNode = 0; iNode < 4; iNode++) {
1117 ID = strtol(tetraPtr, &tetraPtr, 10);
1118 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
1122 // Issue 0020682. Avoid creating nodes and tetras at place where
1123 // volumic elements already exist
1124 if ( !node[1] || !node[0] || !node[2] || !node[3] )
1126 if ( elemSearcher->FindElementsByPoint(( SMESH_MeshEditor::TNodeXYZ(node[0]) +
1127 SMESH_MeshEditor::TNodeXYZ(node[1]) +
1128 SMESH_MeshEditor::TNodeXYZ(node[2]) +
1129 SMESH_MeshEditor::TNodeXYZ(node[3]) ) / 4.,
1130 SMDSAbs_Volume, foundVolumes ))
1133 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
1134 shapeID = theMeshDS->ShapeToIndex( aSolid );
1135 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
1138 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
1140 UnmapViewOfFile(mapPtr);
1141 CloseHandle(hMapObject);
1144 munmap(mapPtr, length);
1155 //=============================================================================
1157 *Here we are going to use the GHS3D mesher
1159 //=============================================================================
1161 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1162 const TopoDS_Shape& theShape)
1165 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1167 // we count the number of shapes
1168 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
1170 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1171 for ( ; expBox.More(); expBox.Next() )
1174 // create bounding box for every shape inside the compound
1177 TopoDS_Shape* tabShape;
1179 tabShape = new TopoDS_Shape[_nbShape];
1180 tabBox = new double*[_nbShape];
1181 for (int i=0; i<_nbShape; i++)
1182 tabBox[i] = new double[6];
1183 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1185 // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh
1186 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
1187 tabShape[iShape] = expBox.Current();
1188 Bnd_Box BoundingBox;
1189 BRepBndLib::Add(expBox.Current(), BoundingBox);
1190 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1191 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1192 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1193 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1197 // a unique working file name
1198 // to avoid access to the same files by eg different users
1199 TCollection_AsciiString aGenericName
1200 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1202 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1203 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1204 aFacesFileName = aGenericName + ".faces"; // in faces
1205 aPointsFileName = aGenericName + ".points"; // in points
1206 aResultFileName = aGenericName + ".noboite";// out points and volumes
1207 aBadResFileName = aGenericName + ".boite"; // out bad result
1208 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1209 aLogFileName = aGenericName + ".log"; // log
1211 // -----------------
1213 // -----------------
1215 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1216 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1219 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1221 INFOS( "Can't write into " << aFacesFileName);
1222 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1224 map <int,int> aSmdsToGhs3dIdMap;
1225 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1226 map<vector<double>,double> enforcedVertices;
1227 int nbEnforcedVertices = 0;
1229 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1230 nbEnforcedVertices = enforcedVertices.size();
1235 SMESH_MesherHelper helper( theMesh );
1236 helper.SetSubShape( theShape );
1238 Ok = writePoints( aPointsFile, helper, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap, enforcedVertices) &&
1239 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1241 // Write aSmdsToGhs3dIdMap to temp file
1242 TCollection_AsciiString aSmdsToGhs3dIdMapFileName;
1243 aSmdsToGhs3dIdMapFileName = aGenericName + ".ids"; // ids relation
1244 ofstream aIdsFile ( aSmdsToGhs3dIdMapFileName.ToCString() , ios::out);
1246 aIdsFile.rdbuf()->is_open();
1248 INFOS( "Can't write into " << aSmdsToGhs3dIdMapFileName);
1249 return error(SMESH_Comment("Can't write into ") << aSmdsToGhs3dIdMapFileName);
1251 aIdsFile << "Smds Ghs3d" << std::endl;
1252 map <int,int>::const_iterator myit;
1253 for (myit=aSmdsToGhs3dIdMap.begin() ; myit != aSmdsToGhs3dIdMap.end() ; ++myit) {
1254 aIdsFile << myit->first << " " << myit->second << std::endl;
1258 aPointsFile.close();
1262 if ( !_keepFiles ) {
1263 removeFile( aFacesFileName );
1264 removeFile( aPointsFileName );
1265 removeFile( aSmdsToGhs3dIdMapFileName );
1267 return error(COMPERR_BAD_INPUT_MESH);
1269 removeFile( aResultFileName ); // needed for boundary recovery module usage
1271 // -----------------
1273 // -----------------
1275 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1276 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1277 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1279 std::cout << std::endl;
1280 std::cout << "Ghs3d execution..." << std::endl;
1281 std::cout << cmd << std::endl;
1283 system( cmd.ToCString() ); // run
1285 std::cout << std::endl;
1286 std::cout << "End of Ghs3d execution !" << std::endl;
1292 // Mapping the result file
1295 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1296 if ( fileOpen < 0 ) {
1297 std::cout << std::endl;
1298 std::cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << std::endl;
1299 std::cout << "Log: " << aLogFileName << std::endl;
1304 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1305 Ok = readResultFile( fileOpen,
1307 aResultFileName.ToCString(),
1309 theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1310 toMeshHoles, nbEnforcedVertices );
1313 // ---------------------
1314 // remove working files
1315 // ---------------------
1320 removeFile( aLogFileName );
1322 else if ( OSD_File( aLogFileName ).Size() > 0 )
1324 // get problem description from the log file
1325 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1326 storeErrorDescription( aLogFileName, conv );
1330 // the log file is empty
1331 removeFile( aLogFileName );
1332 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1333 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1336 if ( !_keepFiles ) {
1337 removeFile( aFacesFileName );
1338 removeFile( aPointsFileName );
1339 removeFile( aResultFileName );
1340 removeFile( aBadResFileName );
1341 removeFile( aBbResFileName );
1342 removeFile( aSmdsToGhs3dIdMapFileName );
1344 std::cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1346 std::cout << "not ";
1347 std::cout << "treated !" << std::endl;
1348 std::cout << std::endl;
1350 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1357 //=============================================================================
1359 *Here we are going to use the GHS3D mesher w/o geometry
1361 //=============================================================================
1362 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1363 SMESH_MesherHelper* aHelper)
1365 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1367 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1368 TopoDS_Shape theShape = aHelper->GetSubShape();
1370 // a unique working file name
1371 // to avoid access to the same files by eg different users
1372 TCollection_AsciiString aGenericName
1373 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1375 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1376 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1377 aFacesFileName = aGenericName + ".faces"; // in faces
1378 aPointsFileName = aGenericName + ".points"; // in points
1379 aResultFileName = aGenericName + ".noboite";// out points and volumes
1380 aBadResFileName = aGenericName + ".boite"; // out bad result
1381 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1382 aLogFileName = aGenericName + ".log"; // log
1384 // -----------------
1386 // -----------------
1388 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1389 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1391 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1394 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1396 GHS3DPlugin_Hypothesis::TEnforcedVertexValues enforcedVertices;
1397 int nbEnforcedVertices = 0;
1399 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1400 nbEnforcedVertices = enforcedVertices.size();
1405 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1407 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1408 writePoints( aPointsFile, &theMesh, aNodeByGhs3dId,enforcedVertices));
1411 aPointsFile.close();
1414 if ( !_keepFiles ) {
1415 removeFile( aFacesFileName );
1416 removeFile( aPointsFileName );
1418 return error(COMPERR_BAD_INPUT_MESH);
1420 removeFile( aResultFileName ); // needed for boundary recovery module usage
1422 // -----------------
1424 // -----------------
1426 TCollection_AsciiString cmd =
1427 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1428 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1429 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1431 system( cmd.ToCString() ); // run
1437 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1438 if ( fileOpen < 0 ) {
1439 std::cout << std::endl;
1440 std::cout << "Error when opening the " << aResultFileName.ToCString() << " file" << std::endl;
1441 std::cout << "Log: " << aLogFileName << std::endl;
1442 std::cout << std::endl;
1446 Ok = readResultFile( fileOpen,
1448 aResultFileName.ToCString(),
1450 theMesh, theShape ,aNodeByGhs3dId, nbEnforcedVertices );
1453 // ---------------------
1454 // remove working files
1455 // ---------------------
1460 removeFile( aLogFileName );
1462 else if ( OSD_File( aLogFileName ).Size() > 0 )
1464 // get problem description from the log file
1465 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1466 storeErrorDescription( aLogFileName, conv );
1469 // the log file is empty
1470 removeFile( aLogFileName );
1471 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1472 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1477 removeFile( aFacesFileName );
1478 removeFile( aPointsFileName );
1479 removeFile( aResultFileName );
1480 removeFile( aBadResFileName );
1481 removeFile( aBbResFileName );
1487 //================================================================================
1489 * \brief Provide human readable text by error code reported by ghs3d
1491 //================================================================================
1493 static string translateError(const int errNum)
1497 return "The surface mesh includes a face of type other than edge, "
1498 "triangle or quadrilateral. This face type is not supported.";
1500 return "Not enough memory for the face table.";
1502 return "Not enough memory.";
1504 return "Not enough memory.";
1506 return "Face is ignored.";
1508 return "End of file. Some data are missing in the file.";
1510 return "Read error on the file. There are wrong data in the file.";
1512 return "the metric file is inadequate (dimension other than 3).";
1514 return "the metric file is inadequate (values not per vertices).";
1516 return "the metric file contains more than one field.";
1518 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1519 "value of number of mesh vertices in the \".noboite\" file.";
1521 return "Too many sub-domains.";
1523 return "the number of vertices is negative or null.";
1525 return "the number of faces is negative or null.";
1527 return "A face has a null vertex.";
1529 return "incompatible data.";
1531 return "the number of vertices is negative or null.";
1533 return "the number of vertices is negative or null (in the \".mesh\" file).";
1535 return "the number of faces is negative or null.";
1537 return "A face appears more than once in the input surface mesh.";
1539 return "An edge appears more than once in the input surface mesh.";
1541 return "A face has a vertex negative or null.";
1543 return "NOT ENOUGH MEMORY.";
1545 return "Not enough available memory.";
1547 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1548 "in terms of quality or the input list of points is wrong.";
1550 return "Some vertices are too close to one another or coincident.";
1552 return "Some vertices are too close to one another or coincident.";
1554 return "A vertex cannot be inserted.";
1556 return "There are at least two points considered as coincident.";
1558 return "Some vertices are too close to one another or coincident.";
1560 return "The surface mesh regeneration step has failed.";
1562 return "Constrained edge cannot be enforced.";
1564 return "Constrained face cannot be enforced.";
1566 return "Missing faces.";
1568 return "No guess to start the definition of the connected component(s).";
1570 return "The surface mesh includes at least one hole. The domain is not well defined.";
1572 return "Impossible to define a component.";
1574 return "The surface edge intersects another surface edge.";
1576 return "The surface edge intersects the surface face.";
1578 return "One boundary point lies within a surface face.";
1580 return "One surface edge intersects a surface face.";
1582 return "One boundary point lies within a surface edge.";
1584 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1585 "to too many swaps.";
1587 return "Edge is unique (i.e., bounds a hole in the surface).";
1589 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1591 return "Too many components, too many sub-domain.";
1593 return "The surface mesh includes at least one hole. "
1594 "Therefore there is no domain properly defined.";
1596 return "Statistics.";
1598 return "Statistics.";
1600 return "Warning, it is dramatically tedious to enforce the boundary items.";
1602 return "Not enough memory at this time, nevertheless, the program continues. "
1603 "The expected mesh will be correct but not really as large as required.";
1605 return "see above error code, resulting quality may be poor.";
1607 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1609 return "Unknown face type.";
1612 return "End of file. Some data are missing in the file.";
1614 return "A too small volume element is detected.";
1616 return "There exists at least a null or negative volume element.";
1618 return "There exist null or negative volume elements.";
1620 return "A too small volume element is detected. A face is considered being degenerated.";
1622 return "Some element is suspected to be very bad shaped or wrong.";
1624 return "A too bad quality face is detected. This face is considered degenerated.";
1626 return "A too bad quality face is detected. This face is degenerated.";
1628 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1630 return "Abnormal error occured, contact hotline.";
1632 return "Not enough memory for the face table.";
1634 return "The algorithm cannot run further. "
1635 "The surface mesh is probably very bad in terms of quality.";
1637 return "Bad vertex number.";
1642 //================================================================================
1644 * \brief Retrieve from a string given number of integers
1646 //================================================================================
1648 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1651 ids.reserve( nbIds );
1654 while ( !isdigit( *ptr )) ++ptr;
1655 if ( ptr[-1] == '-' ) --ptr;
1656 ids.push_back( strtol( ptr, &ptr, 10 ));
1662 //================================================================================
1664 * \brief Retrieve problem description form a log file
1665 * \retval bool - always false
1667 //================================================================================
1669 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1670 const _Ghs2smdsConvertor & toSmdsConvertor )
1674 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1676 int file = ::open (logFile.ToCString(), O_RDONLY);
1679 return error( SMESH_Comment("See ") << logFile << " for problem description");
1682 // struct stat status;
1683 // fstat(file, &status);
1684 // size_t length = status.st_size;
1685 off_t length = lseek( file, 0, SEEK_END);
1686 lseek( file, 0, SEEK_SET);
1689 vector< char > buf( length );
1690 int nBytesRead = ::read (file, & buf[0], length);
1692 char* ptr = & buf[0];
1693 char* bufEnd = ptr + nBytesRead;
1695 SMESH_Comment errDescription;
1697 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1699 // look for errors "ERR #"
1701 set<string> foundErrorStr; // to avoid reporting same error several times
1702 set<int> elemErrorNums; // not to report different types of errors with bad elements
1703 while ( ++ptr < bufEnd )
1705 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1708 list<const SMDS_MeshElement*> badElems;
1709 vector<int> nodeIds;
1713 int errNum = strtol(ptr, &ptr, 10);
1714 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1716 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1717 ptr = getIds(ptr, NODE, nodeIds);
1718 ptr = getIds(ptr, TRIA, nodeIds);
1719 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1721 case 1000: // ERR 1000 : 1 3 2
1722 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1723 ptr = getIds(ptr, TRIA, nodeIds);
1724 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1727 // Edge (e1, e2) appears more than once in the input surface mesh
1728 ptr = getIds(ptr, EDGE, nodeIds);
1729 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1732 // Face (f 1, f 2, f 3) has a vertex negative or null
1733 ptr = getIds(ptr, TRIA, nodeIds);
1734 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1737 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1738 ptr = getIds(ptr, NODE, nodeIds);
1739 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1740 ptr = getIds(ptr, NODE, nodeIds);
1741 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1744 // Vertex v1 cannot be inserted (warning).
1745 ptr = getIds(ptr, NODE, nodeIds);
1746 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1749 // There are at least two points whose distance is dist, i.e., considered as coincident
1750 case 2103: // ERR 2103 : 16 WITH 3
1751 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1752 ptr = getIds(ptr, NODE, nodeIds);
1753 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1754 ptr = getIds(ptr, NODE, nodeIds);
1755 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1758 // Constrained edge (e1, e2) cannot be enforced (warning).
1759 ptr = getIds(ptr, EDGE, nodeIds);
1760 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1763 // Constrained face (f 1, f 2, f 3) cannot be enforced
1764 ptr = getIds(ptr, TRIA, nodeIds);
1765 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1767 case 3103: // ERR 3103 : 1 2 WITH 7 3
1768 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1769 ptr = getIds(ptr, EDGE, nodeIds);
1770 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1771 ptr = getIds(ptr, EDGE, nodeIds);
1772 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1774 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1775 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1776 ptr = getIds(ptr, EDGE, nodeIds);
1777 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1778 ptr = getIds(ptr, TRIA, nodeIds);
1779 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1781 case 3105: // ERR 3105 : 8 IN 2 3 5
1782 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1783 ptr = getIds(ptr, NODE, nodeIds);
1784 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1785 ptr = getIds(ptr, TRIA, nodeIds);
1786 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1789 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1790 ptr = getIds(ptr, EDGE, nodeIds);
1791 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1792 ptr = getIds(ptr, TRIA, nodeIds);
1793 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1795 case 3107: // ERR 3107 : 2 IN 4 1
1796 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1797 ptr = getIds(ptr, NODE, nodeIds);
1798 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1799 ptr = getIds(ptr, EDGE, nodeIds);
1800 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1802 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1803 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1804 ptr = getIds(ptr, EDGE, nodeIds);
1805 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1807 case 9000: // ERR 9000
1808 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1809 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1810 // A too small volume element is detected. Are reported the index of the element,
1811 // its four vertex indices, its volume and the tolerance threshold value
1812 ptr = getIds(ptr, ID, nodeIds);
1813 ptr = getIds(ptr, VOL, nodeIds);
1814 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1815 // even if all nodes found, volume it most probably invisible,
1816 // add its faces to demenstrate it anyhow
1818 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1819 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1820 faceNodes[2] = nodeIds[3]; // 013
1821 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1822 faceNodes[1] = nodeIds[2]; // 023
1823 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1824 faceNodes[0] = nodeIds[1]; // 123
1825 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1828 case 9001: // ERR 9001
1829 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1830 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1831 // %% NUMBER OF NULL VOLUME TETS : 0
1832 // There exists at least a null or negative volume element
1835 // There exist n null or negative volume elements
1838 // A too small volume element is detected
1841 // A too bad quality face is detected. This face is considered degenerated,
1842 // its index, its three vertex indices together with its quality value are reported
1843 break; // same as next
1844 case 9112: // ERR 9112
1845 // FACE 2 WITH VERTICES : 4 2 5
1846 // SMALL INRADIUS : 0.
1847 // A too bad quality face is detected. This face is degenerated,
1848 // its index, its three vertex indices together with its inradius are reported
1849 ptr = getIds(ptr, ID, nodeIds);
1850 ptr = getIds(ptr, TRIA, nodeIds);
1851 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1852 // add triangle edges as it most probably has zero area and hence invisible
1854 vector<int> edgeNodes(2);
1855 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1856 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1857 edgeNodes[1] = nodeIds[2]; // 0-2
1858 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1859 edgeNodes[0] = nodeIds[1]; // 1-2
1860 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1865 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1867 continue; // not to report same error several times
1869 // const SMDS_MeshElement* nullElem = 0;
1870 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1872 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1873 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1874 // if ( oneMoreErrorType )
1875 // continue; // not to report different types of errors with bad elements
1878 // store bad elements
1879 //if ( allElemsOk ) {
1880 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1881 for ( ; elem != badElems.end(); ++elem )
1882 addBadInputElement( *elem );
1886 string text = translateError( errNum );
1887 if ( errDescription.find( text ) == text.npos ) {
1888 if ( !errDescription.empty() )
1889 errDescription << "\n";
1890 errDescription << text;
1895 if ( errDescription.empty() ) { // no errors found
1896 char msg[] = "connection to server failed";
1897 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1898 errDescription << "Licence problems.";
1901 if ( errDescription.empty() )
1902 errDescription << "See " << logFile << " for problem description";
1904 errDescription << "\nSee " << logFile << " for more information";
1906 return error( errDescription );
1909 //================================================================================
1911 * \brief Creates _Ghs2smdsConvertor
1913 //================================================================================
1915 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1916 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1920 //================================================================================
1922 * \brief Creates _Ghs2smdsConvertor
1924 //================================================================================
1926 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1927 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1931 //================================================================================
1933 * \brief Return SMDS element by ids of GHS3D nodes
1935 //================================================================================
1937 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1939 size_t nbNodes = ghsNodes.size();
1940 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1941 for ( size_t i = 0; i < nbNodes; ++i ) {
1942 int ghsNode = ghsNodes[ i ];
1943 if ( _ghs2NodeMap ) {
1944 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1945 if ( in == _ghs2NodeMap->end() )
1947 nodes[ i ] = in->second;
1950 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1952 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1958 if ( nbNodes == 2 ) {
1959 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1961 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1964 if ( nbNodes == 3 ) {
1965 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1967 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1971 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1977 //=============================================================================
1981 //=============================================================================
1982 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
1983 const TopoDS_Shape& aShape,
1984 MapShapeNbElems& aResMap)
1986 int nbtri = 0, nbqua = 0;
1987 double fullArea = 0.0;
1988 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1989 TopoDS_Face F = TopoDS::Face( exp.Current() );
1990 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1991 MapShapeNbElemsItr anIt = aResMap.find(sm);
1992 if( anIt==aResMap.end() ) {
1993 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1994 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1995 "Submesh can not be evaluated",this));
1998 std::vector<int> aVec = (*anIt).second;
1999 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
2000 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
2002 BRepGProp::SurfaceProperties(F,G);
2003 double anArea = G.Mass();
2007 // collect info from edges
2008 int nb0d_e = 0, nb1d_e = 0;
2009 bool IsQuadratic = false;
2010 bool IsFirst = true;
2011 TopTools_MapOfShape tmpMap;
2012 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2013 TopoDS_Edge E = TopoDS::Edge(exp.Current());
2014 if( tmpMap.Contains(E) )
2017 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
2018 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
2019 std::vector<int> aVec = (*anIt).second;
2020 nb0d_e += aVec[SMDSEntity_Node];
2021 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
2023 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
2029 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
2032 BRepGProp::VolumeProperties(aShape,G);
2033 double aVolume = G.Mass();
2034 double tetrVol = 0.1179*ELen*ELen*ELen;
2035 double CoeffQuality = 0.9;
2036 int nbVols = int(aVolume/tetrVol/CoeffQuality);
2037 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2038 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2039 std::vector<int> aVec(SMDSEntity_Last);
2040 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2042 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2043 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2044 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2047 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2048 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2049 aVec[SMDSEntity_Pyramid] = nbqua;
2051 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2052 aResMap.insert(std::make_pair(sm,aVec));