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"
30 #include "SMESH_Gen.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_Comment.hxx"
33 #include "SMESH_MesherHelper.hxx"
34 #include "SMESH_MeshEditor.hxx"
36 #include "SMDS_MeshElement.hxx"
37 #include "SMDS_MeshNode.hxx"
38 #include "SMDS_FaceOfNodes.hxx"
39 #include "SMDS_VolumeOfNodes.hxx"
41 #include <BRepAdaptor_Surface.hxx>
42 #include <BRepBndLib.hxx>
43 #include <BRepClass3d_SolidClassifier.hxx>
44 #include <BRepTools.hxx>
45 #include <BRep_Tool.hxx>
46 #include <Bnd_Box.hxx>
47 #include <GeomAPI_ProjectPointOnSurf.hxx>
48 #include <OSD_File.hxx>
49 #include <Precision.hxx>
50 #include <Quantity_Parameter.hxx>
51 #include <Standard_ProgramError.hxx>
52 #include <Standard_ErrorHandler.hxx>
53 #include <Standard_Failure.hxx>
55 #include <TopExp_Explorer.hxx>
56 #include <TopTools_IndexedMapOfShape.hxx>
57 #include <TopTools_ListIteratorOfListOfShape.hxx>
59 //#include <BRepClass_FaceClassifier.hxx>
60 #include <TopTools_MapOfShape.hxx>
61 #include <BRepGProp.hxx>
62 #include <GProp_GProps.hxx>
64 #include "utilities.h"
69 #include <sys/sysinfo.h>
74 //#include <Standard_Stream.hxx>
77 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
98 static void removeFile( const TCollection_AsciiString& fileName )
101 OSD_File( fileName ).Remove();
103 catch ( Standard_ProgramError ) {
104 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
108 //=============================================================================
112 //=============================================================================
114 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
115 : SMESH_3D_Algo(hypId, studyId, gen)
117 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
119 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
120 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
123 _compatibleHypothesis.push_back("GHS3D_Parameters");
124 _requireShape = false; // can work without shape
127 //=============================================================================
131 //=============================================================================
133 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
135 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
138 //=============================================================================
142 //=============================================================================
144 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
145 const TopoDS_Shape& aShape,
146 Hypothesis_Status& aStatus )
148 aStatus = SMESH_Hypothesis::HYP_OK;
150 // there is only one compatible Hypothesis so far
153 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
155 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
157 _keepFiles = _hyp->GetKeepFiles();
162 //=======================================================================
163 //function : findShape
165 //=======================================================================
167 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
169 const TopoDS_Shape shape[],
172 TopAbs_State * state = 0)
175 int j, iShape, nbNode = 4;
177 for ( j=0; j<nbNode; j++ ) {
178 gp_XYZ p ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
179 if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) {
186 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
187 if (state) *state = SC.State();
188 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
189 for (iShape = 0; iShape < nShape; iShape++) {
190 aShape = shape[iShape];
191 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
192 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
193 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
194 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
195 if (state) *state = SC.State();
196 if (SC.State() == TopAbs_IN)
204 //=======================================================================
205 //function : readMapIntLine
207 //=======================================================================
209 static char* readMapIntLine(char* ptr, int tab[]) {
213 for ( int i=0; i<17; i++ ) {
214 intVal = strtol(ptr, &ptr, 10);
221 //=======================================================================
222 //function : countShape
224 //=======================================================================
226 template < class Mesh, class Shape >
227 static int countShape( Mesh* mesh, Shape shape ) {
228 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
230 for ( ; expShape.More(); expShape.Next() ) {
236 //=======================================================================
237 //function : writeFaces
239 //=======================================================================
241 static bool writeFaces (ofstream & theFile,
242 SMESHDS_Mesh * theMesh,
243 const map <int,int> & theSmdsToGhs3dIdMap)
247 // NB_ELEMS DUMMY_INT
248 // Loop from 1 to NB_ELEMS
249 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
251 int nbShape = countShape( theMesh, TopAbs_FACE );
253 int *tabID; tabID = new int[nbShape];
254 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
256 SMESHDS_SubMesh* theSubMesh;
257 const SMDS_MeshElement* aFace;
258 const char* space = " ";
259 const int dummyint = 0;
260 map<int,int>::const_iterator itOnMap;
261 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
262 int shapeID, nbNodes, aSmdsID;
265 cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
268 theFile << space << theMesh->NbFaces() << space << dummyint << endl;
270 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
271 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
273 aShape = expface.Current();
274 shapeID = theMesh->ShapeToIndex( aShape );
276 for ( int j=0; j<=i; j++) {
277 if ( shapeID == tabID[j] ) {
284 tabShape[i] = aShape;
287 for ( int i =0; i < nbShape; i++ ) {
288 if ( tabID[i] != 0 ) {
289 aShape = tabShape[i];
291 theSubMesh = theMesh->MeshElements( aShape );
292 if ( !theSubMesh ) continue;
293 itOnSubMesh = theSubMesh->GetElements();
294 while ( itOnSubMesh->more() ) {
295 aFace = itOnSubMesh->next();
296 nbNodes = aFace->NbNodes();
298 theFile << space << nbNodes;
300 itOnSubFace = aFace->nodesIterator();
301 while ( itOnSubFace->more() ) {
303 aSmdsID = itOnSubFace->next()->GetID();
304 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
305 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
307 theFile << space << (*itOnMap).second;
310 // (NB_NODES + 1) times: DUMMY_INT
311 for ( int j=0; j<=nbNodes; j++)
312 theFile << space << dummyint;
325 //=======================================================================
326 //function : writeFaces
327 //purpose : Write Faces in case if generate 3D mesh w/o geometry
328 //=======================================================================
330 static bool writeFaces (ofstream & theFile,
331 SMESHDS_Mesh * theMesh,
332 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
336 // NB_ELEMS DUMMY_INT
337 // Loop from 1 to NB_ELEMS
338 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
342 list< const SMDS_MeshElement* > faces;
343 list< const SMDS_MeshElement* >::iterator f;
344 map< const SMDS_MeshNode*,int >::iterator it;
345 SMDS_ElemIteratorPtr nodeIt;
346 const SMDS_MeshElement* elem;
349 const char* space = " ";
350 const int dummyint = 0;
352 //get all faces from mesh
353 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
354 while ( eIt->more() ) {
355 const SMDS_MeshElement* elem = eIt->next();
358 faces.push_back( elem );
365 cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
367 // NB_ELEMS DUMMY_INT
368 theFile << space << nbFaces << space << dummyint << endl;
370 // Loop from 1 to NB_ELEMS
372 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
374 for ( ; f != faces.end(); ++f )
378 nbNodes = elem->NbNodes();
379 theFile << space << nbNodes;
381 // NODE_NB_1 NODE_NB_2 ...
382 nodeIt = elem->nodesIterator();
383 while ( nodeIt->more() )
386 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
387 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
388 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
389 theFile << space << it->second;
392 // (NB_NODES + 1) times: DUMMY_INT
393 for ( int i=0; i<=nbNodes; i++)
394 theFile << space << dummyint;
399 // put nodes to theNodeByGhs3dId vector
400 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
401 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
402 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
404 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
410 //=======================================================================
411 //function : writePoints
413 //=======================================================================
415 static bool writePoints (ofstream & theFile,
416 SMESHDS_Mesh * theMesh,
417 map <int,int> & theSmdsToGhs3dIdMap,
418 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
423 // Loop from 1 to NB_NODES
426 int nbNodes = theMesh->NbNodes();
430 const char* space = " ";
431 const int dummyint = 0;
434 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
435 const SMDS_MeshNode* node;
438 theFile << space << nbNodes << endl;
440 cout << "The initial 2D mesh contains :" << endl;
441 cout << " " << nbNodes << " vertices" << endl;
443 // Loop from 1 to NB_NODES
448 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
449 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
454 << space << node->X()
455 << space << node->Y()
456 << space << node->Z()
457 << space << dummyint;
465 //=======================================================================
466 //function : writePoints
468 //=======================================================================
470 static bool writePoints (ofstream & theFile,
471 SMESHDS_Mesh * theMesh,
472 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
477 // Loop from 1 to NB_NODES
480 //int nbNodes = theMesh->NbNodes();
481 int nbNodes = theNodeByGhs3dId.size();
485 const char* space = " ";
486 const int dummyint = 0;
488 const SMDS_MeshNode* node;
491 theFile << space << nbNodes << endl;
492 cout << nbNodes << " nodes" << endl;
494 // Loop from 1 to NB_NODES
496 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
497 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
498 for ( ; nodeIt != after; ++nodeIt )
504 << space << node->X()
505 << space << node->Y()
506 << space << node->Z()
507 << space << dummyint;
515 //=======================================================================
516 //function : findShapeID
517 //purpose : find the solid corresponding to GHS3D sub-domain following
518 // the technique proposed in GHS3D manual in chapter
519 // "B.4 Subdomain (sub-region) assignment"
520 //=======================================================================
522 static int findShapeID(SMESH_Mesh& mesh,
523 const SMDS_MeshNode* node1,
524 const SMDS_MeshNode* node2,
525 const SMDS_MeshNode* node3,
526 const bool toMeshHoles)
528 const int invalidID = 0;
529 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
531 // face the nodes belong to
532 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
536 // geom face the face assigned to
537 SMESH_MeshEditor editor(&mesh);
538 int geomFaceID = editor.FindShape( face );
541 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
542 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
544 TopoDS_Face geomFace = TopoDS::Face( shape );
546 // solids bounded by geom face
547 TopTools_IndexedMapOfShape solids, shells;
548 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
549 for ( ; ansIt.More(); ansIt.Next() ) {
550 switch ( ansIt.Value().ShapeType() ) {
552 solids.Add( ansIt.Value() ); break;
554 shells.Add( ansIt.Value() ); break;
558 // analyse found solids
559 if ( solids.Extent() == 0 || shells.Extent() == 0)
562 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
563 if ( solids.Extent() == 1 )
566 return meshDS->ShapeToIndex( solid1 );
568 // - Are we at a hole boundary face?
569 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
570 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
572 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
573 // check if any edge of shells(1) belongs to another shell
574 for ( ; eExp.More() && !touch; eExp.Next() ) {
575 ansIt = mesh.GetAncestors( eExp.Current() );
576 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
577 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
578 touch = ( !ansIt.Value().IsSame( shells(1) ));
582 return meshDS->ShapeToIndex( solid1 );
585 // find orientation of geom face within the first solid
586 TopExp_Explorer fExp( solid1, TopAbs_FACE );
587 for ( ; fExp.More(); fExp.Next() )
588 if ( geomFace.IsSame( fExp.Current() )) {
589 geomFace = TopoDS::Face( fExp.Current() );
593 return invalidID; // face not found
595 // normale to triangle
596 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
597 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
598 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
599 gp_Vec vec12( node1Pnt, node2Pnt );
600 gp_Vec vec13( node1Pnt, node3Pnt );
601 gp_Vec meshNormal = vec12 ^ vec13;
602 if ( meshNormal.SquareMagnitude() < DBL_MIN )
605 // find UV of node1 on geomFace
606 SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
607 const SMDS_MeshNode* nNotOnSeamEdge = 0;
608 if ( helper.IsSeamShape( node1->GetPosition()->GetShapeId() ))
609 if ( helper.IsSeamShape( node2->GetPosition()->GetShapeId() ))
610 nNotOnSeamEdge = node3;
612 nNotOnSeamEdge = node2;
614 gp_XY uv = helper.GetNodeUV( geomFace, node1, nNotOnSeamEdge, &checkUV );
615 // check that uv is correct
617 TopoDS_Shape nodeShape = helper.GetSubShapeByNode( node1, meshDS );
618 if ( !nodeShape.IsNull() )
619 switch ( nodeShape.ShapeType() )
621 case TopAbs_FACE: tol = BRep_Tool::Tolerance( TopoDS::Face( nodeShape )); break;
622 case TopAbs_EDGE: tol = BRep_Tool::Tolerance( TopoDS::Edge( nodeShape )); break;
623 case TopAbs_VERTEX: tol = BRep_Tool::Tolerance( TopoDS::Vertex( nodeShape )); break;
626 BRepAdaptor_Surface surface( geomFace );
627 if ( checkUV && node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol )
630 // normale to geomFace at UV
632 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
633 gp_Vec geomNormal = du ^ dv;
634 if ( geomNormal.SquareMagnitude() < DBL_MIN )
635 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
636 if ( geomFace.Orientation() == TopAbs_REVERSED )
637 geomNormal.Reverse();
640 bool isReverse = ( meshNormal * geomNormal ) < 0;
642 return meshDS->ShapeToIndex( solid1 );
644 if ( solids.Extent() == 1 )
645 return HOLE_ID; // we are inside a hole
647 return meshDS->ShapeToIndex( solids(2) );
650 //=======================================================================
651 //function : readResultFile
653 //=======================================================================
655 static bool readResultFile(const int fileOpen,
657 const char* fileName,
660 TopoDS_Shape tabShape[],
663 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
673 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
676 int nbElems, nbNodes, nbInputNodes;
677 int nodeId/*, triangleId*/;
679 int ID, shapeID, ghs3dShapeID;
682 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
684 int *tab, *tabID, *nodeID, *nodeAssigne;
686 const SMDS_MeshNode **node;
689 //tabID = new int[nbShape];
691 coord = new double[3];
692 node = new const SMDS_MeshNode*[4];
695 SMDS_MeshNode * aNewNode;
696 map <int,const SMDS_MeshNode*>::iterator itOnNode;
697 SMDS_MeshElement* aTet;
702 // Read the file state
703 fileStat = fstat(fileOpen, &status);
704 length = status.st_size;
706 // Mapping the result file into memory
708 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
709 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
710 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
711 0, (DWORD)length, NULL);
712 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
714 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
718 ptr = readMapIntLine(ptr, tab);
723 nbInputNodes = tab[2];
725 nodeAssigne = new int[ nbNodes+1 ];
728 aSolid = tabShape[0];
730 // Reading the nodeId
731 for (int i=0; i < 4*nbElems; i++)
732 nodeId = strtol(ptr, &ptr, 10);
734 // Reading the nodeCoor and update the nodeMap
735 for (int iNode=1; iNode <= nbNodes; iNode++) {
736 for (int iCoor=0; iCoor < 3; iCoor++)
737 coord[ iCoor ] = strtod(ptr, &ptr);
738 nodeAssigne[ iNode ] = 1;
739 if ( iNode > nbInputNodes ) {
740 nodeAssigne[ iNode ] = 0;
741 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
742 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
746 // Reading the number of triangles which corresponds to the number of sub-domains
747 nbTriangle = strtol(ptr, &ptr, 10);
749 tabID = new int[nbTriangle];
750 for (int i=0; i < nbTriangle; i++) {
752 // find the solid corresponding to GHS3D sub-domain following
753 // the technique proposed in GHS3D manual in chapter
754 // "B.4 Subdomain (sub-region) assignment"
755 int nodeId1 = strtol(ptr, &ptr, 10);
756 int nodeId2 = strtol(ptr, &ptr, 10);
757 int nodeId3 = strtol(ptr, &ptr, 10);
758 if ( nbTriangle > 1 ) {
759 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
760 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
761 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
764 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
765 // -- 0020330: Pb with ghs3d as a submesh
766 // check that found shape is to be meshed
767 if ( tabID[i] > 0 ) {
768 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] );
769 bool isToBeMeshed = false;
770 for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS )
771 isToBeMeshed = foundShape.IsSame( tabShape[ iS ]);
775 // END -- 0020330: Pb with ghs3d as a submesh
777 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
779 } catch ( Standard_Failure ) {
786 if ( nbTriangle <= nbShape ) // no holes
787 toMeshHoles = true; // not avoid creating tetras in holes
789 // Associating the tetrahedrons to the shapes
790 shapeID = compoundID;
791 for (int iElem = 0; iElem < nbElems; iElem++) {
792 for (int iNode = 0; iNode < 4; iNode++) {
793 ID = strtol(tetraPtr, &tetraPtr, 10);
794 itOnNode = theGhs3dIdToNodeMap.find(ID);
795 node[ iNode ] = itOnNode->second;
796 nodeID[ iNode ] = ID;
798 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
799 // tetras within holes depending on hypo option,
800 // so we first check if aTet is inside a hole and then create it
801 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
802 if ( nbTriangle > 1 ) {
803 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
804 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
805 if ( tabID[ ghs3dShapeID ] == 0 ) {
807 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
808 if ( toMeshHoles || state == TopAbs_IN )
809 shapeID = theMeshDS->ShapeToIndex( aSolid );
810 tabID[ ghs3dShapeID ] = shapeID;
813 shapeID = tabID[ ghs3dShapeID ];
815 else if ( nbShape > 1 ) {
816 // Case where nbTriangle == 1 while nbShape == 2 encountered
817 // with compound of 2 boxes and "To mesh holes"==False,
818 // so there are no subdomains specified for each tetrahedron.
819 // Try to guess a solid by a node already bound to shape
821 for ( int i=0; i<4 && shapeID==0; i++ ) {
822 if ( nodeAssigne[ nodeID[i] ] == 1 &&
823 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
824 node[i]->GetPosition()->GetShapeId() > 1 )
826 shapeID = node[i]->GetPosition()->GetShapeId();
830 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
831 shapeID = theMeshDS->ShapeToIndex( aSolid );
834 // set new nodes and tetrahedron onto the shape
835 for ( int i=0; i<4; i++ ) {
836 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
837 if ( shapeID != HOLE_ID )
838 theMeshDS->SetNodeInVolume( node[i], shapeID );
839 nodeAssigne[ nodeID[i] ] = shapeID;
842 if ( toMeshHoles || shapeID != HOLE_ID ) {
843 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
844 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
847 shapeIDs.insert( shapeID );
850 // Remove nodes of tetras inside holes if !toMeshHoles
851 if ( !toMeshHoles ) {
852 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
853 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
854 ID = itOnNode->first;
855 if ( nodeAssigne[ ID ] == HOLE_ID )
856 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
861 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
863 UnmapViewOfFile(mapPtr);
864 CloseHandle(hMapObject);
867 munmap(mapPtr, length);
876 delete [] nodeAssigne;
879 if ( shapeIDs.size() != nbShape ) {
880 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
881 for (int i=0; i<nbShape; i++) {
882 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
883 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
884 cout << " Solid #" << shapeID << " not found" << endl;
892 //=======================================================================
893 //function : readResultFile
895 //=======================================================================
897 static bool readResultFile(const int fileOpen,
899 const char* fileName,
901 SMESHDS_Mesh* theMeshDS,
903 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
913 int nbElems, nbNodes, nbInputNodes;
914 int nodeId, triangleId;
920 const SMDS_MeshNode **node;
923 coord = new double[3];
924 node = new const SMDS_MeshNode*[4];
926 SMDS_MeshNode * aNewNode;
927 map <int,const SMDS_MeshNode*>::iterator IdNode;
928 SMDS_MeshElement* aTet;
930 // Read the file state
931 fileStat = fstat(fileOpen, &status);
932 length = status.st_size;
934 // Mapping the result file into memory
936 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
937 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
938 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
939 0, (DWORD)length, NULL);
940 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
942 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
946 ptr = readMapIntLine(ptr, tab);
951 nbInputNodes = tab[2];
953 theNodeByGhs3dId.resize( nbNodes );
955 // Reading the nodeId
956 for (int i=0; i < 4*nbElems; i++)
957 nodeId = strtol(ptr, &ptr, 10);
959 // Reading the nodeCoor and update the nodeMap
960 shapeID = theMeshDS->ShapeToIndex( aSolid );
961 for (int iNode=0; iNode < nbNodes; iNode++) {
962 for (int iCoor=0; iCoor < 3; iCoor++)
963 coord[ iCoor ] = strtod(ptr, &ptr);
964 if ((iNode+1) > nbInputNodes) {
965 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
966 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
967 theNodeByGhs3dId[ iNode ] = aNewNode;
971 // Reading the triangles
972 nbTriangle = strtol(ptr, &ptr, 10);
974 for (int i=0; i < 3*nbTriangle; i++)
975 triangleId = strtol(ptr, &ptr, 10);
979 // Associating the tetrahedrons to the shapes
980 for (int iElem = 0; iElem < nbElems; iElem++) {
981 for (int iNode = 0; iNode < 4; iNode++) {
982 ID = strtol(tetraPtr, &tetraPtr, 10);
983 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
985 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
986 shapeID = theMeshDS->ShapeToIndex( aSolid );
987 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
990 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
992 UnmapViewOfFile(mapPtr);
993 CloseHandle(hMapObject);
996 munmap(mapPtr, length);
1007 //=============================================================================
1009 *Here we are going to use the GHS3D mesher
1011 //=============================================================================
1013 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1014 const TopoDS_Shape& theShape)
1017 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1019 // we count the number of shapes
1020 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
1022 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1023 for ( ; expBox.More(); expBox.Next() )
1026 // create bounding box for every shape inside the compound
1029 TopoDS_Shape* tabShape;
1031 tabShape = new TopoDS_Shape[_nbShape];
1032 tabBox = new double*[_nbShape];
1033 for (int i=0; i<_nbShape; i++)
1034 tabBox[i] = new double[6];
1035 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1037 // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh
1038 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
1039 tabShape[iShape] = expBox.Current();
1040 Bnd_Box BoundingBox;
1041 BRepBndLib::Add(expBox.Current(), BoundingBox);
1042 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1043 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1044 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1045 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1049 // a unique working file name
1050 // to avoid access to the same files by eg different users
1051 TCollection_AsciiString aGenericName
1052 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1054 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1055 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1056 aFacesFileName = aGenericName + ".faces"; // in faces
1057 aPointsFileName = aGenericName + ".points"; // in points
1058 aResultFileName = aGenericName + ".noboite";// out points and volumes
1059 aBadResFileName = aGenericName + ".boite"; // out bad result
1060 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1061 aLogFileName = aGenericName + ".log"; // log
1063 // -----------------
1065 // -----------------
1067 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1068 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1071 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1073 INFOS( "Can't write into " << aFacesFileName);
1074 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1076 map <int,int> aSmdsToGhs3dIdMap;
1077 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1079 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1080 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1083 aPointsFile.close();
1086 if ( !_keepFiles ) {
1087 removeFile( aFacesFileName );
1088 removeFile( aPointsFileName );
1090 return error(COMPERR_BAD_INPUT_MESH);
1092 removeFile( aResultFileName ); // needed for boundary recovery module usage
1094 // -----------------
1096 // -----------------
1098 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1099 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1100 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1103 cout << "Ghs3d execution..." << endl;
1104 cout << cmd << endl;
1106 system( cmd.ToCString() ); // run
1109 cout << "End of Ghs3d execution !" << endl;
1115 // Mapping the result file
1118 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1119 if ( fileOpen < 0 ) {
1121 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1122 cout << "Log: " << aLogFileName << endl;
1127 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1128 Ok = readResultFile( fileOpen,
1130 aResultFileName.ToCString(),
1132 theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1136 // ---------------------
1137 // remove working files
1138 // ---------------------
1143 removeFile( aLogFileName );
1145 else if ( OSD_File( aLogFileName ).Size() > 0 )
1147 // get problem description from the log file
1148 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1149 storeErrorDescription( aLogFileName, conv );
1153 // the log file is empty
1154 removeFile( aLogFileName );
1155 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1156 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1159 if ( !_keepFiles ) {
1160 removeFile( aFacesFileName );
1161 removeFile( aPointsFileName );
1162 removeFile( aResultFileName );
1163 removeFile( aBadResFileName );
1164 removeFile( aBbResFileName );
1166 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1169 cout << "treated !" << endl;
1172 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1179 //=============================================================================
1181 *Here we are going to use the GHS3D mesher w/o geometry
1183 //=============================================================================
1184 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1185 SMESH_MesherHelper* aHelper)
1187 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1189 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1190 TopoDS_Shape theShape = aHelper->GetSubShape();
1192 // a unique working file name
1193 // to avoid access to the same files by eg different users
1194 TCollection_AsciiString aGenericName
1195 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1197 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1198 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1199 aFacesFileName = aGenericName + ".faces"; // in faces
1200 aPointsFileName = aGenericName + ".points"; // in points
1201 aResultFileName = aGenericName + ".noboite";// out points and volumes
1202 aBadResFileName = aGenericName + ".boite"; // out bad result
1203 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1204 aLogFileName = aGenericName + ".log"; // log
1206 // -----------------
1208 // -----------------
1210 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1211 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1213 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1216 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1218 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1220 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1221 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1224 aPointsFile.close();
1227 if ( !_keepFiles ) {
1228 removeFile( aFacesFileName );
1229 removeFile( aPointsFileName );
1231 return error(COMPERR_BAD_INPUT_MESH);
1233 removeFile( aResultFileName ); // needed for boundary recovery module usage
1235 // -----------------
1237 // -----------------
1239 TCollection_AsciiString cmd =
1240 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1241 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1242 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1244 system( cmd.ToCString() ); // run
1250 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1251 if ( fileOpen < 0 ) {
1253 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1254 cout << "Log: " << aLogFileName << endl;
1259 Ok = readResultFile( fileOpen,
1261 aResultFileName.ToCString(),
1263 meshDS, theShape ,aNodeByGhs3dId );
1266 // ---------------------
1267 // remove working files
1268 // ---------------------
1273 removeFile( aLogFileName );
1275 else if ( OSD_File( aLogFileName ).Size() > 0 )
1277 // get problem description from the log file
1278 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1279 storeErrorDescription( aLogFileName, conv );
1282 // the log file is empty
1283 removeFile( aLogFileName );
1284 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1285 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1290 removeFile( aFacesFileName );
1291 removeFile( aPointsFileName );
1292 removeFile( aResultFileName );
1293 removeFile( aBadResFileName );
1294 removeFile( aBbResFileName );
1300 //================================================================================
1302 * \brief Provide human readable text by error code reported by ghs3d
1304 //================================================================================
1306 static string translateError(const int errNum)
1310 return "The surface mesh includes a face of type other than edge, "
1311 "triangle or quadrilateral. This face type is not supported.";
1313 return "Not enough memory for the face table.";
1315 return "Not enough memory.";
1317 return "Not enough memory.";
1319 return "Face is ignored.";
1321 return "End of file. Some data are missing in the file.";
1323 return "Read error on the file. There are wrong data in the file.";
1325 return "the metric file is inadequate (dimension other than 3).";
1327 return "the metric file is inadequate (values not per vertices).";
1329 return "the metric file contains more than one field.";
1331 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1332 "value of number of mesh vertices in the \".noboite\" file.";
1334 return "Too many sub-domains.";
1336 return "the number of vertices is negative or null.";
1338 return "the number of faces is negative or null.";
1340 return "A face has a null vertex.";
1342 return "incompatible data.";
1344 return "the number of vertices is negative or null.";
1346 return "the number of vertices is negative or null (in the \".mesh\" file).";
1348 return "the number of faces is negative or null.";
1350 return "A face appears more than once in the input surface mesh.";
1352 return "An edge appears more than once in the input surface mesh.";
1354 return "A face has a vertex negative or null.";
1356 return "NOT ENOUGH MEMORY.";
1358 return "Not enough available memory.";
1360 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1361 "in terms of quality or the input list of points is wrong.";
1363 return "Some vertices are too close to one another or coincident.";
1365 return "Some vertices are too close to one another or coincident.";
1367 return "A vertex cannot be inserted.";
1369 return "There are at least two points considered as coincident.";
1371 return "Some vertices are too close to one another or coincident.";
1373 return "The surface mesh regeneration step has failed.";
1375 return "Constrained edge cannot be enforced.";
1377 return "Constrained face cannot be enforced.";
1379 return "Missing faces.";
1381 return "No guess to start the definition of the connected component(s).";
1383 return "The surface mesh includes at least one hole. The domain is not well defined.";
1385 return "Impossible to define a component.";
1387 return "The surface edge intersects another surface edge.";
1389 return "The surface edge intersects the surface face.";
1391 return "One boundary point lies within a surface face.";
1393 return "One surface edge intersects a surface face.";
1395 return "One boundary point lies within a surface edge.";
1397 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1398 "to too many swaps.";
1400 return "Edge is unique (i.e., bounds a hole in the surface).";
1402 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1404 return "Too many components, too many sub-domain.";
1406 return "The surface mesh includes at least one hole. "
1407 "Therefore there is no domain properly defined.";
1409 return "Statistics.";
1411 return "Statistics.";
1413 return "Warning, it is dramatically tedious to enforce the boundary items.";
1415 return "Not enough memory at this time, nevertheless, the program continues. "
1416 "The expected mesh will be correct but not really as large as required.";
1418 return "see above error code, resulting quality may be poor.";
1420 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1422 return "Unknown face type.";
1425 return "End of file. Some data are missing in the file.";
1427 return "A too small volume element is detected.";
1429 return "There exists at least a null or negative volume element.";
1431 return "There exist null or negative volume elements.";
1433 return "A too small volume element is detected. A face is considered being degenerated.";
1435 return "Some element is suspected to be very bad shaped or wrong.";
1437 return "A too bad quality face is detected. This face is considered degenerated.";
1439 return "A too bad quality face is detected. This face is degenerated.";
1441 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1443 return "Abnormal error occured, contact hotline.";
1445 return "Not enough memory for the face table.";
1447 return "The algorithm cannot run further. "
1448 "The surface mesh is probably very bad in terms of quality.";
1450 return "Bad vertex number.";
1455 //================================================================================
1457 * \brief Retrieve from a string given number of integers
1459 //================================================================================
1461 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1464 ids.reserve( nbIds );
1467 while ( !isdigit( *ptr )) ++ptr;
1468 if ( ptr[-1] == '-' ) --ptr;
1469 ids.push_back( strtol( ptr, &ptr, 10 ));
1475 //================================================================================
1477 * \brief Retrieve problem description form a log file
1478 * \retval bool - always false
1480 //================================================================================
1482 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1483 const _Ghs2smdsConvertor & toSmdsConvertor )
1487 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1489 int file = ::open (logFile.ToCString(), O_RDONLY);
1492 return error( SMESH_Comment("See ") << logFile << " for problem description");
1495 // struct stat status;
1496 // fstat(file, &status);
1497 // size_t length = status.st_size;
1498 off_t length = lseek( file, 0, SEEK_END);
1499 lseek( file, 0, SEEK_SET);
1502 vector< char > buf( length );
1503 int nBytesRead = ::read (file, & buf[0], length);
1505 char* ptr = & buf[0];
1506 char* bufEnd = ptr + nBytesRead;
1508 SMESH_Comment errDescription;
1510 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1512 // look for errors "ERR #"
1514 set<string> foundErrorStr; // to avoid reporting same error several times
1515 set<int> elemErrorNums; // not to report different types of errors with bad elements
1516 while ( ++ptr < bufEnd )
1518 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1521 list<const SMDS_MeshElement*> badElems;
1522 vector<int> nodeIds;
1526 int errNum = strtol(ptr, &ptr, 10);
1527 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1529 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1530 ptr = getIds(ptr, NODE, nodeIds);
1531 ptr = getIds(ptr, TRIA, nodeIds);
1532 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1534 case 1000: // ERR 1000 : 1 3 2
1535 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1536 ptr = getIds(ptr, TRIA, nodeIds);
1537 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1540 // Edge (e1, e2) appears more than once in the input surface mesh
1541 ptr = getIds(ptr, EDGE, nodeIds);
1542 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1545 // Face (f 1, f 2, f 3) has a vertex negative or null
1546 ptr = getIds(ptr, TRIA, nodeIds);
1547 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1550 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1551 ptr = getIds(ptr, NODE, nodeIds);
1552 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1553 ptr = getIds(ptr, NODE, nodeIds);
1554 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1557 // Vertex v1 cannot be inserted (warning).
1558 ptr = getIds(ptr, NODE, nodeIds);
1559 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1562 // There are at least two points whose distance is dist, i.e., considered as coincident
1563 case 2103: // ERR 2103 : 16 WITH 3
1564 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1565 ptr = getIds(ptr, NODE, nodeIds);
1566 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1567 ptr = getIds(ptr, NODE, nodeIds);
1568 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1571 // Constrained edge (e1, e2) cannot be enforced (warning).
1572 ptr = getIds(ptr, EDGE, nodeIds);
1573 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1576 // Constrained face (f 1, f 2, f 3) cannot be enforced
1577 ptr = getIds(ptr, TRIA, nodeIds);
1578 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1580 case 3103: // ERR 3103 : 1 2 WITH 7 3
1581 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1582 ptr = getIds(ptr, EDGE, nodeIds);
1583 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1584 ptr = getIds(ptr, EDGE, nodeIds);
1585 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1587 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1588 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1589 ptr = getIds(ptr, EDGE, nodeIds);
1590 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1591 ptr = getIds(ptr, TRIA, nodeIds);
1592 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1594 case 3105: // ERR 3105 : 8 IN 2 3 5
1595 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1596 ptr = getIds(ptr, NODE, nodeIds);
1597 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1598 ptr = getIds(ptr, TRIA, nodeIds);
1599 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1602 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1603 ptr = getIds(ptr, EDGE, nodeIds);
1604 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1605 ptr = getIds(ptr, TRIA, nodeIds);
1606 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1608 case 3107: // ERR 3107 : 2 IN 4 1
1609 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1610 ptr = getIds(ptr, NODE, nodeIds);
1611 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1612 ptr = getIds(ptr, EDGE, nodeIds);
1613 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1615 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1616 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1617 ptr = getIds(ptr, EDGE, nodeIds);
1618 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1620 case 9000: // ERR 9000
1621 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1622 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1623 // A too small volume element is detected. Are reported the index of the element,
1624 // its four vertex indices, its volume and the tolerance threshold value
1625 ptr = getIds(ptr, ID, nodeIds);
1626 ptr = getIds(ptr, VOL, nodeIds);
1627 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1628 // even if all nodes found, volume it most probably invisible,
1629 // add its faces to demenstrate it anyhow
1631 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1632 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1633 faceNodes[2] = nodeIds[3]; // 013
1634 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1635 faceNodes[1] = nodeIds[2]; // 023
1636 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1637 faceNodes[0] = nodeIds[1]; // 123
1638 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1641 case 9001: // ERR 9001
1642 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1643 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1644 // %% NUMBER OF NULL VOLUME TETS : 0
1645 // There exists at least a null or negative volume element
1648 // There exist n null or negative volume elements
1651 // A too small volume element is detected
1654 // A too bad quality face is detected. This face is considered degenerated,
1655 // its index, its three vertex indices together with its quality value are reported
1656 break; // same as next
1657 case 9112: // ERR 9112
1658 // FACE 2 WITH VERTICES : 4 2 5
1659 // SMALL INRADIUS : 0.
1660 // A too bad quality face is detected. This face is degenerated,
1661 // its index, its three vertex indices together with its inradius are reported
1662 ptr = getIds(ptr, ID, nodeIds);
1663 ptr = getIds(ptr, TRIA, nodeIds);
1664 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1665 // add triangle edges as it most probably has zero area and hence invisible
1667 vector<int> edgeNodes(2);
1668 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1669 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1670 edgeNodes[1] = nodeIds[2]; // 0-2
1671 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1672 edgeNodes[0] = nodeIds[1]; // 1-2
1673 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1678 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1680 continue; // not to report same error several times
1682 // const SMDS_MeshElement* nullElem = 0;
1683 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1685 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1686 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1687 // if ( oneMoreErrorType )
1688 // continue; // not to report different types of errors with bad elements
1691 // store bad elements
1692 //if ( allElemsOk ) {
1693 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1694 for ( ; elem != badElems.end(); ++elem )
1695 addBadInputElement( *elem );
1699 string text = translateError( errNum );
1700 if ( errDescription.find( text ) == text.npos ) {
1701 if ( !errDescription.empty() )
1702 errDescription << "\n";
1703 errDescription << text;
1708 if ( errDescription.empty() ) { // no errors found
1709 char msg[] = "connection to server failed";
1710 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1711 errDescription << "Licence problems.";
1714 if ( errDescription.empty() )
1715 errDescription << "See " << logFile << " for problem description";
1717 errDescription << "\nSee " << logFile << " for more information";
1719 return error( errDescription );
1722 //================================================================================
1724 * \brief Creates _Ghs2smdsConvertor
1726 //================================================================================
1728 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1729 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1733 //================================================================================
1735 * \brief Creates _Ghs2smdsConvertor
1737 //================================================================================
1739 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1740 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1744 //================================================================================
1746 * \brief Return SMDS element by ids of GHS3D nodes
1748 //================================================================================
1750 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1752 size_t nbNodes = ghsNodes.size();
1753 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1754 for ( size_t i = 0; i < nbNodes; ++i ) {
1755 int ghsNode = ghsNodes[ i ];
1756 if ( _ghs2NodeMap ) {
1757 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1758 if ( in == _ghs2NodeMap->end() )
1760 nodes[ i ] = in->second;
1763 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1765 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1771 if ( nbNodes == 2 ) {
1772 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1774 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1777 if ( nbNodes == 3 ) {
1778 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1780 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1784 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1790 //=============================================================================
1794 //=============================================================================
1795 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
1796 const TopoDS_Shape& aShape,
1797 MapShapeNbElems& aResMap)
1799 int nbtri = 0, nbqua = 0;
1800 double fullArea = 0.0;
1801 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1802 TopoDS_Face F = TopoDS::Face( exp.Current() );
1803 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1804 MapShapeNbElemsItr anIt = aResMap.find(sm);
1805 if( anIt==aResMap.end() ) {
1806 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1807 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1808 "Submesh can not be evaluated",this));
1811 std::vector<int> aVec = (*anIt).second;
1812 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
1813 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1815 BRepGProp::SurfaceProperties(F,G);
1816 double anArea = G.Mass();
1820 // collect info from edges
1821 int nb0d_e = 0, nb1d_e = 0;
1822 bool IsQuadratic = false;
1823 bool IsFirst = true;
1824 TopTools_MapOfShape tmpMap;
1825 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1826 TopoDS_Edge E = TopoDS::Edge(exp.Current());
1827 if( tmpMap.Contains(E) )
1830 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
1831 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
1832 std::vector<int> aVec = (*anIt).second;
1833 nb0d_e += aVec[SMDSEntity_Node];
1834 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
1836 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1842 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
1845 BRepGProp::VolumeProperties(aShape,G);
1846 double aVolume = G.Mass();
1847 double tetrVol = 0.1179*ELen*ELen*ELen;
1848 double CoeffQuality = 0.9;
1849 int nbVols = int(aVolume/tetrVol/CoeffQuality);
1850 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
1851 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
1852 std::vector<int> aVec(SMDSEntity_Last);
1853 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1855 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
1856 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
1857 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
1860 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
1861 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
1862 aVec[SMDSEntity_Pyramid] = nbqua;
1864 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1865 aResMap.insert(std::make_pair(sm,aVec));