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 "SMESH_Gen.hxx"
32 #include "SMESH_Mesh.hxx"
33 #include "SMESH_Comment.hxx"
34 #include "SMESH_MesherHelper.hxx"
35 #include "SMESH_MeshEditor.hxx"
37 #include "SMDS_MeshElement.hxx"
38 #include "SMDS_MeshNode.hxx"
39 #include "SMDS_FaceOfNodes.hxx"
40 #include "SMDS_VolumeOfNodes.hxx"
42 #include <BRepAdaptor_Surface.hxx>
43 #include <BRepBndLib.hxx>
44 #include <BRepClass3d_SolidClassifier.hxx>
45 #include <BRepTools.hxx>
46 #include <BRep_Tool.hxx>
47 #include <Bnd_Box.hxx>
48 #include <GeomAPI_ProjectPointOnSurf.hxx>
49 #include <OSD_File.hxx>
50 #include <Precision.hxx>
51 #include <Quantity_Parameter.hxx>
52 #include <Standard_ProgramError.hxx>
53 #include <Standard_ErrorHandler.hxx>
54 #include <Standard_Failure.hxx>
56 #include <TopExp_Explorer.hxx>
57 #include <TopTools_IndexedMapOfShape.hxx>
58 #include <TopTools_ListIteratorOfListOfShape.hxx>
60 //#include <BRepClass_FaceClassifier.hxx>
61 #include <TopTools_MapOfShape.hxx>
62 #include <BRepGProp.hxx>
63 #include <GProp_GProps.hxx>
65 #include "utilities.h"
70 #include <sys/sysinfo.h>
75 //#include <Standard_Stream.hxx>
78 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
99 static void removeFile( const TCollection_AsciiString& fileName )
102 OSD_File( fileName ).Remove();
104 catch ( Standard_ProgramError ) {
105 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
109 //=============================================================================
113 //=============================================================================
115 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
116 : SMESH_3D_Algo(hypId, studyId, gen)
118 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
120 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
121 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
124 _compatibleHypothesis.push_back("GHS3D_Parameters");
125 _requireShape = false; // can work without shape
128 //=============================================================================
132 //=============================================================================
134 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
136 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
139 //=============================================================================
143 //=============================================================================
145 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
146 const TopoDS_Shape& aShape,
147 Hypothesis_Status& aStatus )
149 aStatus = SMESH_Hypothesis::HYP_OK;
151 // there is only one compatible Hypothesis so far
154 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
156 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
158 _keepFiles = _hyp->GetKeepFiles();
163 //=======================================================================
164 //function : findShape
166 //=======================================================================
168 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
170 const TopoDS_Shape shape[],
173 TopAbs_State * state = 0)
176 int j, iShape, nbNode = 4;
178 for ( j=0; j<nbNode; j++ ) {
179 gp_XYZ p ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
180 if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) {
187 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
188 if (state) *state = SC.State();
189 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
190 for (iShape = 0; iShape < nShape; iShape++) {
191 aShape = shape[iShape];
192 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
193 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
194 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
195 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
196 if (state) *state = SC.State();
197 if (SC.State() == TopAbs_IN)
205 //=======================================================================
206 //function : readMapIntLine
208 //=======================================================================
210 static char* readMapIntLine(char* ptr, int tab[]) {
212 std::cout << std::endl;
214 for ( int i=0; i<17; i++ ) {
215 intVal = strtol(ptr, &ptr, 10);
222 //=======================================================================
223 //function : countShape
225 //=======================================================================
227 template < class Mesh, class Shape >
228 static int countShape( Mesh* mesh, Shape shape ) {
229 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
231 for ( ; expShape.More(); expShape.Next() ) {
237 //=======================================================================
238 //function : writeFaces
240 //=======================================================================
242 static bool writeFaces (ofstream & theFile,
243 SMESHDS_Mesh * theMesh,
244 const map <int,int> & theSmdsToGhs3dIdMap)
248 // NB_ELEMS DUMMY_INT
249 // Loop from 1 to NB_ELEMS
250 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
252 int nbShape = countShape( theMesh, TopAbs_FACE );
254 int *tabID; tabID = new int[nbShape];
255 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
257 SMESHDS_SubMesh* theSubMesh;
258 const SMDS_MeshElement* aFace;
259 const char* space = " ";
260 const int dummyint = 0;
261 map<int,int>::const_iterator itOnMap;
262 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
263 int shapeID, nbNodes, aSmdsID;
266 std::cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << std::endl;
267 std::cout << std::endl;
269 theFile << space << theMesh->NbFaces() << space << dummyint << std::endl;
271 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
272 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
274 aShape = expface.Current();
275 shapeID = theMesh->ShapeToIndex( aShape );
277 for ( int j=0; j<=i; j++) {
278 if ( shapeID == tabID[j] ) {
285 tabShape[i] = aShape;
288 for ( int i =0; i < nbShape; i++ ) {
289 if ( tabID[i] != 0 ) {
290 aShape = tabShape[i];
292 theSubMesh = theMesh->MeshElements( aShape );
293 if ( !theSubMesh ) continue;
294 itOnSubMesh = theSubMesh->GetElements();
295 while ( itOnSubMesh->more() ) {
296 aFace = itOnSubMesh->next();
297 nbNodes = aFace->NbNodes();
299 theFile << space << nbNodes;
301 itOnSubFace = aFace->nodesIterator();
302 while ( itOnSubFace->more() ) {
304 aSmdsID = itOnSubFace->next()->GetID();
305 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
306 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
308 theFile << space << (*itOnMap).second;
311 // (NB_NODES + 1) times: DUMMY_INT
312 for ( int j=0; j<=nbNodes; j++)
313 theFile << space << dummyint;
315 theFile << std::endl;
327 //=======================================================================
328 //function : writeFaces
329 //purpose : Write Faces in case if generate 3D mesh w/o geometry
330 //=======================================================================
332 static bool writeFaces (ofstream & theFile,
333 SMESHDS_Mesh * theMesh,
334 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
338 // NB_ELEMS DUMMY_INT
339 // Loop from 1 to NB_ELEMS
340 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
344 list< const SMDS_MeshElement* > faces;
345 list< const SMDS_MeshElement* >::iterator f;
346 map< const SMDS_MeshNode*,int >::iterator it;
347 SMDS_ElemIteratorPtr nodeIt;
348 const SMDS_MeshElement* elem;
351 const char* space = " ";
352 const int dummyint = 0;
354 //get all faces from mesh
355 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
356 while ( eIt->more() ) {
357 const SMDS_MeshElement* elem = eIt->next();
360 faces.push_back( elem );
367 std::cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
369 // NB_ELEMS DUMMY_INT
370 theFile << space << nbFaces << space << dummyint << std::endl;
372 // Loop from 1 to NB_ELEMS
374 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
376 for ( ; f != faces.end(); ++f )
380 nbNodes = elem->NbNodes();
381 theFile << space << nbNodes;
383 // NODE_NB_1 NODE_NB_2 ...
384 nodeIt = elem->nodesIterator();
385 while ( nodeIt->more() )
388 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
389 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
390 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
391 theFile << space << it->second;
394 // (NB_NODES + 1) times: DUMMY_INT
395 for ( int i=0; i<=nbNodes; i++)
396 theFile << space << dummyint;
398 theFile << std::endl;
401 // put nodes to theNodeByGhs3dId vector
402 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
403 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
404 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
406 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
412 //=======================================================================
413 //function : writePoints
415 //=======================================================================
417 static bool writePoints (ofstream & theFile,
418 SMESHDS_Mesh * theMesh,
419 map <int,int> & theSmdsToGhs3dIdMap,
420 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap,
421 map<vector<double>,double> & theEnforcedVertices)
426 // Loop from 1 to NB_NODES
429 int nbNodes = theMesh->NbNodes();
432 int nbEnforcedVertices = theEnforcedVertices.size();
434 const char* space = " ";
435 const int dummyint = 0;
438 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
439 const SMDS_MeshNode* node;
442 std::cout << std::endl;
443 std::cout << "The initial 2D mesh contains :" << std::endl;
444 std::cout << " " << nbNodes << " nodes" << std::endl;
445 if (nbEnforcedVertices > 0) {
446 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
448 std::cout << std::endl;
449 std::cout << "Start writing in 'points' file ..." << std::endl;
450 theFile << space << nbNodes << std::endl;
452 // Loop from 1 to NB_NODES
457 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
458 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
463 << space << node->X()
464 << space << node->Y()
465 << space << node->Z()
466 << space << dummyint;
468 theFile << std::endl;
472 // Iterate over the enforced vertices
473 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
474 const TopoDS_Shape shapeToMesh = theMesh->ShapeToMesh();
475 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
476 double x = vertexIt->first[0];
477 double y = vertexIt->first[1];
478 double z = vertexIt->first[2];
479 // Test if point is inside shape to mesh
480 gp_Pnt myPoint(x,y,z);
481 BRepClass3d_SolidClassifier scl(shapeToMesh);
482 scl.Perform(myPoint, 1e-7);
483 TopAbs_State result = scl.State();
484 if ( result == TopAbs_IN ) {
485 MESSAGE("Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second);
486 // X Y Z PHY_SIZE DUMMY_INT
491 << space << vertexIt->second
492 << space << dummyint;
494 theFile << std::endl;
497 MESSAGE("Enforced vertex (" << x << "," << y <<"," << z << ") is not inside the geometry: it was not added ");
501 std::cout << std::endl;
502 std::cout << "End writing in 'points' file." << std::endl;
507 //=======================================================================
508 //function : writePoints
510 //=======================================================================
512 static bool writePoints (ofstream & theFile,
513 SMESHDS_Mesh * theMesh,
514 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId,
515 map<vector<double>,double> & theEnforcedVertices)
520 // Loop from 1 to NB_NODES
523 //int nbNodes = theMesh->NbNodes();
524 int nbNodes = theNodeByGhs3dId.size();
528 int nbEnforcedVertices = theEnforcedVertices.size();
530 const char* space = " ";
531 const int dummyint = 0;
533 const SMDS_MeshNode* node;
536 std::cout << std::endl;
537 std::cout << "The initial 2D mesh contains :" << std::endl;
538 std::cout << " " << nbNodes << " nodes" << std::endl;
539 std::cout << " " << nbEnforcedVertices << " enforced vertices" << std::endl;
540 std::cout << std::endl;
541 std::cout << "Start writing in 'points' file ..." << std::endl;
542 theFile << space << nbNodes << std::endl;
544 // Loop from 1 to NB_NODES
546 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
547 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
548 for ( ; nodeIt != after; ++nodeIt )
554 << space << node->X()
555 << space << node->Y()
556 << space << node->Z()
557 << space << dummyint;
559 theFile << std::endl;
563 // Iterate over the enforced vertices
564 GHS3DPlugin_Hypothesis::TEnforcedVertexValues::const_iterator vertexIt;
565 const TopoDS_Shape shapeToMesh = theMesh->ShapeToMesh();
566 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
567 double x = vertexIt->first[0];
568 double y = vertexIt->first[1];
569 double z = vertexIt->first[2];
570 // Test if point is inside shape to mesh
571 gp_Pnt myPoint(x,y,z);
572 BRepClass3d_SolidClassifier scl(shapeToMesh);
573 scl.Perform(myPoint, 1e-7);
574 TopAbs_State result = scl.State();
575 if ( result == TopAbs_IN ) {
576 std::cout << "Adding enforced vertex (" << x << "," << y <<"," << z << ") = " << vertexIt->second << std::endl;
578 // X Y Z PHY_SIZE DUMMY_INT
583 << space << vertexIt->second
584 << space << dummyint;
586 theFile << std::endl;
589 std::cout << std::endl;
590 std::cout << "End writing in 'points' file." << std::endl;
595 //=======================================================================
596 //function : findShapeID
597 //purpose : find the solid corresponding to GHS3D sub-domain following
598 // the technique proposed in GHS3D manual in chapter
599 // "B.4 Subdomain (sub-region) assignment"
600 //=======================================================================
602 static int findShapeID(SMESH_Mesh& mesh,
603 const SMDS_MeshNode* node1,
604 const SMDS_MeshNode* node2,
605 const SMDS_MeshNode* node3,
606 const bool toMeshHoles)
608 const int invalidID = 0;
609 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
611 // face the nodes belong to
612 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
616 // geom face the face assigned to
617 SMESH_MeshEditor editor(&mesh);
618 int geomFaceID = editor.FindShape( face );
621 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
622 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
624 TopoDS_Face geomFace = TopoDS::Face( shape );
626 // solids bounded by geom face
627 TopTools_IndexedMapOfShape solids, shells;
628 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
629 for ( ; ansIt.More(); ansIt.Next() ) {
630 switch ( ansIt.Value().ShapeType() ) {
632 solids.Add( ansIt.Value() ); break;
634 shells.Add( ansIt.Value() ); break;
638 // analyse found solids
639 if ( solids.Extent() == 0 || shells.Extent() == 0)
642 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
643 if ( solids.Extent() == 1 )
646 return meshDS->ShapeToIndex( solid1 );
648 // - Are we at a hole boundary face?
649 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
650 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
652 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
653 // check if any edge of shells(1) belongs to another shell
654 for ( ; eExp.More() && !touch; eExp.Next() ) {
655 ansIt = mesh.GetAncestors( eExp.Current() );
656 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
657 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
658 touch = ( !ansIt.Value().IsSame( shells(1) ));
662 return meshDS->ShapeToIndex( solid1 );
665 // find orientation of geom face within the first solid
666 TopExp_Explorer fExp( solid1, TopAbs_FACE );
667 for ( ; fExp.More(); fExp.Next() )
668 if ( geomFace.IsSame( fExp.Current() )) {
669 geomFace = TopoDS::Face( fExp.Current() );
673 return invalidID; // face not found
675 // normale to triangle
676 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
677 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
678 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
679 gp_Vec vec12( node1Pnt, node2Pnt );
680 gp_Vec vec13( node1Pnt, node3Pnt );
681 gp_Vec meshNormal = vec12 ^ vec13;
682 if ( meshNormal.SquareMagnitude() < DBL_MIN )
685 // find UV of node1 on geomFace
686 SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
687 const SMDS_MeshNode* nNotOnSeamEdge = 0;
688 if ( helper.IsSeamShape( node1->GetPosition()->GetShapeId() ))
689 if ( helper.IsSeamShape( node2->GetPosition()->GetShapeId() ))
690 nNotOnSeamEdge = node3;
692 nNotOnSeamEdge = node2;
694 gp_XY uv = helper.GetNodeUV( geomFace, node1, nNotOnSeamEdge, &uvOK );
695 // check that uv is correct
697 TopoDS_Shape nodeShape = helper.GetSubShapeByNode( node1, meshDS );
698 if ( !nodeShape.IsNull() )
699 switch ( nodeShape.ShapeType() )
701 case TopAbs_FACE: tol = BRep_Tool::Tolerance( TopoDS::Face( nodeShape )); break;
702 case TopAbs_EDGE: tol = BRep_Tool::Tolerance( TopoDS::Edge( nodeShape )); break;
703 case TopAbs_VERTEX: tol = BRep_Tool::Tolerance( TopoDS::Vertex( nodeShape )); break;
706 BRepAdaptor_Surface surface( geomFace );
707 if ( !uvOK || node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol )
710 // normale to geomFace at UV
712 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
713 gp_Vec geomNormal = du ^ dv;
714 if ( geomNormal.SquareMagnitude() < DBL_MIN )
715 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
716 if ( geomFace.Orientation() == TopAbs_REVERSED )
717 geomNormal.Reverse();
720 bool isReverse = ( meshNormal * geomNormal ) < 0;
722 return meshDS->ShapeToIndex( solid1 );
724 if ( solids.Extent() == 1 )
725 return HOLE_ID; // we are inside a hole
727 return meshDS->ShapeToIndex( solids(2) );
730 //=======================================================================
731 //function : readResultFile
733 //=======================================================================
735 static bool readResultFile(const int fileOpen,
737 const char* fileName,
740 TopoDS_Shape tabShape[],
743 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
745 int nbEnforcedVertices)
747 MESSAGE("GHS3DPlugin_GHS3D::readResultFile()");
755 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
758 int nbElems, nbNodes, nbInputNodes;
759 int nodeId/*, triangleId*/;
761 int ID, shapeID, ghs3dShapeID;
764 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
766 int *tab, *tabID, *nodeID, *nodeAssigne;
768 const SMDS_MeshNode **node;
771 //tabID = new int[nbShape];
773 coord = new double[3];
774 node = new const SMDS_MeshNode*[4];
777 SMDS_MeshNode * aNewNode;
778 map <int,const SMDS_MeshNode*>::iterator itOnNode;
779 SMDS_MeshElement* aTet;
784 // Read the file state
785 fileStat = fstat(fileOpen, &status);
786 length = status.st_size;
788 // Mapping the result file into memory
790 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
791 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
792 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
793 0, (DWORD)length, NULL);
794 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
796 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
800 ptr = readMapIntLine(ptr, tab);
805 nbInputNodes = tab[2];
807 nodeAssigne = new int[ nbNodes+1 ];
810 aSolid = tabShape[0];
812 // Reading the nodeId
813 for (int i=0; i < 4*nbElems; i++)
814 nodeId = strtol(ptr, &ptr, 10);
816 MESSAGE("nbInputNodes: "<<nbInputNodes);
817 MESSAGE("nbEnforcedVertices: "<<nbEnforcedVertices);
818 // Reading the nodeCoor and update the nodeMap
819 for (int iNode=1; iNode <= nbNodes; iNode++) {
820 for (int iCoor=0; iCoor < 3; iCoor++)
821 coord[ iCoor ] = strtod(ptr, &ptr);
822 nodeAssigne[ iNode ] = 1;
823 if ( iNode > (nbInputNodes-nbEnforcedVertices) ) {
824 // Creating SMESH nodes
825 // - for enforced vertices
826 // - for vertices of forced edges
828 nodeAssigne[ iNode ] = 0;
829 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
830 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
834 // Reading the number of triangles which corresponds to the number of sub-domains
835 nbTriangle = strtol(ptr, &ptr, 10);
837 tabID = new int[nbTriangle];
838 for (int i=0; i < nbTriangle; i++) {
840 // find the solid corresponding to GHS3D sub-domain following
841 // the technique proposed in GHS3D manual in chapter
842 // "B.4 Subdomain (sub-region) assignment"
843 int nodeId1 = strtol(ptr, &ptr, 10);
844 int nodeId2 = strtol(ptr, &ptr, 10);
845 int nodeId3 = strtol(ptr, &ptr, 10);
846 if ( nbTriangle > 1 ) {
847 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
848 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
849 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
852 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
853 // -- 0020330: Pb with ghs3d as a submesh
854 // check that found shape is to be meshed
855 if ( tabID[i] > 0 ) {
856 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] );
857 bool isToBeMeshed = false;
858 for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS )
859 isToBeMeshed = foundShape.IsSame( tabShape[ iS ]);
863 // END -- 0020330: Pb with ghs3d as a submesh
865 std::cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << std::endl;
867 } catch ( Standard_Failure ) {
874 if ( nbTriangle <= nbShape ) // no holes
875 toMeshHoles = true; // not avoid creating tetras in holes
877 // Associating the tetrahedrons to the shapes
878 shapeID = compoundID;
879 for (int iElem = 0; iElem < nbElems; iElem++) {
880 for (int iNode = 0; iNode < 4; iNode++) {
881 ID = strtol(tetraPtr, &tetraPtr, 10);
882 itOnNode = theGhs3dIdToNodeMap.find(ID);
883 node[ iNode ] = itOnNode->second;
884 nodeID[ iNode ] = ID;
886 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
887 // tetras within holes depending on hypo option,
888 // so we first check if aTet is inside a hole and then create it
889 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
890 if ( nbTriangle > 1 ) {
891 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
892 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
893 if ( tabID[ ghs3dShapeID ] == 0 ) {
895 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
896 if ( toMeshHoles || state == TopAbs_IN )
897 shapeID = theMeshDS->ShapeToIndex( aSolid );
898 tabID[ ghs3dShapeID ] = shapeID;
901 shapeID = tabID[ ghs3dShapeID ];
903 else if ( nbShape > 1 ) {
904 // Case where nbTriangle == 1 while nbShape == 2 encountered
905 // with compound of 2 boxes and "To mesh holes"==False,
906 // so there are no subdomains specified for each tetrahedron.
907 // Try to guess a solid by a node already bound to shape
909 for ( int i=0; i<4 && shapeID==0; i++ ) {
910 if ( nodeAssigne[ nodeID[i] ] == 1 &&
911 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
912 node[i]->GetPosition()->GetShapeId() > 1 )
914 shapeID = node[i]->GetPosition()->GetShapeId();
918 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
919 shapeID = theMeshDS->ShapeToIndex( aSolid );
922 // set new nodes and tetrahedron onto the shape
923 for ( int i=0; i<4; i++ ) {
924 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
925 if ( shapeID != HOLE_ID )
926 theMeshDS->SetNodeInVolume( node[i], shapeID );
927 nodeAssigne[ nodeID[i] ] = shapeID;
930 if ( toMeshHoles || shapeID != HOLE_ID ) {
931 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
932 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
935 shapeIDs.insert( shapeID );
938 // Remove nodes of tetras inside holes if !toMeshHoles
939 if ( !toMeshHoles ) {
940 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
941 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
942 ID = itOnNode->first;
943 if ( nodeAssigne[ ID ] == HOLE_ID )
944 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
949 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
951 UnmapViewOfFile(mapPtr);
952 CloseHandle(hMapObject);
955 munmap(mapPtr, length);
964 delete [] nodeAssigne;
967 if ( shapeIDs.size() != nbShape ) {
968 std::cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << std::endl;
969 for (int i=0; i<nbShape; i++) {
970 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
971 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
972 std::cout << " Solid #" << shapeID << " not found" << std::endl;
980 //=======================================================================
981 //function : readResultFile
983 //=======================================================================
985 static bool readResultFile(const int fileOpen,
987 const char* fileName,
989 SMESHDS_Mesh* theMeshDS,
991 vector <const SMDS_MeshNode*>& theNodeByGhs3dId,
992 int nbEnforcedVertices) {
1002 int nbElems, nbNodes, nbInputNodes;
1003 int nodeId, triangleId;
1009 const SMDS_MeshNode **node;
1012 coord = new double[3];
1013 node = new const SMDS_MeshNode*[4];
1015 SMDS_MeshNode * aNewNode;
1016 map <int,const SMDS_MeshNode*>::iterator IdNode;
1017 SMDS_MeshElement* aTet;
1019 // Read the file state
1020 fileStat = fstat(fileOpen, &status);
1021 length = status.st_size;
1023 // Mapping the result file into memory
1025 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
1026 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
1027 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
1028 0, (DWORD)length, NULL);
1029 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
1031 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
1035 ptr = readMapIntLine(ptr, tab);
1040 nbInputNodes = tab[2];
1042 theNodeByGhs3dId.resize( nbNodes );
1044 // Reading the nodeId
1045 for (int i=0; i < 4*nbElems; i++)
1046 nodeId = strtol(ptr, &ptr, 10);
1048 // Reading the nodeCoor and update the nodeMap
1049 shapeID = theMeshDS->ShapeToIndex( aSolid );
1050 for (int iNode=0; iNode < nbNodes; iNode++) {
1051 for (int iCoor=0; iCoor < 3; iCoor++)
1052 coord[ iCoor ] = strtod(ptr, &ptr);
1053 if ((iNode+1) > (nbInputNodes-nbEnforcedVertices)) {
1054 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
1055 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
1056 theNodeByGhs3dId[ iNode ] = aNewNode;
1060 // Reading the triangles
1061 nbTriangle = strtol(ptr, &ptr, 10);
1063 for (int i=0; i < 3*nbTriangle; i++)
1064 triangleId = strtol(ptr, &ptr, 10);
1068 // Associating the tetrahedrons to the shapes
1069 for (int iElem = 0; iElem < nbElems; iElem++) {
1070 for (int iNode = 0; iNode < 4; iNode++) {
1071 ID = strtol(tetraPtr, &tetraPtr, 10);
1072 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
1074 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
1075 shapeID = theMeshDS->ShapeToIndex( aSolid );
1076 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
1079 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
1081 UnmapViewOfFile(mapPtr);
1082 CloseHandle(hMapObject);
1085 munmap(mapPtr, length);
1096 //=============================================================================
1098 *Here we are going to use the GHS3D mesher
1100 //=============================================================================
1102 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1103 const TopoDS_Shape& theShape)
1106 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1108 // we count the number of shapes
1109 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
1111 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1112 for ( ; expBox.More(); expBox.Next() )
1115 // create bounding box for every shape inside the compound
1118 TopoDS_Shape* tabShape;
1120 tabShape = new TopoDS_Shape[_nbShape];
1121 tabBox = new double*[_nbShape];
1122 for (int i=0; i<_nbShape; i++)
1123 tabBox[i] = new double[6];
1124 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1126 // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh
1127 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
1128 tabShape[iShape] = expBox.Current();
1129 Bnd_Box BoundingBox;
1130 BRepBndLib::Add(expBox.Current(), BoundingBox);
1131 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1132 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1133 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1134 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1138 // a unique working file name
1139 // to avoid access to the same files by eg different users
1140 TCollection_AsciiString aGenericName
1141 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1143 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1144 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1145 aFacesFileName = aGenericName + ".faces"; // in faces
1146 aPointsFileName = aGenericName + ".points"; // in points
1147 aResultFileName = aGenericName + ".noboite";// out points and volumes
1148 aBadResFileName = aGenericName + ".boite"; // out bad result
1149 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1150 aLogFileName = aGenericName + ".log"; // log
1152 // -----------------
1154 // -----------------
1156 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1157 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1160 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1162 INFOS( "Can't write into " << aFacesFileName);
1163 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1165 map <int,int> aSmdsToGhs3dIdMap;
1166 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1167 map<vector<double>,double> enforcedVertices;
1168 int nbEnforcedVertices = 0;
1170 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1171 nbEnforcedVertices = enforcedVertices.size();
1176 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap, enforcedVertices) &&
1177 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1180 aPointsFile.close();
1183 if ( !_keepFiles ) {
1184 removeFile( aFacesFileName );
1185 removeFile( aPointsFileName );
1187 return error(COMPERR_BAD_INPUT_MESH);
1189 removeFile( aResultFileName ); // needed for boundary recovery module usage
1191 // -----------------
1193 // -----------------
1195 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1196 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1197 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1199 std::cout << std::endl;
1200 std::cout << "Ghs3d execution..." << std::endl;
1201 std::cout << cmd << std::endl;
1203 system( cmd.ToCString() ); // run
1205 std::cout << std::endl;
1206 std::cout << "End of Ghs3d execution !" << std::endl;
1212 // Mapping the result file
1215 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1216 if ( fileOpen < 0 ) {
1217 std::cout << std::endl;
1218 std::cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << std::endl;
1219 std::cout << "Log: " << aLogFileName << std::endl;
1224 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1225 Ok = readResultFile( fileOpen,
1227 aResultFileName.ToCString(),
1229 theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1230 toMeshHoles, nbEnforcedVertices );
1233 // ---------------------
1234 // remove working files
1235 // ---------------------
1240 removeFile( aLogFileName );
1242 else if ( OSD_File( aLogFileName ).Size() > 0 )
1244 // get problem description from the log file
1245 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1246 storeErrorDescription( aLogFileName, conv );
1250 // the log file is empty
1251 removeFile( aLogFileName );
1252 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1253 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1256 if ( !_keepFiles ) {
1257 removeFile( aFacesFileName );
1258 removeFile( aPointsFileName );
1259 removeFile( aResultFileName );
1260 removeFile( aBadResFileName );
1261 removeFile( aBbResFileName );
1263 std::cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1265 std::cout << "not ";
1266 std::cout << "treated !" << std::endl;
1267 std::cout << std::endl;
1269 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1276 //=============================================================================
1278 *Here we are going to use the GHS3D mesher w/o geometry
1280 //=============================================================================
1281 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1282 SMESH_MesherHelper* aHelper)
1284 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1286 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1287 TopoDS_Shape theShape = aHelper->GetSubShape();
1289 // a unique working file name
1290 // to avoid access to the same files by eg different users
1291 TCollection_AsciiString aGenericName
1292 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1294 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1295 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1296 aFacesFileName = aGenericName + ".faces"; // in faces
1297 aPointsFileName = aGenericName + ".points"; // in points
1298 aResultFileName = aGenericName + ".noboite";// out points and volumes
1299 aBadResFileName = aGenericName + ".boite"; // out bad result
1300 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1301 aLogFileName = aGenericName + ".log"; // log
1303 // -----------------
1305 // -----------------
1307 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1308 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1310 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1313 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1315 GHS3DPlugin_Hypothesis::TEnforcedVertexValues enforcedVertices;
1316 int nbEnforcedVertices = 0;
1318 enforcedVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1319 nbEnforcedVertices = enforcedVertices.size();
1324 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1326 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1327 writePoints( aPointsFile, meshDS, aNodeByGhs3dId,enforcedVertices));
1330 aPointsFile.close();
1333 if ( !_keepFiles ) {
1334 removeFile( aFacesFileName );
1335 removeFile( aPointsFileName );
1337 return error(COMPERR_BAD_INPUT_MESH);
1339 removeFile( aResultFileName ); // needed for boundary recovery module usage
1341 // -----------------
1343 // -----------------
1345 TCollection_AsciiString cmd =
1346 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1347 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1348 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1350 system( cmd.ToCString() ); // run
1356 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1357 if ( fileOpen < 0 ) {
1358 std::cout << std::endl;
1359 std::cout << "Error when opening the " << aResultFileName.ToCString() << " file" << std::endl;
1360 std::cout << "Log: " << aLogFileName << std::endl;
1361 std::cout << std::endl;
1365 Ok = readResultFile( fileOpen,
1367 aResultFileName.ToCString(),
1369 meshDS, theShape ,aNodeByGhs3dId, nbEnforcedVertices );
1372 // ---------------------
1373 // remove working files
1374 // ---------------------
1379 removeFile( aLogFileName );
1381 else if ( OSD_File( aLogFileName ).Size() > 0 )
1383 // get problem description from the log file
1384 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1385 storeErrorDescription( aLogFileName, conv );
1388 // the log file is empty
1389 removeFile( aLogFileName );
1390 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1391 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1396 removeFile( aFacesFileName );
1397 removeFile( aPointsFileName );
1398 removeFile( aResultFileName );
1399 removeFile( aBadResFileName );
1400 removeFile( aBbResFileName );
1406 //================================================================================
1408 * \brief Provide human readable text by error code reported by ghs3d
1410 //================================================================================
1412 static string translateError(const int errNum)
1416 return "The surface mesh includes a face of type other than edge, "
1417 "triangle or quadrilateral. This face type is not supported.";
1419 return "Not enough memory for the face table.";
1421 return "Not enough memory.";
1423 return "Not enough memory.";
1425 return "Face is ignored.";
1427 return "End of file. Some data are missing in the file.";
1429 return "Read error on the file. There are wrong data in the file.";
1431 return "the metric file is inadequate (dimension other than 3).";
1433 return "the metric file is inadequate (values not per vertices).";
1435 return "the metric file contains more than one field.";
1437 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1438 "value of number of mesh vertices in the \".noboite\" file.";
1440 return "Too many sub-domains.";
1442 return "the number of vertices is negative or null.";
1444 return "the number of faces is negative or null.";
1446 return "A face has a null vertex.";
1448 return "incompatible data.";
1450 return "the number of vertices is negative or null.";
1452 return "the number of vertices is negative or null (in the \".mesh\" file).";
1454 return "the number of faces is negative or null.";
1456 return "A face appears more than once in the input surface mesh.";
1458 return "An edge appears more than once in the input surface mesh.";
1460 return "A face has a vertex negative or null.";
1462 return "NOT ENOUGH MEMORY.";
1464 return "Not enough available memory.";
1466 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1467 "in terms of quality or the input list of points is wrong.";
1469 return "Some vertices are too close to one another or coincident.";
1471 return "Some vertices are too close to one another or coincident.";
1473 return "A vertex cannot be inserted.";
1475 return "There are at least two points considered as coincident.";
1477 return "Some vertices are too close to one another or coincident.";
1479 return "The surface mesh regeneration step has failed.";
1481 return "Constrained edge cannot be enforced.";
1483 return "Constrained face cannot be enforced.";
1485 return "Missing faces.";
1487 return "No guess to start the definition of the connected component(s).";
1489 return "The surface mesh includes at least one hole. The domain is not well defined.";
1491 return "Impossible to define a component.";
1493 return "The surface edge intersects another surface edge.";
1495 return "The surface edge intersects the surface face.";
1497 return "One boundary point lies within a surface face.";
1499 return "One surface edge intersects a surface face.";
1501 return "One boundary point lies within a surface edge.";
1503 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1504 "to too many swaps.";
1506 return "Edge is unique (i.e., bounds a hole in the surface).";
1508 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1510 return "Too many components, too many sub-domain.";
1512 return "The surface mesh includes at least one hole. "
1513 "Therefore there is no domain properly defined.";
1515 return "Statistics.";
1517 return "Statistics.";
1519 return "Warning, it is dramatically tedious to enforce the boundary items.";
1521 return "Not enough memory at this time, nevertheless, the program continues. "
1522 "The expected mesh will be correct but not really as large as required.";
1524 return "see above error code, resulting quality may be poor.";
1526 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1528 return "Unknown face type.";
1531 return "End of file. Some data are missing in the file.";
1533 return "A too small volume element is detected.";
1535 return "There exists at least a null or negative volume element.";
1537 return "There exist null or negative volume elements.";
1539 return "A too small volume element is detected. A face is considered being degenerated.";
1541 return "Some element is suspected to be very bad shaped or wrong.";
1543 return "A too bad quality face is detected. This face is considered degenerated.";
1545 return "A too bad quality face is detected. This face is degenerated.";
1547 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1549 return "Abnormal error occured, contact hotline.";
1551 return "Not enough memory for the face table.";
1553 return "The algorithm cannot run further. "
1554 "The surface mesh is probably very bad in terms of quality.";
1556 return "Bad vertex number.";
1561 //================================================================================
1563 * \brief Retrieve from a string given number of integers
1565 //================================================================================
1567 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1570 ids.reserve( nbIds );
1573 while ( !isdigit( *ptr )) ++ptr;
1574 if ( ptr[-1] == '-' ) --ptr;
1575 ids.push_back( strtol( ptr, &ptr, 10 ));
1581 //================================================================================
1583 * \brief Retrieve problem description form a log file
1584 * \retval bool - always false
1586 //================================================================================
1588 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1589 const _Ghs2smdsConvertor & toSmdsConvertor )
1593 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1595 int file = ::open (logFile.ToCString(), O_RDONLY);
1598 return error( SMESH_Comment("See ") << logFile << " for problem description");
1601 // struct stat status;
1602 // fstat(file, &status);
1603 // size_t length = status.st_size;
1604 off_t length = lseek( file, 0, SEEK_END);
1605 lseek( file, 0, SEEK_SET);
1608 vector< char > buf( length );
1609 int nBytesRead = ::read (file, & buf[0], length);
1611 char* ptr = & buf[0];
1612 char* bufEnd = ptr + nBytesRead;
1614 SMESH_Comment errDescription;
1616 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1618 // look for errors "ERR #"
1620 set<string> foundErrorStr; // to avoid reporting same error several times
1621 set<int> elemErrorNums; // not to report different types of errors with bad elements
1622 while ( ++ptr < bufEnd )
1624 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1627 list<const SMDS_MeshElement*> badElems;
1628 vector<int> nodeIds;
1632 int errNum = strtol(ptr, &ptr, 10);
1633 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1635 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1636 ptr = getIds(ptr, NODE, nodeIds);
1637 ptr = getIds(ptr, TRIA, nodeIds);
1638 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1640 case 1000: // ERR 1000 : 1 3 2
1641 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1642 ptr = getIds(ptr, TRIA, nodeIds);
1643 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1646 // Edge (e1, e2) appears more than once in the input surface mesh
1647 ptr = getIds(ptr, EDGE, nodeIds);
1648 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1651 // Face (f 1, f 2, f 3) has a vertex negative or null
1652 ptr = getIds(ptr, TRIA, nodeIds);
1653 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1656 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1657 ptr = getIds(ptr, NODE, nodeIds);
1658 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1659 ptr = getIds(ptr, NODE, nodeIds);
1660 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1663 // Vertex v1 cannot be inserted (warning).
1664 ptr = getIds(ptr, NODE, nodeIds);
1665 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1668 // There are at least two points whose distance is dist, i.e., considered as coincident
1669 case 2103: // ERR 2103 : 16 WITH 3
1670 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1671 ptr = getIds(ptr, NODE, nodeIds);
1672 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1673 ptr = getIds(ptr, NODE, nodeIds);
1674 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1677 // Constrained edge (e1, e2) cannot be enforced (warning).
1678 ptr = getIds(ptr, EDGE, nodeIds);
1679 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1682 // Constrained face (f 1, f 2, f 3) cannot be enforced
1683 ptr = getIds(ptr, TRIA, nodeIds);
1684 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1686 case 3103: // ERR 3103 : 1 2 WITH 7 3
1687 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1688 ptr = getIds(ptr, EDGE, nodeIds);
1689 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1690 ptr = getIds(ptr, EDGE, nodeIds);
1691 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1693 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1694 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1695 ptr = getIds(ptr, EDGE, nodeIds);
1696 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1697 ptr = getIds(ptr, TRIA, nodeIds);
1698 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1700 case 3105: // ERR 3105 : 8 IN 2 3 5
1701 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1702 ptr = getIds(ptr, NODE, nodeIds);
1703 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1704 ptr = getIds(ptr, TRIA, nodeIds);
1705 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1708 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1709 ptr = getIds(ptr, EDGE, nodeIds);
1710 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1711 ptr = getIds(ptr, TRIA, nodeIds);
1712 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1714 case 3107: // ERR 3107 : 2 IN 4 1
1715 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1716 ptr = getIds(ptr, NODE, nodeIds);
1717 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1718 ptr = getIds(ptr, EDGE, nodeIds);
1719 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1721 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1722 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1723 ptr = getIds(ptr, EDGE, nodeIds);
1724 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1726 case 9000: // ERR 9000
1727 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1728 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1729 // A too small volume element is detected. Are reported the index of the element,
1730 // its four vertex indices, its volume and the tolerance threshold value
1731 ptr = getIds(ptr, ID, nodeIds);
1732 ptr = getIds(ptr, VOL, nodeIds);
1733 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1734 // even if all nodes found, volume it most probably invisible,
1735 // add its faces to demenstrate it anyhow
1737 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1738 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1739 faceNodes[2] = nodeIds[3]; // 013
1740 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1741 faceNodes[1] = nodeIds[2]; // 023
1742 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1743 faceNodes[0] = nodeIds[1]; // 123
1744 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1747 case 9001: // ERR 9001
1748 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1749 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1750 // %% NUMBER OF NULL VOLUME TETS : 0
1751 // There exists at least a null or negative volume element
1754 // There exist n null or negative volume elements
1757 // A too small volume element is detected
1760 // A too bad quality face is detected. This face is considered degenerated,
1761 // its index, its three vertex indices together with its quality value are reported
1762 break; // same as next
1763 case 9112: // ERR 9112
1764 // FACE 2 WITH VERTICES : 4 2 5
1765 // SMALL INRADIUS : 0.
1766 // A too bad quality face is detected. This face is degenerated,
1767 // its index, its three vertex indices together with its inradius are reported
1768 ptr = getIds(ptr, ID, nodeIds);
1769 ptr = getIds(ptr, TRIA, nodeIds);
1770 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1771 // add triangle edges as it most probably has zero area and hence invisible
1773 vector<int> edgeNodes(2);
1774 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1775 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1776 edgeNodes[1] = nodeIds[2]; // 0-2
1777 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1778 edgeNodes[0] = nodeIds[1]; // 1-2
1779 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1784 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1786 continue; // not to report same error several times
1788 // const SMDS_MeshElement* nullElem = 0;
1789 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1791 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1792 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1793 // if ( oneMoreErrorType )
1794 // continue; // not to report different types of errors with bad elements
1797 // store bad elements
1798 //if ( allElemsOk ) {
1799 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1800 for ( ; elem != badElems.end(); ++elem )
1801 addBadInputElement( *elem );
1805 string text = translateError( errNum );
1806 if ( errDescription.find( text ) == text.npos ) {
1807 if ( !errDescription.empty() )
1808 errDescription << "\n";
1809 errDescription << text;
1814 if ( errDescription.empty() ) { // no errors found
1815 char msg[] = "connection to server failed";
1816 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1817 errDescription << "Licence problems.";
1820 if ( errDescription.empty() )
1821 errDescription << "See " << logFile << " for problem description";
1823 errDescription << "\nSee " << logFile << " for more information";
1825 return error( errDescription );
1828 //================================================================================
1830 * \brief Creates _Ghs2smdsConvertor
1832 //================================================================================
1834 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1835 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1839 //================================================================================
1841 * \brief Creates _Ghs2smdsConvertor
1843 //================================================================================
1845 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1846 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1850 //================================================================================
1852 * \brief Return SMDS element by ids of GHS3D nodes
1854 //================================================================================
1856 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1858 size_t nbNodes = ghsNodes.size();
1859 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1860 for ( size_t i = 0; i < nbNodes; ++i ) {
1861 int ghsNode = ghsNodes[ i ];
1862 if ( _ghs2NodeMap ) {
1863 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1864 if ( in == _ghs2NodeMap->end() )
1866 nodes[ i ] = in->second;
1869 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1871 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1877 if ( nbNodes == 2 ) {
1878 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1880 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1883 if ( nbNodes == 3 ) {
1884 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1886 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1890 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1896 //=============================================================================
1900 //=============================================================================
1901 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
1902 const TopoDS_Shape& aShape,
1903 MapShapeNbElems& aResMap)
1905 int nbtri = 0, nbqua = 0;
1906 double fullArea = 0.0;
1907 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1908 TopoDS_Face F = TopoDS::Face( exp.Current() );
1909 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1910 MapShapeNbElemsItr anIt = aResMap.find(sm);
1911 if( anIt==aResMap.end() ) {
1912 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1913 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1914 "Submesh can not be evaluated",this));
1917 std::vector<int> aVec = (*anIt).second;
1918 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
1919 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1921 BRepGProp::SurfaceProperties(F,G);
1922 double anArea = G.Mass();
1926 // collect info from edges
1927 int nb0d_e = 0, nb1d_e = 0;
1928 bool IsQuadratic = false;
1929 bool IsFirst = true;
1930 TopTools_MapOfShape tmpMap;
1931 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1932 TopoDS_Edge E = TopoDS::Edge(exp.Current());
1933 if( tmpMap.Contains(E) )
1936 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
1937 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
1938 std::vector<int> aVec = (*anIt).second;
1939 nb0d_e += aVec[SMDSEntity_Node];
1940 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
1942 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1948 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
1951 BRepGProp::VolumeProperties(aShape,G);
1952 double aVolume = G.Mass();
1953 double tetrVol = 0.1179*ELen*ELen*ELen;
1954 double CoeffQuality = 0.9;
1955 int nbVols = int(aVolume/tetrVol/CoeffQuality);
1956 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
1957 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
1958 std::vector<int> aVec(SMDSEntity_Last);
1959 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1961 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
1962 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
1963 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
1966 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
1967 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
1968 aVec[SMDSEntity_Pyramid] = nbqua;
1970 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1971 aResMap.insert(std::make_pair(sm,aVec));