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 //=============================================================================
29 #include "GHS3DPlugin_GHS3D.hxx"
30 #include "GHS3DPlugin_Hypothesis.hxx"
32 #include "SMESH_Gen.hxx"
33 #include "SMESH_Mesh.hxx"
34 #include "SMESH_Comment.hxx"
35 #include "SMESH_MesherHelper.hxx"
36 #include "SMESH_MeshEditor.hxx"
38 #include "SMDS_MeshElement.hxx"
39 #include "SMDS_MeshNode.hxx"
40 #include "SMDS_FaceOfNodes.hxx"
41 #include "SMDS_VolumeOfNodes.hxx"
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepClass3d_SolidClassifier.hxx>
46 #include <BRepTools.hxx>
47 #include <BRep_Tool.hxx>
48 #include <Bnd_Box.hxx>
49 #include <GeomAPI_ProjectPointOnSurf.hxx>
50 #include <OSD_File.hxx>
51 #include <Precision.hxx>
52 #include <Quantity_Parameter.hxx>
53 #include <Standard_ErrorHandler.hxx>
54 #include <Standard_Failure.hxx>
56 #include <TopExp_Explorer.hxx>
57 #include <TopTools_IndexedMapOfShape.hxx>
58 #include <TopTools_ListIteratorOfListOfShape.hxx>
60 //#include <BRepClass_FaceClassifier.hxx>
61 //#include <BRepGProp.hxx>
62 //#include <GProp_GProps.hxx>
64 #include "utilities.h"
67 #include <sys/sysinfo.h>
70 //#include <Standard_Stream.hxx>
73 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
94 //=============================================================================
98 //=============================================================================
100 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
101 : SMESH_3D_Algo(hypId, studyId, gen)
103 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
105 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
106 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
109 _compatibleHypothesis.push_back("GHS3D_Parameters");
110 _requireShape = false; // can work without shape
113 //=============================================================================
117 //=============================================================================
119 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
121 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
124 //=============================================================================
128 //=============================================================================
130 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
131 const TopoDS_Shape& aShape,
132 Hypothesis_Status& aStatus )
134 aStatus = SMESH_Hypothesis::HYP_OK;
136 // there is only one compatible Hypothesis so far
139 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
141 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
143 _keepFiles = _hyp->GetKeepFiles();
148 //=======================================================================
149 //function : findShape
151 //=======================================================================
153 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
155 const TopoDS_Shape shape[],
158 TopAbs_State * state = 0)
161 int j, iShape, nbNode = 4;
163 for ( j=0; j<nbNode; j++ )
164 aPnt += gp_XYZ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
167 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
168 if (state) *state = SC.State();
169 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
170 for (iShape = 0; iShape < nShape; iShape++) {
171 aShape = shape[iShape];
172 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
173 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
174 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
175 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
176 if (state) *state = SC.State();
177 if (SC.State() == TopAbs_IN)
185 //=======================================================================
186 //function : readMapIntLine
188 //=======================================================================
190 static char* readMapIntLine(char* ptr, int tab[]) {
194 for ( int i=0; i<17; i++ ) {
195 intVal = strtol(ptr, &ptr, 10);
202 //=======================================================================
203 //function : countShape
205 //=======================================================================
207 template < class Mesh, class Shape >
208 static int countShape( Mesh* mesh, Shape shape ) {
209 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
211 for ( ; expShape.More(); expShape.Next() ) {
217 //=======================================================================
218 //function : writeFaces
220 //=======================================================================
222 static bool writeFaces (ofstream & theFile,
223 SMESHDS_Mesh * theMesh,
224 const map <int,int> & theSmdsToGhs3dIdMap)
228 // NB_ELEMS DUMMY_INT
229 // Loop from 1 to NB_ELEMS
230 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
232 int nbShape = countShape( theMesh, TopAbs_FACE );
234 int *tabID; tabID = new int[nbShape];
235 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
237 SMESHDS_SubMesh* theSubMesh;
238 const SMDS_MeshElement* aFace;
239 const char* space = " ";
240 const int dummyint = 0;
241 map<int,int>::const_iterator itOnMap;
242 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
243 int shapeID, nbNodes, aSmdsID;
246 cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
249 theFile << space << theMesh->NbFaces() << space << dummyint << endl;
251 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
252 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
254 aShape = expface.Current();
255 shapeID = theMesh->ShapeToIndex( aShape );
257 for ( int j=0; j<=i; j++) {
258 if ( shapeID == tabID[j] ) {
265 tabShape[i] = aShape;
268 for ( int i =0; i < nbShape; i++ ) {
269 if ( tabID[i] != 0 ) {
270 aShape = tabShape[i];
272 theSubMesh = theMesh->MeshElements( aShape );
273 if ( !theSubMesh ) continue;
274 itOnSubMesh = theSubMesh->GetElements();
275 while ( itOnSubMesh->more() ) {
276 aFace = itOnSubMesh->next();
277 nbNodes = aFace->NbNodes();
279 theFile << space << nbNodes;
281 itOnSubFace = aFace->nodesIterator();
282 while ( itOnSubFace->more() ) {
284 aSmdsID = itOnSubFace->next()->GetID();
285 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
286 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
288 theFile << space << (*itOnMap).second;
291 // (NB_NODES + 1) times: DUMMY_INT
292 for ( int j=0; j<=nbNodes; j++)
293 theFile << space << dummyint;
306 //=======================================================================
307 //function : writeFaces
308 //purpose : Write Faces in case if generate 3D mesh w/o geometry
309 //=======================================================================
311 static bool writeFaces (ofstream & theFile,
312 SMESHDS_Mesh * theMesh,
313 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
317 // NB_ELEMS DUMMY_INT
318 // Loop from 1 to NB_ELEMS
319 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
323 list< const SMDS_MeshElement* > faces;
324 list< const SMDS_MeshElement* >::iterator f;
325 map< const SMDS_MeshNode*,int >::iterator it;
326 SMDS_ElemIteratorPtr nodeIt;
327 const SMDS_MeshElement* elem;
330 const char* space = " ";
331 const int dummyint = 0;
333 //get all faces from mesh
334 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
335 while ( eIt->more() ) {
336 const SMDS_MeshElement* elem = eIt->next();
339 faces.push_back( elem );
346 cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
348 // NB_ELEMS DUMMY_INT
349 theFile << space << nbFaces << space << dummyint << endl;
351 // Loop from 1 to NB_ELEMS
353 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
355 for ( ; f != faces.end(); ++f )
359 nbNodes = elem->NbNodes();
360 theFile << space << nbNodes;
362 // NODE_NB_1 NODE_NB_2 ...
363 nodeIt = elem->nodesIterator();
364 while ( nodeIt->more() )
367 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
368 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
369 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
370 theFile << space << it->second;
373 // (NB_NODES + 1) times: DUMMY_INT
374 for ( int i=0; i<=nbNodes; i++)
375 theFile << space << dummyint;
380 // put nodes to theNodeByGhs3dId vector
381 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
382 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
383 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
385 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
391 //=======================================================================
392 //function : writePoints
394 //=======================================================================
396 static bool writePoints (ofstream & theFile,
397 SMESHDS_Mesh * theMesh,
398 map <int,int> & theSmdsToGhs3dIdMap,
399 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
404 // Loop from 1 to NB_NODES
407 int nbNodes = theMesh->NbNodes();
411 const char* space = " ";
412 const int dummyint = 0;
415 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
416 const SMDS_MeshNode* node;
419 theFile << space << nbNodes << endl;
421 cout << "The initial 2D mesh contains :" << endl;
422 cout << " " << nbNodes << " vertices" << endl;
424 // Loop from 1 to NB_NODES
429 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
430 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
435 << space << node->X()
436 << space << node->Y()
437 << space << node->Z()
438 << space << dummyint;
446 //=======================================================================
447 //function : writePoints
449 //=======================================================================
451 static bool writePoints (ofstream & theFile,
452 SMESHDS_Mesh * theMesh,
453 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
458 // Loop from 1 to NB_NODES
461 //int nbNodes = theMesh->NbNodes();
462 int nbNodes = theNodeByGhs3dId.size();
466 const char* space = " ";
467 const int dummyint = 0;
469 const SMDS_MeshNode* node;
472 theFile << space << nbNodes << endl;
473 cout << nbNodes << " nodes" << endl;
475 // Loop from 1 to NB_NODES
477 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
478 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
479 for ( ; nodeIt != after; ++nodeIt )
485 << space << node->X()
486 << space << node->Y()
487 << space << node->Z()
488 << space << dummyint;
496 //=======================================================================
497 //function : findShapeID
498 //purpose : find the solid corresponding to GHS3D sub-domain following
499 // the technique proposed in GHS3D manual in chapter
500 // "B.4 Subdomain (sub-region) assignment"
501 //=======================================================================
503 static int findShapeID(SMESH_Mesh& mesh,
504 const SMDS_MeshNode* node1,
505 const SMDS_MeshNode* node2,
506 const SMDS_MeshNode* node3,
507 const bool toMeshHoles)
509 const int invalidID = 0;
510 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
512 // face th enodes belong to
513 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
517 // geom face the face assigned to
518 SMESH_MeshEditor editor(&mesh);
519 int geomFaceID = editor.FindShape( face );
522 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
523 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
525 TopoDS_Face geomFace = TopoDS::Face( shape );
527 // solids bounded by geom face
528 TopTools_IndexedMapOfShape solids, shells;
529 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
530 for ( ; ansIt.More(); ansIt.Next() ) {
531 switch ( ansIt.Value().ShapeType() ) {
533 solids.Add( ansIt.Value() ); break;
535 shells.Add( ansIt.Value() ); break;
539 // analyse found solids
540 if ( solids.Extent() == 0 || shells.Extent() == 0)
543 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
544 if ( solids.Extent() == 1 )
547 return meshDS->ShapeToIndex( solid1 );
549 // - Are we at a hole boundary face?
550 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
551 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
553 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
554 // check if any edge of shells(1) belongs to another shell
555 for ( ; eExp.More() && !touch; eExp.Next() ) {
556 ansIt = mesh.GetAncestors( eExp.Current() );
557 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
558 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
559 touch = ( !ansIt.Value().IsSame( shells(1) ));
563 return meshDS->ShapeToIndex( solid1 );
566 // find orientation of geom face within the first solid
567 TopExp_Explorer fExp( solid1, TopAbs_FACE );
568 for ( ; fExp.More(); fExp.Next() )
569 if ( geomFace.IsSame( fExp.Current() )) {
570 geomFace = TopoDS::Face( fExp.Current() );
574 return invalidID; // face not found
576 // find UV of the center of triangle on geomFace
577 SMESH_MesherHelper helper( mesh );
578 BRepAdaptor_Surface surface( geomFace );
580 const SMDS_MeshNode* nodes[3] = { node1, node2, node3 };
581 const SMDS_MeshNode* nNotOnSeamEdge = 0;
582 for ( int n = 0; !nNotOnSeamEdge && n < 3; ++n )
583 if ( !helper.IsSeamShape( nodes[n]->GetPosition()->GetShapeId() ))
584 nNotOnSeamEdge = nodes[n];
585 for ( int n = 0; n < 3; ++n )
587 const SMDS_MeshNode* node = nodes[n];
588 gp_XY uv = helper.GetNodeUV( geomFace, node, nNotOnSeamEdge );
590 // check that uv is correct
591 gp_Pnt nodePnt ( node->X(), node->Y(), node->Z() );
592 double tol = BRep_Tool::Tolerance( geomFace );
593 if ( nodePnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
594 // project node onto geomFace to get right UV
595 GeomAPI_ProjectPointOnSurf projector( nodePnt, surface.Surface().Surface() );
596 if ( !projector.IsDone() || projector.NbPoints() < 1 )
598 Quantity_Parameter U,V;
599 projector.LowerDistanceParameters(U,V);
604 // normale to face at node1
605 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
606 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
607 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
608 gp_Vec vec12( node1Pnt, node2Pnt );
609 gp_Vec vec13( node1Pnt, node3Pnt );
610 gp_Vec meshNormal = vec12 ^ vec13;
611 if ( meshNormal.SquareMagnitude() < DBL_MIN )
614 // normale to geomFace at UV
616 surface.D1( UV.X(), UV.Y(), node1Pnt, du, dv );
617 gp_Vec geomNormal = du ^ dv;
618 if ( geomNormal.SquareMagnitude() < DBL_MIN )
619 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
620 if ( geomFace.Orientation() == TopAbs_REVERSED )
621 geomNormal.Reverse();
624 bool isReverse = ( meshNormal * geomNormal ) < 0;
626 return meshDS->ShapeToIndex( solid1 );
628 if ( solids.Extent() == 1 )
629 return HOLE_ID; // we are inside a hole
631 return meshDS->ShapeToIndex( solids(2) );
634 //=======================================================================
635 //function : readResultFile
637 //=======================================================================
639 static bool readResultFile(const int fileOpen,
641 TopoDS_Shape tabShape[],
644 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
654 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
657 int nbElems, nbNodes, nbInputNodes;
658 int nodeId/*, triangleId*/;
660 int ID, shapeID, ghs3dShapeID;
663 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
665 int *tab, *tabID, *nodeID, *nodeAssigne;
667 const SMDS_MeshNode **node;
670 //tabID = new int[nbShape];
672 coord = new double[3];
673 node = new const SMDS_MeshNode*[4];
676 SMDS_MeshNode * aNewNode;
677 map <int,const SMDS_MeshNode*>::iterator itOnNode;
678 SMDS_MeshElement* aTet;
683 // Read the file state
684 fileStat = fstat(fileOpen, &status);
685 length = status.st_size;
687 // Mapping the result file into memory
688 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
691 ptr = readMapIntLine(ptr, tab);
696 nbInputNodes = tab[2];
698 nodeAssigne = new int[ nbNodes+1 ];
701 aSolid = tabShape[0];
703 // Reading the nodeId
704 for (int i=0; i < 4*nbElems; i++)
705 nodeId = strtol(ptr, &ptr, 10);
707 // Reading the nodeCoor and update the nodeMap
708 for (int iNode=1; iNode <= nbNodes; iNode++) {
709 for (int iCoor=0; iCoor < 3; iCoor++)
710 coord[ iCoor ] = strtod(ptr, &ptr);
711 nodeAssigne[ iNode ] = 1;
712 if ( iNode > nbInputNodes ) {
713 nodeAssigne[ iNode ] = 0;
714 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
715 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
719 // Reading the number of triangles which corresponds to the number of sub-domains
720 nbTriangle = strtol(ptr, &ptr, 10);
722 tabID = new int[nbTriangle];
723 for (int i=0; i < nbTriangle; i++) {
725 // find the solid corresponding to GHS3D sub-domain following
726 // the technique proposed in GHS3D manual in chapter
727 // "B.4 Subdomain (sub-region) assignment"
728 int nodeId1 = strtol(ptr, &ptr, 10);
729 int nodeId2 = strtol(ptr, &ptr, 10);
730 int nodeId3 = strtol(ptr, &ptr, 10);
731 if ( nbTriangle > 1 ) {
732 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
733 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
734 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
737 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
739 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
741 } catch ( Standard_Failure ) {
748 if ( nbTriangle <= nbShape ) // no holes
749 toMeshHoles = true; // not avoid creating tetras in holes
751 // Associating the tetrahedrons to the shapes
752 shapeID = compoundID;
753 for (int iElem = 0; iElem < nbElems; iElem++) {
754 for (int iNode = 0; iNode < 4; iNode++) {
755 ID = strtol(tetraPtr, &tetraPtr, 10);
756 itOnNode = theGhs3dIdToNodeMap.find(ID);
757 node[ iNode ] = itOnNode->second;
758 nodeID[ iNode ] = ID;
760 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
761 // tetras within holes depending on hypo option,
762 // so we first check if aTet is inside a hole and then create it
763 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
764 if ( nbTriangle > 1 ) {
765 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
766 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
767 if ( tabID[ ghs3dShapeID ] == 0 ) {
769 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
770 if ( toMeshHoles || state == TopAbs_IN )
771 shapeID = theMeshDS->ShapeToIndex( aSolid );
772 tabID[ ghs3dShapeID ] = shapeID;
775 shapeID = tabID[ ghs3dShapeID ];
777 else if ( nbShape > 1 ) {
778 // Case where nbTriangle == 1 while nbShape == 2 encountered
779 // with compound of 2 boxes and "To mesh holes"==False,
780 // so there are no subdomains specified for each tetrahedron.
781 // Try to guess a solid by a node already bound to shape
783 for ( int i=0; i<4 && shapeID==0; i++ ) {
784 if ( nodeAssigne[ nodeID[i] ] == 1 &&
785 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
786 node[i]->GetPosition()->GetShapeId() > 1 )
788 shapeID = node[i]->GetPosition()->GetShapeId();
792 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
793 shapeID = theMeshDS->ShapeToIndex( aSolid );
796 // set new nodes and tetrahedron onto the shape
797 for ( int i=0; i<4; i++ ) {
798 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
799 if ( shapeID != HOLE_ID )
800 theMeshDS->SetNodeInVolume( node[i], shapeID );
801 nodeAssigne[ nodeID[i] ] = shapeID;
804 if ( toMeshHoles || shapeID != HOLE_ID ) {
805 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
806 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
809 shapeIDs.insert( shapeID );
812 // Remove nodes of tetras inside holes if !toMeshHoles
813 if ( !toMeshHoles ) {
814 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
815 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
816 ID = itOnNode->first;
817 if ( nodeAssigne[ ID ] == HOLE_ID )
818 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
823 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
824 munmap(mapPtr, length);
832 delete [] nodeAssigne;
835 if ( shapeIDs.size() != nbShape ) {
836 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
837 for (int i=0; i<nbShape; i++) {
838 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
839 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
840 cout << " Solid #" << shapeID << " not found" << endl;
848 //=======================================================================
849 //function : readResultFile
851 //=======================================================================
853 static bool readResultFile(const int fileOpen,
854 SMESHDS_Mesh* theMeshDS,
856 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
866 int nbElems, nbNodes, nbInputNodes;
867 int nodeId, triangleId;
873 const SMDS_MeshNode **node;
876 coord = new double[3];
877 node = new const SMDS_MeshNode*[4];
879 SMDS_MeshNode * aNewNode;
880 map <int,const SMDS_MeshNode*>::iterator IdNode;
881 SMDS_MeshElement* aTet;
883 // Read the file state
884 fileStat = fstat(fileOpen, &status);
885 length = status.st_size;
887 // Mapping the result file into memory
888 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
891 ptr = readMapIntLine(ptr, tab);
896 nbInputNodes = tab[2];
898 theNodeByGhs3dId.resize( nbNodes );
900 // Reading the nodeId
901 for (int i=0; i < 4*nbElems; i++)
902 nodeId = strtol(ptr, &ptr, 10);
904 // Reading the nodeCoor and update the nodeMap
905 shapeID = theMeshDS->ShapeToIndex( aSolid );
906 for (int iNode=0; iNode < nbNodes; iNode++) {
907 for (int iCoor=0; iCoor < 3; iCoor++)
908 coord[ iCoor ] = strtod(ptr, &ptr);
909 if ((iNode+1) > nbInputNodes) {
910 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
911 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
912 theNodeByGhs3dId[ iNode ] = aNewNode;
916 // Reading the triangles
917 nbTriangle = strtol(ptr, &ptr, 10);
919 for (int i=0; i < 3*nbTriangle; i++)
920 triangleId = strtol(ptr, &ptr, 10);
924 // Associating the tetrahedrons to the shapes
925 for (int iElem = 0; iElem < nbElems; iElem++) {
926 for (int iNode = 0; iNode < 4; iNode++) {
927 ID = strtol(tetraPtr, &tetraPtr, 10);
928 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
930 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
931 shapeID = theMeshDS->ShapeToIndex( aSolid );
932 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
935 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
936 munmap(mapPtr, length);
946 //=============================================================================
948 *Here we are going to use the GHS3D mesher
950 //=============================================================================
952 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
953 const TopoDS_Shape& theShape)
956 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
958 _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
960 // create bounding box for every shape inside the compound
963 TopoDS_Shape* tabShape;
965 tabShape = new TopoDS_Shape[_nbShape];
966 tabBox = new double*[_nbShape];
967 for (int i=0; i<_nbShape; i++)
968 tabBox[i] = new double[6];
969 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
971 TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
972 for (; expBox.More(); expBox.Next()) {
973 tabShape[iShape] = expBox.Current();
975 BRepBndLib::Add(expBox.Current(), BoundingBox);
976 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
977 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
978 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
979 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
983 // a unique working file name
984 // to avoid access to the same files by eg different users
985 TCollection_AsciiString aGenericName
986 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
988 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
989 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
990 aFacesFileName = aGenericName + ".faces"; // in faces
991 aPointsFileName = aGenericName + ".points"; // in points
992 aResultFileName = aGenericName + ".noboite";// out points and volumes
993 aBadResFileName = aGenericName + ".boite"; // out bad result
994 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
995 aLogFileName = aGenericName + ".log"; // log
1001 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1002 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1006 aFacesFile->is_open() && aPointsFile->is_open();
1008 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1011 INFOS( "Can't write into " << aFacesFileName);
1012 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1014 map <int,int> aSmdsToGhs3dIdMap;
1015 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1017 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1018 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1021 aPointsFile.close();
1024 if ( !_keepFiles ) {
1025 OSD_File( aFacesFileName ).Remove();
1026 OSD_File( aPointsFileName ).Remove();
1028 return error(COMPERR_BAD_INPUT_MESH);
1030 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1032 // -----------------
1034 // -----------------
1036 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1037 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1038 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1041 cout << "Ghs3d execution..." << endl;
1042 cout << cmd << endl;
1044 system( cmd.ToCString() ); // run
1047 cout << "End of Ghs3d execution !" << endl;
1053 // Mapping the result file
1056 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1057 if ( fileOpen < 0 ) {
1059 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1060 cout << "Log: " << aLogFileName << endl;
1065 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1066 Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1070 // ---------------------
1071 // remove working files
1072 // ---------------------
1077 OSD_File( aLogFileName ).Remove();
1079 else if ( OSD_File( aLogFileName ).Size() > 0 )
1081 // get problem description from the log file
1082 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1083 storeErrorDescription( aLogFileName, conv );
1087 // the log file is empty
1088 OSD_File( aLogFileName ).Remove();
1089 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1090 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1093 if ( !_keepFiles ) {
1094 OSD_File( aFacesFileName ).Remove();
1095 OSD_File( aPointsFileName ).Remove();
1096 OSD_File( aResultFileName ).Remove();
1097 OSD_File( aBadResFileName ).Remove();
1098 OSD_File( aBbResFileName ).Remove();
1100 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1103 cout << "treated !" << endl;
1106 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1113 //=============================================================================
1115 *Here we are going to use the GHS3D mesher w/o geometry
1117 //=============================================================================
1118 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1119 SMESH_MesherHelper* aHelper)
1121 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1123 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1124 TopoDS_Shape theShape = aHelper->GetSubShape();
1126 // a unique working file name
1127 // to avoid access to the same files by eg different users
1128 TCollection_AsciiString aGenericName
1129 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1131 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1132 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1133 aFacesFileName = aGenericName + ".faces"; // in faces
1134 aPointsFileName = aGenericName + ".points"; // in points
1135 aResultFileName = aGenericName + ".noboite";// out points and volumes
1136 aBadResFileName = aGenericName + ".boite"; // out bad result
1137 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1138 aLogFileName = aGenericName + ".log"; // log
1140 // -----------------
1142 // -----------------
1144 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1145 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1148 aFacesFile->is_open() && aPointsFile->is_open();
1150 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1154 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1156 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1158 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1159 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1162 aPointsFile.close();
1165 if ( !_keepFiles ) {
1166 OSD_File( aFacesFileName ).Remove();
1167 OSD_File( aPointsFileName ).Remove();
1169 return error(COMPERR_BAD_INPUT_MESH);
1171 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1173 // -----------------
1175 // -----------------
1177 TCollection_AsciiString cmd =
1178 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1179 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1180 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1182 system( cmd.ToCString() ); // run
1188 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1189 if ( fileOpen < 0 ) {
1191 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1192 cout << "Log: " << aLogFileName << endl;
1197 Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1200 // ---------------------
1201 // remove working files
1202 // ---------------------
1207 OSD_File( aLogFileName ).Remove();
1209 else if ( OSD_File( aLogFileName ).Size() > 0 )
1211 // get problem description from the log file
1212 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1213 storeErrorDescription( aLogFileName, conv );
1216 // the log file is empty
1217 OSD_File( aLogFileName ).Remove();
1218 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1219 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1224 OSD_File( aFacesFileName ).Remove();
1225 OSD_File( aPointsFileName ).Remove();
1226 OSD_File( aResultFileName ).Remove();
1227 OSD_File( aBadResFileName ).Remove();
1228 OSD_File( aBbResFileName ).Remove();
1234 //================================================================================
1236 * \brief Provide human readable text by error code reported by ghs3d
1238 //================================================================================
1240 static string translateError(const int errNum)
1244 return "The surface mesh includes a face of type other than edge, "
1245 "triangle or quadrilateral. This face type is not supported.";
1247 return "Not enough memory for the face table.";
1249 return "Not enough memory.";
1251 return "Not enough memory.";
1253 return "Face is ignored.";
1255 return "End of file. Some data are missing in the file.";
1257 return "Read error on the file. There are wrong data in the file.";
1259 return "the metric file is inadequate (dimension other than 3).";
1261 return "the metric file is inadequate (values not per vertices).";
1263 return "the metric file contains more than one field.";
1265 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1266 "value of number of mesh vertices in the \".noboite\" file.";
1268 return "Too many sub-domains.";
1270 return "the number of vertices is negative or null.";
1272 return "the number of faces is negative or null.";
1274 return "A face has a null vertex.";
1276 return "incompatible data.";
1278 return "the number of vertices is negative or null.";
1280 return "the number of vertices is negative or null (in the \".mesh\" file).";
1282 return "the number of faces is negative or null.";
1284 return "A face appears more than once in the input surface mesh.";
1286 return "An edge appears more than once in the input surface mesh.";
1288 return "A face has a vertex negative or null.";
1290 return "NOT ENOUGH MEMORY.";
1292 return "Not enough available memory.";
1294 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1295 "in terms of quality or the input list of points is wrong.";
1297 return "Some vertices are too close to one another or coincident.";
1299 return "Some vertices are too close to one another or coincident.";
1301 return "A vertex cannot be inserted.";
1303 return "There are at least two points considered as coincident.";
1305 return "Some vertices are too close to one another or coincident.";
1307 return "The surface mesh regeneration step has failed.";
1309 return "Constrained edge cannot be enforced.";
1311 return "Constrained face cannot be enforced.";
1313 return "Missing faces.";
1315 return "No guess to start the definition of the connected component(s).";
1317 return "The surface mesh includes at least one hole. The domain is not well defined.";
1319 return "Impossible to define a component.";
1321 return "The surface edge intersects another surface edge.";
1323 return "The surface edge intersects the surface face.";
1325 return "One boundary point lies within a surface face.";
1327 return "One surface edge intersects a surface face.";
1329 return "One boundary point lies within a surface edge.";
1331 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1332 "to too many swaps.";
1334 return "Edge is unique (i.e., bounds a hole in the surface).";
1336 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1338 return "Too many components, too many sub-domain.";
1340 return "The surface mesh includes at least one hole. "
1341 "Therefore there is no domain properly defined.";
1343 return "Statistics.";
1345 return "Statistics.";
1347 return "Warning, it is dramatically tedious to enforce the boundary items.";
1349 return "Not enough memory at this time, nevertheless, the program continues. "
1350 "The expected mesh will be correct but not really as large as required.";
1352 return "see above error code, resulting quality may be poor.";
1354 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1356 return "Unknown face type.";
1359 return "End of file. Some data are missing in the file.";
1361 return "A too small volume element is detected.";
1363 return "There exists at least a null or negative volume element.";
1365 return "There exist null or negative volume elements.";
1367 return "A too small volume element is detected. A face is considered being degenerated.";
1369 return "Some element is suspected to be very bad shaped or wrong.";
1371 return "A too bad quality face is detected. This face is considered degenerated.";
1373 return "A too bad quality face is detected. This face is degenerated.";
1375 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1377 return "Abnormal error occured, contact hotline.";
1379 return "Not enough memory for the face table.";
1381 return "The algorithm cannot run further. "
1382 "The surface mesh is probably very bad in terms of quality.";
1384 return "Bad vertex number.";
1389 //================================================================================
1391 * \brief Retrieve from a string given number of integers
1393 //================================================================================
1395 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1398 ids.reserve( nbIds );
1401 while ( !isdigit( *ptr )) ++ptr;
1402 if ( ptr[-1] == '-' ) --ptr;
1403 ids.push_back( strtol( ptr, &ptr, 10 ));
1409 //================================================================================
1411 * \brief Retrieve problem description form a log file
1412 * \retval bool - always false
1414 //================================================================================
1416 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1417 const _Ghs2smdsConvertor & toSmdsConvertor )
1421 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1423 int file = ::open (logFile.ToCString(), O_RDONLY);
1426 return error( SMESH_Comment("See ") << logFile << " for problem description");
1429 // struct stat status;
1430 // fstat(file, &status);
1431 // size_t length = status.st_size;
1432 off_t length = lseek( file, 0, SEEK_END);
1433 lseek( file, 0, SEEK_SET);
1436 vector< char > buf( length );
1437 int nBytesRead = ::read (file, & buf[0], length);
1439 char* ptr = & buf[0];
1440 char* bufEnd = ptr + nBytesRead;
1442 SMESH_Comment errDescription;
1444 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1446 // look for errors "ERR #"
1448 set<string> foundErrorStr; // to avoid reporting same error several times
1449 set<int> elemErrorNums; // not to report different types of errors with bad elements
1450 while ( ++ptr < bufEnd )
1452 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1455 list<const SMDS_MeshElement*> badElems;
1456 vector<int> nodeIds;
1460 int errNum = strtol(ptr, &ptr, 10);
1461 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1463 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1464 ptr = getIds(ptr, NODE, nodeIds);
1465 ptr = getIds(ptr, TRIA, nodeIds);
1466 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1468 case 1000: // ERR 1000 : 1 3 2
1469 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1470 ptr = getIds(ptr, TRIA, nodeIds);
1471 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1474 // Edge (e1, e2) appears more than once in the input surface mesh
1475 ptr = getIds(ptr, EDGE, nodeIds);
1476 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1479 // Face (f 1, f 2, f 3) has a vertex negative or null
1480 ptr = getIds(ptr, TRIA, nodeIds);
1481 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1484 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1485 ptr = getIds(ptr, NODE, nodeIds);
1486 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1487 ptr = getIds(ptr, NODE, nodeIds);
1488 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1491 // Vertex v1 cannot be inserted (warning).
1492 ptr = getIds(ptr, NODE, nodeIds);
1493 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1496 // There are at least two points whose distance is dist, i.e., considered as coincident
1497 case 2103: // ERR 2103 : 16 WITH 3
1498 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1499 ptr = getIds(ptr, NODE, nodeIds);
1500 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1501 ptr = getIds(ptr, NODE, nodeIds);
1502 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1505 // Constrained edge (e1, e2) cannot be enforced (warning).
1506 ptr = getIds(ptr, EDGE, nodeIds);
1507 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1510 // Constrained face (f 1, f 2, f 3) cannot be enforced
1511 ptr = getIds(ptr, TRIA, nodeIds);
1512 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1514 case 3103: // ERR 3103 : 1 2 WITH 7 3
1515 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1516 ptr = getIds(ptr, EDGE, nodeIds);
1517 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1518 ptr = getIds(ptr, EDGE, nodeIds);
1519 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1521 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1522 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1523 ptr = getIds(ptr, EDGE, nodeIds);
1524 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1525 ptr = getIds(ptr, TRIA, nodeIds);
1526 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1528 case 3105: // ERR 3105 : 8 IN 2 3 5
1529 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1530 ptr = getIds(ptr, NODE, nodeIds);
1531 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1532 ptr = getIds(ptr, TRIA, nodeIds);
1533 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1536 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1537 ptr = getIds(ptr, EDGE, nodeIds);
1538 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1539 ptr = getIds(ptr, TRIA, nodeIds);
1540 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1542 case 3107: // ERR 3107 : 2 IN 4 1
1543 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1544 ptr = getIds(ptr, NODE, nodeIds);
1545 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1546 ptr = getIds(ptr, EDGE, nodeIds);
1547 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1549 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1550 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1551 ptr = getIds(ptr, EDGE, nodeIds);
1552 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1554 case 9000: // ERR 9000
1555 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1556 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1557 // A too small volume element is detected. Are reported the index of the element,
1558 // its four vertex indices, its volume and the tolerance threshold value
1559 ptr = getIds(ptr, ID, nodeIds);
1560 ptr = getIds(ptr, VOL, nodeIds);
1561 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1562 // even if all nodes found, volume it most probably invisible,
1563 // add its faces to demenstrate it anyhow
1565 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1566 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1567 faceNodes[2] = nodeIds[3]; // 013
1568 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1569 faceNodes[1] = nodeIds[2]; // 023
1570 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1571 faceNodes[0] = nodeIds[1]; // 123
1572 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1575 case 9001: // ERR 9001
1576 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1577 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1578 // %% NUMBER OF NULL VOLUME TETS : 0
1579 // There exists at least a null or negative volume element
1582 // There exist n null or negative volume elements
1585 // A too small volume element is detected
1588 // A too bad quality face is detected. This face is considered degenerated,
1589 // its index, its three vertex indices together with its quality value are reported
1590 break; // same as next
1591 case 9112: // ERR 9112
1592 // FACE 2 WITH VERTICES : 4 2 5
1593 // SMALL INRADIUS : 0.
1594 // A too bad quality face is detected. This face is degenerated,
1595 // its index, its three vertex indices together with its inradius are reported
1596 ptr = getIds(ptr, ID, nodeIds);
1597 ptr = getIds(ptr, TRIA, nodeIds);
1598 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1599 // add triangle edges as it most probably has zero area and hence invisible
1601 vector<int> edgeNodes(2);
1602 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1603 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1604 edgeNodes[1] = nodeIds[2]; // 0-2
1605 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1606 edgeNodes[0] = nodeIds[1]; // 1-2
1607 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1612 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1614 continue; // not to report same error several times
1616 // const SMDS_MeshElement* nullElem = 0;
1617 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1619 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1620 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1621 // if ( oneMoreErrorType )
1622 // continue; // not to report different types of errors with bad elements
1625 // store bad elements
1626 //if ( allElemsOk ) {
1627 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1628 for ( ; elem != badElems.end(); ++elem )
1629 addBadInputElement( *elem );
1633 string text = translateError( errNum );
1634 if ( errDescription.find( text ) == text.npos ) {
1635 if ( !errDescription.empty() )
1636 errDescription << "\n";
1637 errDescription << text;
1642 if ( errDescription.empty() ) { // no errors found
1643 char msg[] = "connection to server failed";
1644 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1645 errDescription << "Licence problems.";
1648 if ( errDescription.empty() )
1649 errDescription << "See " << logFile << " for problem description";
1651 errDescription << "\nSee " << logFile << " for more information";
1653 return error( errDescription );
1656 //================================================================================
1658 * \brief Creates _Ghs2smdsConvertor
1660 //================================================================================
1662 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1663 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1667 //================================================================================
1669 * \brief Creates _Ghs2smdsConvertor
1671 //================================================================================
1673 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1674 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1678 //================================================================================
1680 * \brief Return SMDS element by ids of GHS3D nodes
1682 //================================================================================
1684 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1686 size_t nbNodes = ghsNodes.size();
1687 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1688 for ( size_t i = 0; i < nbNodes; ++i ) {
1689 int ghsNode = ghsNodes[ i ];
1690 if ( _ghs2NodeMap ) {
1691 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1692 if ( in == _ghs2NodeMap->end() )
1694 nodes[ i ] = in->second;
1697 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1699 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1705 if ( nbNodes == 2 ) {
1706 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1708 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1711 if ( nbNodes == 3 ) {
1712 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1714 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1718 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );