1 // Copyright (C) 2004-2008 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
19 //=============================================================================
20 // File : GHS3DPlugin_GHS3D.cxx
22 // 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 <BRepAdaptor_Surface.hxx>
45 #include <BRepBndLib.hxx>
46 #include <BRepClass3d_SolidClassifier.hxx>
47 #include <BRepTools.hxx>
48 #include <BRep_Tool.hxx>
49 #include <Bnd_Box.hxx>
50 #include <GeomAPI_ProjectPointOnSurf.hxx>
51 #include <OSD_File.hxx>
52 #include <Precision.hxx>
53 #include <Quantity_Parameter.hxx>
54 #include <Standard_ProgramError.hxx>
55 #include <Standard_ErrorHandler.hxx>
56 #include <Standard_Failure.hxx>
58 #include <TopExp_Explorer.hxx>
59 #include <TopTools_IndexedMapOfShape.hxx>
60 #include <TopTools_ListIteratorOfListOfShape.hxx>
62 //#include <BRepClass_FaceClassifier.hxx>
63 #include <TopTools_MapOfShape.hxx>
64 #include <BRepGProp.hxx>
65 #include <GProp_GProps.hxx>
67 #include "utilities.h"
72 #include <sys/sysinfo.h>
77 //#include <Standard_Stream.hxx>
80 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
101 static void removeFile( const TCollection_AsciiString& fileName )
104 OSD_File( fileName ).Remove();
106 catch ( Standard_ProgramError ) {
107 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
111 //=============================================================================
115 //=============================================================================
117 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
118 : SMESH_3D_Algo(hypId, studyId, gen)
120 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
122 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
123 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
126 _compatibleHypothesis.push_back("GHS3D_Parameters");
127 _requireShape = false; // can work without shape
130 //=============================================================================
134 //=============================================================================
136 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
138 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
141 //=============================================================================
145 //=============================================================================
147 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
148 const TopoDS_Shape& aShape,
149 Hypothesis_Status& aStatus )
151 aStatus = SMESH_Hypothesis::HYP_OK;
153 // there is only one compatible Hypothesis so far
156 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
158 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
160 _keepFiles = _hyp->GetKeepFiles();
165 //=======================================================================
166 //function : findShape
168 //=======================================================================
170 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
172 const TopoDS_Shape shape[],
175 TopAbs_State * state = 0)
178 int j, iShape, nbNode = 4;
180 for ( j=0; j<nbNode; j++ ) {
181 gp_XYZ p ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
182 if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) {
189 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
190 if (state) *state = SC.State();
191 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
192 for (iShape = 0; iShape < nShape; iShape++) {
193 aShape = shape[iShape];
194 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
195 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
196 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
197 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
198 if (state) *state = SC.State();
199 if (SC.State() == TopAbs_IN)
207 //=======================================================================
208 //function : readMapIntLine
210 //=======================================================================
212 static char* readMapIntLine(char* ptr, int tab[]) {
214 std::cout << std::endl;
216 for ( int i=0; i<17; i++ ) {
217 intVal = strtol(ptr, &ptr, 10);
224 //=======================================================================
225 //function : countShape
227 //=======================================================================
229 template < class Mesh, class Shape >
230 static int countShape( Mesh* mesh, Shape shape ) {
231 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
233 for ( ; expShape.More(); expShape.Next() ) {
239 //=======================================================================
240 //function : writeFaces
242 //=======================================================================
244 static bool writeFaces (ofstream & theFile,
245 SMESHDS_Mesh * theMesh,
246 const map <int,int> & theSmdsToGhs3dIdMap)
250 // NB_ELEMS DUMMY_INT
251 // Loop from 1 to NB_ELEMS
252 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
254 int nbShape = countShape( theMesh, TopAbs_FACE );
256 int *tabID; tabID = new int[nbShape];
257 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
259 SMESHDS_SubMesh* theSubMesh;
260 const SMDS_MeshElement* aFace;
261 const char* space = " ";
262 const int dummyint = 0;
263 map<int,int>::const_iterator itOnMap;
264 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
265 int shapeID, nbNodes, aSmdsID;
268 std::cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << std::endl;
269 std::cout << std::endl;
271 theFile << space << theMesh->NbFaces() << space << dummyint << std::endl;
273 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
274 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
276 aShape = expface.Current();
277 shapeID = theMesh->ShapeToIndex( aShape );
279 for ( int j=0; j<=i; j++) {
280 if ( shapeID == tabID[j] ) {
287 tabShape[i] = aShape;
290 for ( int i =0; i < nbShape; i++ ) {
291 if ( tabID[i] != 0 ) {
292 aShape = tabShape[i];
294 theSubMesh = theMesh->MeshElements( aShape );
295 if ( !theSubMesh ) continue;
296 itOnSubMesh = theSubMesh->GetElements();
297 while ( itOnSubMesh->more() ) {
298 aFace = itOnSubMesh->next();
299 nbNodes = aFace->NbNodes();
301 theFile << space << nbNodes;
303 itOnSubFace = aFace->nodesIterator();
304 while ( itOnSubFace->more() ) {
306 aSmdsID = itOnSubFace->next()->GetID();
307 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
308 // if ( itOnMap == theSmdsToGhs3dIdMap.end() ) {
309 // cout << "not found node: " << aSmdsID << endl;
312 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
314 theFile << space << (*itOnMap).second;
317 // (NB_NODES + 1) times: DUMMY_INT
318 for ( int j=0; j<=nbNodes; j++)
319 theFile << space << dummyint;
321 theFile << std::endl;
333 //=======================================================================
334 //function : writeFaces
335 //purpose : Write Faces in case if generate 3D mesh w/o geometry
336 //=======================================================================
338 static bool writeFaces (ofstream & theFile,
339 SMESHDS_Mesh * theMesh,
340 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
344 // NB_ELEMS DUMMY_INT
345 // Loop from 1 to NB_ELEMS
346 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
350 list< const SMDS_MeshElement* > faces;
351 list< const SMDS_MeshElement* >::iterator f;
352 map< const SMDS_MeshNode*,int >::iterator it;
353 SMDS_ElemIteratorPtr nodeIt;
354 const SMDS_MeshElement* elem;
357 const char* space = " ";
358 const int dummyint = 0;
360 //get all faces from mesh
361 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
362 while ( eIt->more() ) {
363 const SMDS_MeshElement* elem = eIt->next();
366 faces.push_back( elem );
373 std::cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
375 // NB_ELEMS DUMMY_INT
376 theFile << space << nbFaces << space << dummyint << std::endl;
378 // Loop from 1 to NB_ELEMS
380 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
382 for ( ; f != faces.end(); ++f )
386 nbNodes = elem->NbNodes();
387 theFile << space << nbNodes;
389 // NODE_NB_1 NODE_NB_2 ...
390 nodeIt = elem->nodesIterator();
391 while ( nodeIt->more() )
394 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
395 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
396 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
397 theFile << space << it->second;
400 // (NB_NODES + 1) times: DUMMY_INT
401 for ( int i=0; i<=nbNodes; i++)
402 theFile << space << dummyint;
404 theFile << std::endl;
407 // put nodes to theNodeByGhs3dId vector
408 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
409 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
410 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
412 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
418 //=======================================================================
419 //function : writePoints
421 //=======================================================================
423 static bool writePoints (ofstream & theFile,
424 SMESH_MesherHelper& theHelper,
425 map <int,int> & theSmdsToGhs3dIdMap,
426 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap,
427 map<vector<double>,double> & theEnforcedVertices)
432 // Loop from 1 to NB_NODES
435 SMESHDS_Mesh * theMesh = theHelper.GetMeshDS();
436 int nbNodes = theMesh->NbNodes();
439 int nbEnforcedVertices = theEnforcedVertices.size();
441 // Issue 020674: EDF 870 SMESH: Mesh generated by Netgen not usable by GHS3D
442 // The problem is in nodes on degenerated edges, we need to skip them
443 if ( theHelper.HasDegeneratedEdges() )
445 // here we decrease total nb of nodes by nb of nodes on degenerated edges
447 for (TopExp_Explorer e(theMesh->ShapeToMesh(), TopAbs_EDGE ); e.More(); e.Next())
449 SMESH_subMesh* sm = theHelper.GetMesh()->GetSubMesh( e.Current() );
450 if ( checkedSM.insert( sm->GetId() ).second && theHelper.IsDegenShape(sm->GetId() ))
451 nbNodes -= sm->GetSubMeshDS()->NbNodes();
454 const char* space = " ";
455 const int dummyint = 0;
458 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
459 const SMDS_MeshNode* node;
462 std::cout << std::endl;
463 std::cout << "The initial 2D mesh contains :" << std::endl;
464 std::cout << " " << nbNodes << " nodes" << std::endl;
465 if (nbEnforcedVertices > 0) {
466 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
468 std::cout << std::endl;
469 std::cout << "Start writing in 'points' file ..." << std::endl;
470 theFile << space << nbNodes << std::endl;
472 // Loop from 1 to NB_NODES
477 if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE &&
478 theHelper.IsDegenShape( node->GetPosition()->GetShapeId() )) // Issue 020674
481 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
482 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
487 << space << node->X()
488 << space << node->Y()
489 << space << node->Z()
490 << space << dummyint;
492 theFile << std::endl;
496 // Iterate over the enforced vertices
497 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
498 const TopoDS_Shape shapeToMesh = theMesh->ShapeToMesh();
499 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
500 double x = vertexIt->first[0];
501 double y = vertexIt->first[1];
502 double z = vertexIt->first[2];
503 // Test if point is inside shape to mesh
504 gp_Pnt myPoint(x,y,z);
505 BRepClass3d_SolidClassifier scl(shapeToMesh);
506 scl.Perform(myPoint, 1e-7);
507 TopAbs_State result = scl.State();
508 if ( result == TopAbs_IN ) {
509 MESSAGE("Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second);
510 // X Y Z PHY_SIZE DUMMY_INT
515 << space << vertexIt->second
516 << space << dummyint;
518 theFile << std::endl;
521 MESSAGE("Enforced vertex (" << x << "," << y <<"," << z << ") is not inside the geometry: it was not added ");
525 std::cout << std::endl;
526 std::cout << "End writing in 'points' file." << std::endl;
531 //=======================================================================
532 //function : writePoints
534 //=======================================================================
536 static bool writePoints (ofstream & theFile,
537 SMESHDS_Mesh * theMesh,
538 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId,
539 map<vector<double>,double> & theEnforcedVertices)
544 // Loop from 1 to NB_NODES
547 //int nbNodes = theMesh->NbNodes();
548 int nbNodes = theNodeByGhs3dId.size();
552 int nbEnforcedVertices = theEnforcedVertices.size();
554 const char* space = " ";
555 const int dummyint = 0;
557 const SMDS_MeshNode* node;
560 std::cout << std::endl;
561 std::cout << "The initial 2D mesh contains :" << std::endl;
562 std::cout << " " << nbNodes << " nodes" << std::endl;
563 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
564 std::cout << std::endl;
565 std::cout << "Start writing in 'points' file ..." << std::endl;
566 theFile << space << nbNodes << std::endl;
568 // Loop from 1 to NB_NODES
570 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
571 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
572 for ( ; nodeIt != after; ++nodeIt )
578 << space << node->X()
579 << space << node->Y()
580 << space << node->Z()
581 << space << dummyint;
583 theFile << std::endl;
587 // Iterate over the enforced vertices
588 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
589 const TopoDS_Shape shapeToMesh = theMesh->ShapeToMesh();
590 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
591 double x = vertexIt->first[0];
592 double y = vertexIt->first[1];
593 double z = vertexIt->first[2];
594 // Test if point is inside shape to mesh
595 gp_Pnt myPoint(x,y,z);
596 BRepClass3d_SolidClassifier scl(shapeToMesh);
597 scl.Perform(myPoint, 1e-7);
598 TopAbs_State result = scl.State();
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 SMESHDS_Mesh* theMeshDS,
1015 TopoDS_Shape aSolid,
1016 vector <const SMDS_MeshNode*>& theNodeByGhs3dId,
1017 int nbEnforcedVertices) {
1019 Kernel_Utils::Localizer loc;
1028 int nbElems, nbNodes, nbInputNodes;
1029 int nodeId, triangleId;
1035 const SMDS_MeshNode **node;
1038 coord = new double[3];
1039 node = new const SMDS_MeshNode*[4];
1041 SMDS_MeshNode * aNewNode;
1042 map <int,const SMDS_MeshNode*>::iterator IdNode;
1043 SMDS_MeshElement* aTet;
1045 // Read the file state
1046 fileStat = fstat(fileOpen, &status);
1047 length = status.st_size;
1049 // Mapping the result file into memory
1051 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
1052 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
1053 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
1054 0, (DWORD)length, NULL);
1055 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
1057 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
1061 ptr = readMapIntLine(ptr, tab);
1066 nbInputNodes = tab[2];
1068 theNodeByGhs3dId.resize( nbNodes );
1070 // Reading the nodeId
1071 for (int i=0; i < 4*nbElems; i++)
1072 nodeId = strtol(ptr, &ptr, 10);
1074 // Reading the nodeCoor and update the nodeMap
1075 shapeID = theMeshDS->ShapeToIndex( aSolid );
1076 for (int iNode=0; iNode < nbNodes; iNode++) {
1077 for (int iCoor=0; iCoor < 3; iCoor++)
1078 coord[ iCoor ] = strtod(ptr, &ptr);
1079 if ((iNode+1) > (nbInputNodes-nbEnforcedVertices)) {
1080 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
1081 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
1082 theNodeByGhs3dId[ iNode ] = aNewNode;
1086 // Reading the triangles
1087 nbTriangle = strtol(ptr, &ptr, 10);
1089 for (int i=0; i < 3*nbTriangle; i++)
1090 triangleId = strtol(ptr, &ptr, 10);
1094 // Associating the tetrahedrons to the shapes
1095 for (int iElem = 0; iElem < nbElems; iElem++) {
1096 for (int iNode = 0; iNode < 4; iNode++) {
1097 ID = strtol(tetraPtr, &tetraPtr, 10);
1098 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
1100 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
1101 shapeID = theMeshDS->ShapeToIndex( aSolid );
1102 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
1105 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
1107 UnmapViewOfFile(mapPtr);
1108 CloseHandle(hMapObject);
1111 munmap(mapPtr, length);
1122 //=============================================================================
1124 *Here we are going to use the GHS3D mesher
1126 //=============================================================================
1128 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1129 const TopoDS_Shape& theShape)
1132 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1134 // we count the number of shapes
1135 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
1137 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1138 for ( ; expBox.More(); expBox.Next() )
1141 // create bounding box for every shape inside the compound
1144 TopoDS_Shape* tabShape;
1146 tabShape = new TopoDS_Shape[_nbShape];
1147 tabBox = new double*[_nbShape];
1148 for (int i=0; i<_nbShape; i++)
1149 tabBox[i] = new double[6];
1150 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1152 // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh
1153 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
1154 tabShape[iShape] = expBox.Current();
1155 Bnd_Box BoundingBox;
1156 BRepBndLib::Add(expBox.Current(), BoundingBox);
1157 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1158 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1159 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1160 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1164 // a unique working file name
1165 // to avoid access to the same files by eg different users
1166 TCollection_AsciiString aGenericName
1167 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1169 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1170 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1171 aFacesFileName = aGenericName + ".faces"; // in faces
1172 aPointsFileName = aGenericName + ".points"; // in points
1173 aResultFileName = aGenericName + ".noboite";// out points and volumes
1174 aBadResFileName = aGenericName + ".boite"; // out bad result
1175 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1176 aLogFileName = aGenericName + ".log"; // log
1178 // -----------------
1180 // -----------------
1182 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1183 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1186 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1188 INFOS( "Can't write into " << aFacesFileName);
1189 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1191 map <int,int> aSmdsToGhs3dIdMap;
1192 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1193 map<vector<double>,double> enforcedVertices;
1194 int nbEnforcedVertices = 0;
1196 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1197 nbEnforcedVertices = enforcedVertices.size();
1202 SMESH_MesherHelper helper( theMesh );
1203 helper.SetSubShape( theShape );
1205 Ok = writePoints( aPointsFile, helper, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap, enforcedVertices) &&
1206 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1208 // Write aSmdsToGhs3dIdMap to temp file
1209 TCollection_AsciiString aSmdsToGhs3dIdMapFileName;
1210 aSmdsToGhs3dIdMapFileName = aGenericName + ".ids"; // ids relation
1211 ofstream aIdsFile ( aSmdsToGhs3dIdMapFileName.ToCString() , ios::out);
1213 aIdsFile.rdbuf()->is_open();
1215 INFOS( "Can't write into " << aSmdsToGhs3dIdMapFileName);
1216 return error(SMESH_Comment("Can't write into ") << aSmdsToGhs3dIdMapFileName);
1218 aIdsFile << "Smds Ghs3d" << std::endl;
1219 map <int,int>::const_iterator myit;
1220 for (myit=aSmdsToGhs3dIdMap.begin() ; myit != aSmdsToGhs3dIdMap.end() ; ++myit) {
1221 aIdsFile << myit->first << " " << myit->second << std::endl;
1225 aPointsFile.close();
1229 if ( !_keepFiles ) {
1230 removeFile( aFacesFileName );
1231 removeFile( aPointsFileName );
1232 removeFile( aSmdsToGhs3dIdMapFileName );
1234 return error(COMPERR_BAD_INPUT_MESH);
1236 removeFile( aResultFileName ); // needed for boundary recovery module usage
1238 // -----------------
1240 // -----------------
1242 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1243 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1244 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1246 std::cout << std::endl;
1247 std::cout << "Ghs3d execution..." << std::endl;
1248 std::cout << cmd << std::endl;
1250 system( cmd.ToCString() ); // run
1252 std::cout << std::endl;
1253 std::cout << "End of Ghs3d execution !" << std::endl;
1259 // Mapping the result file
1262 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1263 if ( fileOpen < 0 ) {
1264 std::cout << std::endl;
1265 std::cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << std::endl;
1266 std::cout << "Log: " << aLogFileName << std::endl;
1271 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1272 Ok = readResultFile( fileOpen,
1274 aResultFileName.ToCString(),
1276 theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1277 toMeshHoles, nbEnforcedVertices );
1280 // ---------------------
1281 // remove working files
1282 // ---------------------
1287 removeFile( aLogFileName );
1289 else if ( OSD_File( aLogFileName ).Size() > 0 )
1291 // get problem description from the log file
1292 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1293 storeErrorDescription( aLogFileName, conv );
1297 // the log file is empty
1298 removeFile( aLogFileName );
1299 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1300 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1303 if ( !_keepFiles ) {
1304 removeFile( aFacesFileName );
1305 removeFile( aPointsFileName );
1306 removeFile( aResultFileName );
1307 removeFile( aBadResFileName );
1308 removeFile( aBbResFileName );
1309 removeFile( aSmdsToGhs3dIdMapFileName );
1311 std::cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1313 std::cout << "not ";
1314 std::cout << "treated !" << std::endl;
1315 std::cout << std::endl;
1317 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1324 //=============================================================================
1326 *Here we are going to use the GHS3D mesher w/o geometry
1328 //=============================================================================
1329 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1330 SMESH_MesherHelper* aHelper)
1332 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1334 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1335 TopoDS_Shape theShape = aHelper->GetSubShape();
1337 // a unique working file name
1338 // to avoid access to the same files by eg different users
1339 TCollection_AsciiString aGenericName
1340 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1342 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1343 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1344 aFacesFileName = aGenericName + ".faces"; // in faces
1345 aPointsFileName = aGenericName + ".points"; // in points
1346 aResultFileName = aGenericName + ".noboite";// out points and volumes
1347 aBadResFileName = aGenericName + ".boite"; // out bad result
1348 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1349 aLogFileName = aGenericName + ".log"; // log
1351 // -----------------
1353 // -----------------
1355 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1356 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1358 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1361 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1363 GHS3DPlugin_Hypothesis::TEnforcedVertexValues enforcedVertices;
1364 int nbEnforcedVertices = 0;
1366 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1367 nbEnforcedVertices = enforcedVertices.size();
1372 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1374 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1375 writePoints( aPointsFile, meshDS, aNodeByGhs3dId,enforcedVertices));
1378 aPointsFile.close();
1381 if ( !_keepFiles ) {
1382 removeFile( aFacesFileName );
1383 removeFile( aPointsFileName );
1385 return error(COMPERR_BAD_INPUT_MESH);
1387 removeFile( aResultFileName ); // needed for boundary recovery module usage
1389 // -----------------
1391 // -----------------
1393 TCollection_AsciiString cmd =
1394 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1395 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1396 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1398 system( cmd.ToCString() ); // run
1404 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1405 if ( fileOpen < 0 ) {
1406 std::cout << std::endl;
1407 std::cout << "Error when opening the " << aResultFileName.ToCString() << " file" << std::endl;
1408 std::cout << "Log: " << aLogFileName << std::endl;
1409 std::cout << std::endl;
1413 Ok = readResultFile( fileOpen,
1415 aResultFileName.ToCString(),
1417 meshDS, theShape ,aNodeByGhs3dId, 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( aNodeByGhs3dId );
1433 storeErrorDescription( aLogFileName, conv );
1436 // the log file is empty
1437 removeFile( aLogFileName );
1438 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1439 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1444 removeFile( aFacesFileName );
1445 removeFile( aPointsFileName );
1446 removeFile( aResultFileName );
1447 removeFile( aBadResFileName );
1448 removeFile( aBbResFileName );
1454 //================================================================================
1456 * \brief Provide human readable text by error code reported by ghs3d
1458 //================================================================================
1460 static string translateError(const int errNum)
1464 return "The surface mesh includes a face of type other than edge, "
1465 "triangle or quadrilateral. This face type is not supported.";
1467 return "Not enough memory for the face table.";
1469 return "Not enough memory.";
1471 return "Not enough memory.";
1473 return "Face is ignored.";
1475 return "End of file. Some data are missing in the file.";
1477 return "Read error on the file. There are wrong data in the file.";
1479 return "the metric file is inadequate (dimension other than 3).";
1481 return "the metric file is inadequate (values not per vertices).";
1483 return "the metric file contains more than one field.";
1485 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1486 "value of number of mesh vertices in the \".noboite\" file.";
1488 return "Too many sub-domains.";
1490 return "the number of vertices is negative or null.";
1492 return "the number of faces is negative or null.";
1494 return "A face has a null vertex.";
1496 return "incompatible data.";
1498 return "the number of vertices is negative or null.";
1500 return "the number of vertices is negative or null (in the \".mesh\" file).";
1502 return "the number of faces is negative or null.";
1504 return "A face appears more than once in the input surface mesh.";
1506 return "An edge appears more than once in the input surface mesh.";
1508 return "A face has a vertex negative or null.";
1510 return "NOT ENOUGH MEMORY.";
1512 return "Not enough available memory.";
1514 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1515 "in terms of quality or the input list of points is wrong.";
1517 return "Some vertices are too close to one another or coincident.";
1519 return "Some vertices are too close to one another or coincident.";
1521 return "A vertex cannot be inserted.";
1523 return "There are at least two points considered as coincident.";
1525 return "Some vertices are too close to one another or coincident.";
1527 return "The surface mesh regeneration step has failed.";
1529 return "Constrained edge cannot be enforced.";
1531 return "Constrained face cannot be enforced.";
1533 return "Missing faces.";
1535 return "No guess to start the definition of the connected component(s).";
1537 return "The surface mesh includes at least one hole. The domain is not well defined.";
1539 return "Impossible to define a component.";
1541 return "The surface edge intersects another surface edge.";
1543 return "The surface edge intersects the surface face.";
1545 return "One boundary point lies within a surface face.";
1547 return "One surface edge intersects a surface face.";
1549 return "One boundary point lies within a surface edge.";
1551 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1552 "to too many swaps.";
1554 return "Edge is unique (i.e., bounds a hole in the surface).";
1556 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1558 return "Too many components, too many sub-domain.";
1560 return "The surface mesh includes at least one hole. "
1561 "Therefore there is no domain properly defined.";
1563 return "Statistics.";
1565 return "Statistics.";
1567 return "Warning, it is dramatically tedious to enforce the boundary items.";
1569 return "Not enough memory at this time, nevertheless, the program continues. "
1570 "The expected mesh will be correct but not really as large as required.";
1572 return "see above error code, resulting quality may be poor.";
1574 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1576 return "Unknown face type.";
1579 return "End of file. Some data are missing in the file.";
1581 return "A too small volume element is detected.";
1583 return "There exists at least a null or negative volume element.";
1585 return "There exist null or negative volume elements.";
1587 return "A too small volume element is detected. A face is considered being degenerated.";
1589 return "Some element is suspected to be very bad shaped or wrong.";
1591 return "A too bad quality face is detected. This face is considered degenerated.";
1593 return "A too bad quality face is detected. This face is degenerated.";
1595 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1597 return "Abnormal error occured, contact hotline.";
1599 return "Not enough memory for the face table.";
1601 return "The algorithm cannot run further. "
1602 "The surface mesh is probably very bad in terms of quality.";
1604 return "Bad vertex number.";
1609 //================================================================================
1611 * \brief Retrieve from a string given number of integers
1613 //================================================================================
1615 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1618 ids.reserve( nbIds );
1621 while ( !isdigit( *ptr )) ++ptr;
1622 if ( ptr[-1] == '-' ) --ptr;
1623 ids.push_back( strtol( ptr, &ptr, 10 ));
1629 //================================================================================
1631 * \brief Retrieve problem description form a log file
1632 * \retval bool - always false
1634 //================================================================================
1636 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1637 const _Ghs2smdsConvertor & toSmdsConvertor )
1641 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1643 int file = ::open (logFile.ToCString(), O_RDONLY);
1646 return error( SMESH_Comment("See ") << logFile << " for problem description");
1649 // struct stat status;
1650 // fstat(file, &status);
1651 // size_t length = status.st_size;
1652 off_t length = lseek( file, 0, SEEK_END);
1653 lseek( file, 0, SEEK_SET);
1656 vector< char > buf( length );
1657 int nBytesRead = ::read (file, & buf[0], length);
1659 char* ptr = & buf[0];
1660 char* bufEnd = ptr + nBytesRead;
1662 SMESH_Comment errDescription;
1664 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1666 // look for errors "ERR #"
1668 set<string> foundErrorStr; // to avoid reporting same error several times
1669 set<int> elemErrorNums; // not to report different types of errors with bad elements
1670 while ( ++ptr < bufEnd )
1672 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1675 list<const SMDS_MeshElement*> badElems;
1676 vector<int> nodeIds;
1680 int errNum = strtol(ptr, &ptr, 10);
1681 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1683 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1684 ptr = getIds(ptr, NODE, nodeIds);
1685 ptr = getIds(ptr, TRIA, nodeIds);
1686 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1688 case 1000: // ERR 1000 : 1 3 2
1689 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1690 ptr = getIds(ptr, TRIA, nodeIds);
1691 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1694 // Edge (e1, e2) appears more than once in the input surface mesh
1695 ptr = getIds(ptr, EDGE, nodeIds);
1696 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1699 // Face (f 1, f 2, f 3) has a vertex negative or null
1700 ptr = getIds(ptr, TRIA, nodeIds);
1701 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1704 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1705 ptr = getIds(ptr, NODE, nodeIds);
1706 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1707 ptr = getIds(ptr, NODE, nodeIds);
1708 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1711 // Vertex v1 cannot be inserted (warning).
1712 ptr = getIds(ptr, NODE, nodeIds);
1713 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1716 // There are at least two points whose distance is dist, i.e., considered as coincident
1717 case 2103: // ERR 2103 : 16 WITH 3
1718 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1719 ptr = getIds(ptr, NODE, nodeIds);
1720 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1721 ptr = getIds(ptr, NODE, nodeIds);
1722 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1725 // Constrained edge (e1, e2) cannot be enforced (warning).
1726 ptr = getIds(ptr, EDGE, nodeIds);
1727 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1730 // Constrained face (f 1, f 2, f 3) cannot be enforced
1731 ptr = getIds(ptr, TRIA, nodeIds);
1732 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1734 case 3103: // ERR 3103 : 1 2 WITH 7 3
1735 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1736 ptr = getIds(ptr, EDGE, nodeIds);
1737 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1738 ptr = getIds(ptr, EDGE, nodeIds);
1739 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1741 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1742 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1743 ptr = getIds(ptr, EDGE, nodeIds);
1744 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1745 ptr = getIds(ptr, TRIA, nodeIds);
1746 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1748 case 3105: // ERR 3105 : 8 IN 2 3 5
1749 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1750 ptr = getIds(ptr, NODE, nodeIds);
1751 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1752 ptr = getIds(ptr, TRIA, nodeIds);
1753 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1756 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1757 ptr = getIds(ptr, EDGE, nodeIds);
1758 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1759 ptr = getIds(ptr, TRIA, nodeIds);
1760 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1762 case 3107: // ERR 3107 : 2 IN 4 1
1763 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1764 ptr = getIds(ptr, NODE, nodeIds);
1765 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1766 ptr = getIds(ptr, EDGE, nodeIds);
1767 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1769 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1770 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1771 ptr = getIds(ptr, EDGE, nodeIds);
1772 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1774 case 9000: // ERR 9000
1775 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1776 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1777 // A too small volume element is detected. Are reported the index of the element,
1778 // its four vertex indices, its volume and the tolerance threshold value
1779 ptr = getIds(ptr, ID, nodeIds);
1780 ptr = getIds(ptr, VOL, nodeIds);
1781 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1782 // even if all nodes found, volume it most probably invisible,
1783 // add its faces to demenstrate it anyhow
1785 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1786 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1787 faceNodes[2] = nodeIds[3]; // 013
1788 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1789 faceNodes[1] = nodeIds[2]; // 023
1790 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1791 faceNodes[0] = nodeIds[1]; // 123
1792 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1795 case 9001: // ERR 9001
1796 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1797 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1798 // %% NUMBER OF NULL VOLUME TETS : 0
1799 // There exists at least a null or negative volume element
1802 // There exist n null or negative volume elements
1805 // A too small volume element is detected
1808 // A too bad quality face is detected. This face is considered degenerated,
1809 // its index, its three vertex indices together with its quality value are reported
1810 break; // same as next
1811 case 9112: // ERR 9112
1812 // FACE 2 WITH VERTICES : 4 2 5
1813 // SMALL INRADIUS : 0.
1814 // A too bad quality face is detected. This face is degenerated,
1815 // its index, its three vertex indices together with its inradius are reported
1816 ptr = getIds(ptr, ID, nodeIds);
1817 ptr = getIds(ptr, TRIA, nodeIds);
1818 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1819 // add triangle edges as it most probably has zero area and hence invisible
1821 vector<int> edgeNodes(2);
1822 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1823 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1824 edgeNodes[1] = nodeIds[2]; // 0-2
1825 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1826 edgeNodes[0] = nodeIds[1]; // 1-2
1827 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1832 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1834 continue; // not to report same error several times
1836 // const SMDS_MeshElement* nullElem = 0;
1837 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1839 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1840 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1841 // if ( oneMoreErrorType )
1842 // continue; // not to report different types of errors with bad elements
1845 // store bad elements
1846 //if ( allElemsOk ) {
1847 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1848 for ( ; elem != badElems.end(); ++elem )
1849 addBadInputElement( *elem );
1853 string text = translateError( errNum );
1854 if ( errDescription.find( text ) == text.npos ) {
1855 if ( !errDescription.empty() )
1856 errDescription << "\n";
1857 errDescription << text;
1862 if ( errDescription.empty() ) { // no errors found
1863 char msg[] = "connection to server failed";
1864 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1865 errDescription << "Licence problems.";
1868 if ( errDescription.empty() )
1869 errDescription << "See " << logFile << " for problem description";
1871 errDescription << "\nSee " << logFile << " for more information";
1873 return error( errDescription );
1876 //================================================================================
1878 * \brief Creates _Ghs2smdsConvertor
1880 //================================================================================
1882 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1883 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1887 //================================================================================
1889 * \brief Creates _Ghs2smdsConvertor
1891 //================================================================================
1893 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1894 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1898 //================================================================================
1900 * \brief Return SMDS element by ids of GHS3D nodes
1902 //================================================================================
1904 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1906 size_t nbNodes = ghsNodes.size();
1907 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1908 for ( size_t i = 0; i < nbNodes; ++i ) {
1909 int ghsNode = ghsNodes[ i ];
1910 if ( _ghs2NodeMap ) {
1911 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1912 if ( in == _ghs2NodeMap->end() )
1914 nodes[ i ] = in->second;
1917 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1919 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1925 if ( nbNodes == 2 ) {
1926 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1928 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1931 if ( nbNodes == 3 ) {
1932 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1934 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1938 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1944 //=============================================================================
1948 //=============================================================================
1949 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
1950 const TopoDS_Shape& aShape,
1951 MapShapeNbElems& aResMap)
1953 int nbtri = 0, nbqua = 0;
1954 double fullArea = 0.0;
1955 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1956 TopoDS_Face F = TopoDS::Face( exp.Current() );
1957 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1958 MapShapeNbElemsItr anIt = aResMap.find(sm);
1959 if( anIt==aResMap.end() ) {
1960 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1961 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1962 "Submesh can not be evaluated",this));
1965 std::vector<int> aVec = (*anIt).second;
1966 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
1967 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1969 BRepGProp::SurfaceProperties(F,G);
1970 double anArea = G.Mass();
1974 // collect info from edges
1975 int nb0d_e = 0, nb1d_e = 0;
1976 bool IsQuadratic = false;
1977 bool IsFirst = true;
1978 TopTools_MapOfShape tmpMap;
1979 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1980 TopoDS_Edge E = TopoDS::Edge(exp.Current());
1981 if( tmpMap.Contains(E) )
1984 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
1985 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
1986 std::vector<int> aVec = (*anIt).second;
1987 nb0d_e += aVec[SMDSEntity_Node];
1988 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
1990 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1996 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
1999 BRepGProp::VolumeProperties(aShape,G);
2000 double aVolume = G.Mass();
2001 double tetrVol = 0.1179*ELen*ELen*ELen;
2002 double CoeffQuality = 0.9;
2003 int nbVols = int(aVolume/tetrVol/CoeffQuality);
2004 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2005 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2006 std::vector<int> aVec(SMDSEntity_Last);
2007 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2009 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2010 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2011 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2014 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2015 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2016 aVec[SMDSEntity_Pyramid] = nbqua;
2018 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2019 aResMap.insert(std::make_pair(sm,aVec));