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
616 BRepAdaptor_Surface surface( geomFace );
617 double tol = BRep_Tool::Tolerance( geomFace );
618 if ( checkUV && node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
619 // GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface(), 2 * tol );
620 // if ( !projector.IsDone() || projector.NbPoints() < 1 || projector.LowerDistance() > 2 * tol)
622 // Quantity_Parameter U,V;
623 // projector.LowerDistanceParameters(U,V);
624 // if ( node1Pnt.Distance( surface.Value( U, V )) > 2 * tol )
626 // uv.SetCoord( U,V );
628 // normale to geomFace at UV
630 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
631 gp_Vec geomNormal = du ^ dv;
632 if ( geomNormal.SquareMagnitude() < DBL_MIN )
633 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
634 if ( geomFace.Orientation() == TopAbs_REVERSED )
635 geomNormal.Reverse();
638 bool isReverse = ( meshNormal * geomNormal ) < 0;
640 return meshDS->ShapeToIndex( solid1 );
642 if ( solids.Extent() == 1 )
643 return HOLE_ID; // we are inside a hole
645 return meshDS->ShapeToIndex( solids(2) );
648 //=======================================================================
649 //function : readResultFile
651 //=======================================================================
653 static bool readResultFile(const int fileOpen,
655 const char* fileName,
658 TopoDS_Shape tabShape[],
661 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
671 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
674 int nbElems, nbNodes, nbInputNodes;
675 int nodeId/*, triangleId*/;
677 int ID, shapeID, ghs3dShapeID;
680 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
682 int *tab, *tabID, *nodeID, *nodeAssigne;
684 const SMDS_MeshNode **node;
687 //tabID = new int[nbShape];
689 coord = new double[3];
690 node = new const SMDS_MeshNode*[4];
693 SMDS_MeshNode * aNewNode;
694 map <int,const SMDS_MeshNode*>::iterator itOnNode;
695 SMDS_MeshElement* aTet;
700 // Read the file state
701 fileStat = fstat(fileOpen, &status);
702 length = status.st_size;
704 // Mapping the result file into memory
706 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
707 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
708 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
709 0, (DWORD)length, NULL);
710 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
712 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
716 ptr = readMapIntLine(ptr, tab);
721 nbInputNodes = tab[2];
723 nodeAssigne = new int[ nbNodes+1 ];
726 aSolid = tabShape[0];
728 // Reading the nodeId
729 for (int i=0; i < 4*nbElems; i++)
730 nodeId = strtol(ptr, &ptr, 10);
732 // Reading the nodeCoor and update the nodeMap
733 for (int iNode=1; iNode <= nbNodes; iNode++) {
734 for (int iCoor=0; iCoor < 3; iCoor++)
735 coord[ iCoor ] = strtod(ptr, &ptr);
736 nodeAssigne[ iNode ] = 1;
737 if ( iNode > nbInputNodes ) {
738 nodeAssigne[ iNode ] = 0;
739 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
740 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
744 // Reading the number of triangles which corresponds to the number of sub-domains
745 nbTriangle = strtol(ptr, &ptr, 10);
747 tabID = new int[nbTriangle];
748 for (int i=0; i < nbTriangle; i++) {
750 // find the solid corresponding to GHS3D sub-domain following
751 // the technique proposed in GHS3D manual in chapter
752 // "B.4 Subdomain (sub-region) assignment"
753 int nodeId1 = strtol(ptr, &ptr, 10);
754 int nodeId2 = strtol(ptr, &ptr, 10);
755 int nodeId3 = strtol(ptr, &ptr, 10);
756 if ( nbTriangle > 1 ) {
757 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
758 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
759 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
762 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
763 // -- 0020330: Pb with ghs3d as a submesh
764 // check that found shape is to be meshed
765 if ( tabID[i] > 0 ) {
766 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] );
767 bool isToBeMeshed = false;
768 for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS )
769 isToBeMeshed = foundShape.IsSame( tabShape[ iS ]);
773 // END -- 0020330: Pb with ghs3d as a submesh
775 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
777 } catch ( Standard_Failure ) {
784 if ( nbTriangle <= nbShape ) // no holes
785 toMeshHoles = true; // not avoid creating tetras in holes
787 // Associating the tetrahedrons to the shapes
788 shapeID = compoundID;
789 for (int iElem = 0; iElem < nbElems; iElem++) {
790 for (int iNode = 0; iNode < 4; iNode++) {
791 ID = strtol(tetraPtr, &tetraPtr, 10);
792 itOnNode = theGhs3dIdToNodeMap.find(ID);
793 node[ iNode ] = itOnNode->second;
794 nodeID[ iNode ] = ID;
796 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
797 // tetras within holes depending on hypo option,
798 // so we first check if aTet is inside a hole and then create it
799 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
800 if ( nbTriangle > 1 ) {
801 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
802 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
803 if ( tabID[ ghs3dShapeID ] == 0 ) {
805 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
806 if ( toMeshHoles || state == TopAbs_IN )
807 shapeID = theMeshDS->ShapeToIndex( aSolid );
808 tabID[ ghs3dShapeID ] = shapeID;
811 shapeID = tabID[ ghs3dShapeID ];
813 else if ( nbShape > 1 ) {
814 // Case where nbTriangle == 1 while nbShape == 2 encountered
815 // with compound of 2 boxes and "To mesh holes"==False,
816 // so there are no subdomains specified for each tetrahedron.
817 // Try to guess a solid by a node already bound to shape
819 for ( int i=0; i<4 && shapeID==0; i++ ) {
820 if ( nodeAssigne[ nodeID[i] ] == 1 &&
821 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
822 node[i]->GetPosition()->GetShapeId() > 1 )
824 shapeID = node[i]->GetPosition()->GetShapeId();
828 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
829 shapeID = theMeshDS->ShapeToIndex( aSolid );
832 // set new nodes and tetrahedron onto the shape
833 for ( int i=0; i<4; i++ ) {
834 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
835 if ( shapeID != HOLE_ID )
836 theMeshDS->SetNodeInVolume( node[i], shapeID );
837 nodeAssigne[ nodeID[i] ] = shapeID;
840 if ( toMeshHoles || shapeID != HOLE_ID ) {
841 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
842 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
845 shapeIDs.insert( shapeID );
848 // Remove nodes of tetras inside holes if !toMeshHoles
849 if ( !toMeshHoles ) {
850 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
851 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
852 ID = itOnNode->first;
853 if ( nodeAssigne[ ID ] == HOLE_ID )
854 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
859 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
861 UnmapViewOfFile(mapPtr);
862 CloseHandle(hMapObject);
865 munmap(mapPtr, length);
874 delete [] nodeAssigne;
877 if ( shapeIDs.size() != nbShape ) {
878 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
879 for (int i=0; i<nbShape; i++) {
880 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
881 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
882 cout << " Solid #" << shapeID << " not found" << endl;
890 //=======================================================================
891 //function : readResultFile
893 //=======================================================================
895 static bool readResultFile(const int fileOpen,
897 const char* fileName,
899 SMESHDS_Mesh* theMeshDS,
901 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
911 int nbElems, nbNodes, nbInputNodes;
912 int nodeId, triangleId;
918 const SMDS_MeshNode **node;
921 coord = new double[3];
922 node = new const SMDS_MeshNode*[4];
924 SMDS_MeshNode * aNewNode;
925 map <int,const SMDS_MeshNode*>::iterator IdNode;
926 SMDS_MeshElement* aTet;
928 // Read the file state
929 fileStat = fstat(fileOpen, &status);
930 length = status.st_size;
932 // Mapping the result file into memory
934 HANDLE fd = CreateFile(fileName, GENERIC_READ, FILE_SHARE_READ,
935 NULL, OPEN_EXISTING, FILE_ATTRIBUTE_NORMAL, NULL);
936 HANDLE hMapObject = CreateFileMapping(fd, NULL, PAGE_READONLY,
937 0, (DWORD)length, NULL);
938 ptr = ( char* ) MapViewOfFile(hMapObject, FILE_MAP_READ, 0, 0, 0 );
940 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
944 ptr = readMapIntLine(ptr, tab);
949 nbInputNodes = tab[2];
951 theNodeByGhs3dId.resize( nbNodes );
953 // Reading the nodeId
954 for (int i=0; i < 4*nbElems; i++)
955 nodeId = strtol(ptr, &ptr, 10);
957 // Reading the nodeCoor and update the nodeMap
958 shapeID = theMeshDS->ShapeToIndex( aSolid );
959 for (int iNode=0; iNode < nbNodes; iNode++) {
960 for (int iCoor=0; iCoor < 3; iCoor++)
961 coord[ iCoor ] = strtod(ptr, &ptr);
962 if ((iNode+1) > nbInputNodes) {
963 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
964 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
965 theNodeByGhs3dId[ iNode ] = aNewNode;
969 // Reading the triangles
970 nbTriangle = strtol(ptr, &ptr, 10);
972 for (int i=0; i < 3*nbTriangle; i++)
973 triangleId = strtol(ptr, &ptr, 10);
977 // Associating the tetrahedrons to the shapes
978 for (int iElem = 0; iElem < nbElems; iElem++) {
979 for (int iNode = 0; iNode < 4; iNode++) {
980 ID = strtol(tetraPtr, &tetraPtr, 10);
981 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
983 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
984 shapeID = theMeshDS->ShapeToIndex( aSolid );
985 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
988 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
990 UnmapViewOfFile(mapPtr);
991 CloseHandle(hMapObject);
994 munmap(mapPtr, length);
1005 //=============================================================================
1007 *Here we are going to use the GHS3D mesher
1009 //=============================================================================
1011 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1012 const TopoDS_Shape& theShape)
1015 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1017 // we count the number of shapes
1018 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
1020 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1021 for ( ; expBox.More(); expBox.Next() )
1024 // create bounding box for every shape inside the compound
1027 TopoDS_Shape* tabShape;
1029 tabShape = new TopoDS_Shape[_nbShape];
1030 tabBox = new double*[_nbShape];
1031 for (int i=0; i<_nbShape; i++)
1032 tabBox[i] = new double[6];
1033 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1035 // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh
1036 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
1037 tabShape[iShape] = expBox.Current();
1038 Bnd_Box BoundingBox;
1039 BRepBndLib::Add(expBox.Current(), BoundingBox);
1040 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1041 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1042 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1043 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1047 // a unique working file name
1048 // to avoid access to the same files by eg different users
1049 TCollection_AsciiString aGenericName
1050 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1052 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1053 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1054 aFacesFileName = aGenericName + ".faces"; // in faces
1055 aPointsFileName = aGenericName + ".points"; // in points
1056 aResultFileName = aGenericName + ".noboite";// out points and volumes
1057 aBadResFileName = aGenericName + ".boite"; // out bad result
1058 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1059 aLogFileName = aGenericName + ".log"; // log
1061 // -----------------
1063 // -----------------
1065 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1066 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1069 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1071 INFOS( "Can't write into " << aFacesFileName);
1072 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1074 map <int,int> aSmdsToGhs3dIdMap;
1075 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1077 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1078 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1081 aPointsFile.close();
1084 if ( !_keepFiles ) {
1085 removeFile( aFacesFileName );
1086 removeFile( aPointsFileName );
1088 return error(COMPERR_BAD_INPUT_MESH);
1090 removeFile( aResultFileName ); // needed for boundary recovery module usage
1092 // -----------------
1094 // -----------------
1096 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1097 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1098 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1101 cout << "Ghs3d execution..." << endl;
1102 cout << cmd << endl;
1104 system( cmd.ToCString() ); // run
1107 cout << "End of Ghs3d execution !" << endl;
1113 // Mapping the result file
1116 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1117 if ( fileOpen < 0 ) {
1119 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1120 cout << "Log: " << aLogFileName << endl;
1125 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1126 Ok = readResultFile( fileOpen,
1128 aResultFileName.ToCString(),
1130 theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1134 // ---------------------
1135 // remove working files
1136 // ---------------------
1141 removeFile( aLogFileName );
1143 else if ( OSD_File( aLogFileName ).Size() > 0 )
1145 // get problem description from the log file
1146 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1147 storeErrorDescription( aLogFileName, conv );
1151 // the log file is empty
1152 removeFile( aLogFileName );
1153 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1154 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1157 if ( !_keepFiles ) {
1158 removeFile( aFacesFileName );
1159 removeFile( aPointsFileName );
1160 removeFile( aResultFileName );
1161 removeFile( aBadResFileName );
1162 removeFile( aBbResFileName );
1164 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1167 cout << "treated !" << endl;
1170 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1177 //=============================================================================
1179 *Here we are going to use the GHS3D mesher w/o geometry
1181 //=============================================================================
1182 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1183 SMESH_MesherHelper* aHelper)
1185 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1187 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1188 TopoDS_Shape theShape = aHelper->GetSubShape();
1190 // a unique working file name
1191 // to avoid access to the same files by eg different users
1192 TCollection_AsciiString aGenericName
1193 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1195 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1196 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1197 aFacesFileName = aGenericName + ".faces"; // in faces
1198 aPointsFileName = aGenericName + ".points"; // in points
1199 aResultFileName = aGenericName + ".noboite";// out points and volumes
1200 aBadResFileName = aGenericName + ".boite"; // out bad result
1201 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1202 aLogFileName = aGenericName + ".log"; // log
1204 // -----------------
1206 // -----------------
1208 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1209 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1211 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1214 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1216 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1218 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1219 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1222 aPointsFile.close();
1225 if ( !_keepFiles ) {
1226 removeFile( aFacesFileName );
1227 removeFile( aPointsFileName );
1229 return error(COMPERR_BAD_INPUT_MESH);
1231 removeFile( aResultFileName ); // needed for boundary recovery module usage
1233 // -----------------
1235 // -----------------
1237 TCollection_AsciiString cmd =
1238 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1239 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1240 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1242 system( cmd.ToCString() ); // run
1248 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1249 if ( fileOpen < 0 ) {
1251 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1252 cout << "Log: " << aLogFileName << endl;
1257 Ok = readResultFile( fileOpen,
1259 aResultFileName.ToCString(),
1261 meshDS, theShape ,aNodeByGhs3dId );
1264 // ---------------------
1265 // remove working files
1266 // ---------------------
1271 removeFile( aLogFileName );
1273 else if ( OSD_File( aLogFileName ).Size() > 0 )
1275 // get problem description from the log file
1276 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1277 storeErrorDescription( aLogFileName, conv );
1280 // the log file is empty
1281 removeFile( aLogFileName );
1282 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1283 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1288 removeFile( aFacesFileName );
1289 removeFile( aPointsFileName );
1290 removeFile( aResultFileName );
1291 removeFile( aBadResFileName );
1292 removeFile( aBbResFileName );
1298 //================================================================================
1300 * \brief Provide human readable text by error code reported by ghs3d
1302 //================================================================================
1304 static string translateError(const int errNum)
1308 return "The surface mesh includes a face of type other than edge, "
1309 "triangle or quadrilateral. This face type is not supported.";
1311 return "Not enough memory for the face table.";
1313 return "Not enough memory.";
1315 return "Not enough memory.";
1317 return "Face is ignored.";
1319 return "End of file. Some data are missing in the file.";
1321 return "Read error on the file. There are wrong data in the file.";
1323 return "the metric file is inadequate (dimension other than 3).";
1325 return "the metric file is inadequate (values not per vertices).";
1327 return "the metric file contains more than one field.";
1329 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1330 "value of number of mesh vertices in the \".noboite\" file.";
1332 return "Too many sub-domains.";
1334 return "the number of vertices is negative or null.";
1336 return "the number of faces is negative or null.";
1338 return "A face has a null vertex.";
1340 return "incompatible data.";
1342 return "the number of vertices is negative or null.";
1344 return "the number of vertices is negative or null (in the \".mesh\" file).";
1346 return "the number of faces is negative or null.";
1348 return "A face appears more than once in the input surface mesh.";
1350 return "An edge appears more than once in the input surface mesh.";
1352 return "A face has a vertex negative or null.";
1354 return "NOT ENOUGH MEMORY.";
1356 return "Not enough available memory.";
1358 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1359 "in terms of quality or the input list of points is wrong.";
1361 return "Some vertices are too close to one another or coincident.";
1363 return "Some vertices are too close to one another or coincident.";
1365 return "A vertex cannot be inserted.";
1367 return "There are at least two points considered as coincident.";
1369 return "Some vertices are too close to one another or coincident.";
1371 return "The surface mesh regeneration step has failed.";
1373 return "Constrained edge cannot be enforced.";
1375 return "Constrained face cannot be enforced.";
1377 return "Missing faces.";
1379 return "No guess to start the definition of the connected component(s).";
1381 return "The surface mesh includes at least one hole. The domain is not well defined.";
1383 return "Impossible to define a component.";
1385 return "The surface edge intersects another surface edge.";
1387 return "The surface edge intersects the surface face.";
1389 return "One boundary point lies within a surface face.";
1391 return "One surface edge intersects a surface face.";
1393 return "One boundary point lies within a surface edge.";
1395 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1396 "to too many swaps.";
1398 return "Edge is unique (i.e., bounds a hole in the surface).";
1400 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1402 return "Too many components, too many sub-domain.";
1404 return "The surface mesh includes at least one hole. "
1405 "Therefore there is no domain properly defined.";
1407 return "Statistics.";
1409 return "Statistics.";
1411 return "Warning, it is dramatically tedious to enforce the boundary items.";
1413 return "Not enough memory at this time, nevertheless, the program continues. "
1414 "The expected mesh will be correct but not really as large as required.";
1416 return "see above error code, resulting quality may be poor.";
1418 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1420 return "Unknown face type.";
1423 return "End of file. Some data are missing in the file.";
1425 return "A too small volume element is detected.";
1427 return "There exists at least a null or negative volume element.";
1429 return "There exist null or negative volume elements.";
1431 return "A too small volume element is detected. A face is considered being degenerated.";
1433 return "Some element is suspected to be very bad shaped or wrong.";
1435 return "A too bad quality face is detected. This face is considered degenerated.";
1437 return "A too bad quality face is detected. This face is degenerated.";
1439 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1441 return "Abnormal error occured, contact hotline.";
1443 return "Not enough memory for the face table.";
1445 return "The algorithm cannot run further. "
1446 "The surface mesh is probably very bad in terms of quality.";
1448 return "Bad vertex number.";
1453 //================================================================================
1455 * \brief Retrieve from a string given number of integers
1457 //================================================================================
1459 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1462 ids.reserve( nbIds );
1465 while ( !isdigit( *ptr )) ++ptr;
1466 if ( ptr[-1] == '-' ) --ptr;
1467 ids.push_back( strtol( ptr, &ptr, 10 ));
1473 //================================================================================
1475 * \brief Retrieve problem description form a log file
1476 * \retval bool - always false
1478 //================================================================================
1480 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1481 const _Ghs2smdsConvertor & toSmdsConvertor )
1485 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1487 int file = ::open (logFile.ToCString(), O_RDONLY);
1490 return error( SMESH_Comment("See ") << logFile << " for problem description");
1493 // struct stat status;
1494 // fstat(file, &status);
1495 // size_t length = status.st_size;
1496 off_t length = lseek( file, 0, SEEK_END);
1497 lseek( file, 0, SEEK_SET);
1500 vector< char > buf( length );
1501 int nBytesRead = ::read (file, & buf[0], length);
1503 char* ptr = & buf[0];
1504 char* bufEnd = ptr + nBytesRead;
1506 SMESH_Comment errDescription;
1508 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1510 // look for errors "ERR #"
1512 set<string> foundErrorStr; // to avoid reporting same error several times
1513 set<int> elemErrorNums; // not to report different types of errors with bad elements
1514 while ( ++ptr < bufEnd )
1516 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1519 list<const SMDS_MeshElement*> badElems;
1520 vector<int> nodeIds;
1524 int errNum = strtol(ptr, &ptr, 10);
1525 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1527 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1528 ptr = getIds(ptr, NODE, nodeIds);
1529 ptr = getIds(ptr, TRIA, nodeIds);
1530 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1532 case 1000: // ERR 1000 : 1 3 2
1533 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1534 ptr = getIds(ptr, TRIA, nodeIds);
1535 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1538 // Edge (e1, e2) appears more than once in the input surface mesh
1539 ptr = getIds(ptr, EDGE, nodeIds);
1540 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1543 // Face (f 1, f 2, f 3) has a vertex negative or null
1544 ptr = getIds(ptr, TRIA, nodeIds);
1545 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1548 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1549 ptr = getIds(ptr, NODE, nodeIds);
1550 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1551 ptr = getIds(ptr, NODE, nodeIds);
1552 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1555 // Vertex v1 cannot be inserted (warning).
1556 ptr = getIds(ptr, NODE, nodeIds);
1557 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1560 // There are at least two points whose distance is dist, i.e., considered as coincident
1561 case 2103: // ERR 2103 : 16 WITH 3
1562 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1563 ptr = getIds(ptr, NODE, nodeIds);
1564 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1565 ptr = getIds(ptr, NODE, nodeIds);
1566 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1569 // Constrained edge (e1, e2) cannot be enforced (warning).
1570 ptr = getIds(ptr, EDGE, nodeIds);
1571 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1574 // Constrained face (f 1, f 2, f 3) cannot be enforced
1575 ptr = getIds(ptr, TRIA, nodeIds);
1576 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1578 case 3103: // ERR 3103 : 1 2 WITH 7 3
1579 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1580 ptr = getIds(ptr, EDGE, nodeIds);
1581 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1582 ptr = getIds(ptr, EDGE, nodeIds);
1583 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1585 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1586 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1587 ptr = getIds(ptr, EDGE, nodeIds);
1588 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1589 ptr = getIds(ptr, TRIA, nodeIds);
1590 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1592 case 3105: // ERR 3105 : 8 IN 2 3 5
1593 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1594 ptr = getIds(ptr, NODE, nodeIds);
1595 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1596 ptr = getIds(ptr, TRIA, nodeIds);
1597 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1600 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1601 ptr = getIds(ptr, EDGE, nodeIds);
1602 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1603 ptr = getIds(ptr, TRIA, nodeIds);
1604 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1606 case 3107: // ERR 3107 : 2 IN 4 1
1607 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1608 ptr = getIds(ptr, NODE, nodeIds);
1609 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1610 ptr = getIds(ptr, EDGE, nodeIds);
1611 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1613 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1614 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1615 ptr = getIds(ptr, EDGE, nodeIds);
1616 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1618 case 9000: // ERR 9000
1619 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1620 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1621 // A too small volume element is detected. Are reported the index of the element,
1622 // its four vertex indices, its volume and the tolerance threshold value
1623 ptr = getIds(ptr, ID, nodeIds);
1624 ptr = getIds(ptr, VOL, nodeIds);
1625 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1626 // even if all nodes found, volume it most probably invisible,
1627 // add its faces to demenstrate it anyhow
1629 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1630 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1631 faceNodes[2] = nodeIds[3]; // 013
1632 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1633 faceNodes[1] = nodeIds[2]; // 023
1634 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1635 faceNodes[0] = nodeIds[1]; // 123
1636 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1639 case 9001: // ERR 9001
1640 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1641 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1642 // %% NUMBER OF NULL VOLUME TETS : 0
1643 // There exists at least a null or negative volume element
1646 // There exist n null or negative volume elements
1649 // A too small volume element is detected
1652 // A too bad quality face is detected. This face is considered degenerated,
1653 // its index, its three vertex indices together with its quality value are reported
1654 break; // same as next
1655 case 9112: // ERR 9112
1656 // FACE 2 WITH VERTICES : 4 2 5
1657 // SMALL INRADIUS : 0.
1658 // A too bad quality face is detected. This face is degenerated,
1659 // its index, its three vertex indices together with its inradius are reported
1660 ptr = getIds(ptr, ID, nodeIds);
1661 ptr = getIds(ptr, TRIA, nodeIds);
1662 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1663 // add triangle edges as it most probably has zero area and hence invisible
1665 vector<int> edgeNodes(2);
1666 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1667 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1668 edgeNodes[1] = nodeIds[2]; // 0-2
1669 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1670 edgeNodes[0] = nodeIds[1]; // 1-2
1671 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1676 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1678 continue; // not to report same error several times
1680 // const SMDS_MeshElement* nullElem = 0;
1681 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1683 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1684 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1685 // if ( oneMoreErrorType )
1686 // continue; // not to report different types of errors with bad elements
1689 // store bad elements
1690 //if ( allElemsOk ) {
1691 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1692 for ( ; elem != badElems.end(); ++elem )
1693 addBadInputElement( *elem );
1697 string text = translateError( errNum );
1698 if ( errDescription.find( text ) == text.npos ) {
1699 if ( !errDescription.empty() )
1700 errDescription << "\n";
1701 errDescription << text;
1706 if ( errDescription.empty() ) { // no errors found
1707 char msg[] = "connection to server failed";
1708 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1709 errDescription << "Licence problems.";
1712 if ( errDescription.empty() )
1713 errDescription << "See " << logFile << " for problem description";
1715 errDescription << "\nSee " << logFile << " for more information";
1717 return error( errDescription );
1720 //================================================================================
1722 * \brief Creates _Ghs2smdsConvertor
1724 //================================================================================
1726 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1727 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1731 //================================================================================
1733 * \brief Creates _Ghs2smdsConvertor
1735 //================================================================================
1737 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1738 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1742 //================================================================================
1744 * \brief Return SMDS element by ids of GHS3D nodes
1746 //================================================================================
1748 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1750 size_t nbNodes = ghsNodes.size();
1751 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1752 for ( size_t i = 0; i < nbNodes; ++i ) {
1753 int ghsNode = ghsNodes[ i ];
1754 if ( _ghs2NodeMap ) {
1755 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1756 if ( in == _ghs2NodeMap->end() )
1758 nodes[ i ] = in->second;
1761 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1763 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1769 if ( nbNodes == 2 ) {
1770 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1772 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1775 if ( nbNodes == 3 ) {
1776 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1778 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1782 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1788 //=============================================================================
1792 //=============================================================================
1793 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
1794 const TopoDS_Shape& aShape,
1795 MapShapeNbElems& aResMap)
1797 int nbtri = 0, nbqua = 0;
1798 double fullArea = 0.0;
1799 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1800 TopoDS_Face F = TopoDS::Face( exp.Current() );
1801 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1802 MapShapeNbElemsItr anIt = aResMap.find(sm);
1803 if( anIt==aResMap.end() ) {
1804 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1805 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1806 "Submesh can not be evaluated",this));
1809 std::vector<int> aVec = (*anIt).second;
1810 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
1811 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1813 BRepGProp::SurfaceProperties(F,G);
1814 double anArea = G.Mass();
1818 // collect info from edges
1819 int nb0d_e = 0, nb1d_e = 0;
1820 bool IsQuadratic = false;
1821 bool IsFirst = true;
1822 TopTools_MapOfShape tmpMap;
1823 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1824 TopoDS_Edge E = TopoDS::Edge(exp.Current());
1825 if( tmpMap.Contains(E) )
1828 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
1829 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
1830 std::vector<int> aVec = (*anIt).second;
1831 nb0d_e += aVec[SMDSEntity_Node];
1832 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
1834 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1840 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
1843 BRepGProp::VolumeProperties(aShape,G);
1844 double aVolume = G.Mass();
1845 double tetrVol = 0.1179*ELen*ELen*ELen;
1846 double CoeffQuality = 0.9;
1847 int nbVols = int(aVolume/tetrVol/CoeffQuality);
1848 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
1849 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
1850 std::vector<int> aVec(SMDSEntity_Last);
1851 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1853 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
1854 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
1855 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
1858 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
1859 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
1860 aVec[SMDSEntity_Pyramid] = nbqua;
1862 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1863 aResMap.insert(std::make_pair(sm,aVec));