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 node1 on geomFace
577 SMESH_MesherHelper helper( mesh );
578 gp_XY uv = helper.GetNodeUV( geomFace, node1 );
580 // check that uv is correct
581 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
582 double tol = BRep_Tool::Tolerance( geomFace );
583 BRepAdaptor_Surface surface( geomFace );
584 if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
585 // project node1 onto geomFace to get right UV
586 GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface() );
587 if ( !projector.IsDone() || projector.NbPoints() < 1 )
589 Quantity_Parameter U,V;
590 projector.LowerDistanceParameters(U,V);
593 // normale to face at node1
594 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
595 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
596 gp_Vec vec12( node1Pnt, node2Pnt );
597 gp_Vec vec13( node1Pnt, node3Pnt );
598 gp_Vec meshNormal = vec12 ^ vec13;
599 if ( meshNormal.SquareMagnitude() < DBL_MIN )
602 // normale to geomFace at UV
604 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
605 gp_Vec geomNormal = du ^ dv;
606 if ( geomNormal.SquareMagnitude() < DBL_MIN )
607 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
608 if ( geomFace.Orientation() == TopAbs_REVERSED )
609 geomNormal.Reverse();
612 bool isReverse = ( meshNormal * geomNormal ) < 0;
614 return meshDS->ShapeToIndex( solid1 );
616 if ( solids.Extent() == 1 )
617 return HOLE_ID; // we are inside a hole
619 return meshDS->ShapeToIndex( solids(2) );
622 //=======================================================================
623 //function : readResultFile
625 //=======================================================================
627 static bool readResultFile(const int fileOpen,
629 TopoDS_Shape tabShape[],
632 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
642 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
645 int nbElems, nbNodes, nbInputNodes;
646 int nodeId/*, triangleId*/;
648 int ID, shapeID, ghs3dShapeID;
651 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
653 int *tab, *tabID, *nodeID, *nodeAssigne;
655 const SMDS_MeshNode **node;
658 //tabID = new int[nbShape];
660 coord = new double[3];
661 node = new const SMDS_MeshNode*[4];
664 SMDS_MeshNode * aNewNode;
665 map <int,const SMDS_MeshNode*>::iterator itOnNode;
666 SMDS_MeshElement* aTet;
671 // Read the file state
672 fileStat = fstat(fileOpen, &status);
673 length = status.st_size;
675 // Mapping the result file into memory
676 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
679 ptr = readMapIntLine(ptr, tab);
684 nbInputNodes = tab[2];
686 nodeAssigne = new int[ nbNodes+1 ];
689 aSolid = tabShape[0];
691 // Reading the nodeId
692 for (int i=0; i < 4*nbElems; i++)
693 nodeId = strtol(ptr, &ptr, 10);
695 // Reading the nodeCoor and update the nodeMap
696 for (int iNode=1; iNode <= nbNodes; iNode++) {
697 for (int iCoor=0; iCoor < 3; iCoor++)
698 coord[ iCoor ] = strtod(ptr, &ptr);
699 nodeAssigne[ iNode ] = 1;
700 if ( iNode > nbInputNodes ) {
701 nodeAssigne[ iNode ] = 0;
702 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
703 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
707 // Reading the number of triangles which corresponds to the number of sub-domains
708 nbTriangle = strtol(ptr, &ptr, 10);
710 tabID = new int[nbTriangle];
711 for (int i=0; i < nbTriangle; i++) {
713 // find the solid corresponding to GHS3D sub-domain following
714 // the technique proposed in GHS3D manual in chapter
715 // "B.4 Subdomain (sub-region) assignment"
716 int nodeId1 = strtol(ptr, &ptr, 10);
717 int nodeId2 = strtol(ptr, &ptr, 10);
718 int nodeId3 = strtol(ptr, &ptr, 10);
719 if ( nbTriangle > 1 ) {
720 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
721 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
722 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
725 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
727 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
729 } catch ( Standard_Failure ) {
736 if ( nbTriangle <= nbShape ) // no holes
737 toMeshHoles = true; // not avoid creating tetras in holes
739 // Associating the tetrahedrons to the shapes
740 shapeID = compoundID;
741 for (int iElem = 0; iElem < nbElems; iElem++) {
742 for (int iNode = 0; iNode < 4; iNode++) {
743 ID = strtol(tetraPtr, &tetraPtr, 10);
744 itOnNode = theGhs3dIdToNodeMap.find(ID);
745 node[ iNode ] = itOnNode->second;
746 nodeID[ iNode ] = ID;
748 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
749 // tetras within holes depending on hypo option,
750 // so we first check if aTet is inside a hole and then create it
751 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
752 if ( nbTriangle > 1 ) {
753 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
754 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
755 if ( tabID[ ghs3dShapeID ] == 0 ) {
757 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
758 if ( toMeshHoles || state == TopAbs_IN )
759 shapeID = theMeshDS->ShapeToIndex( aSolid );
760 tabID[ ghs3dShapeID ] = shapeID;
763 shapeID = tabID[ ghs3dShapeID ];
765 else if ( nbShape > 1 ) {
766 // Case where nbTriangle == 1 while nbShape == 2 encountered
767 // with compound of 2 boxes and "To mesh holes"==False,
768 // so there are no subdomains specified for each tetrahedron.
769 // Try to guess a solid by a node already bound to shape
771 for ( int i=0; i<4 && shapeID==0; i++ ) {
772 if ( nodeAssigne[ nodeID[i] ] == 1 &&
773 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
774 node[i]->GetPosition()->GetShapeId() > 1 )
776 shapeID = node[i]->GetPosition()->GetShapeId();
780 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
781 shapeID = theMeshDS->ShapeToIndex( aSolid );
784 // set new nodes and tetrahedron onto the shape
785 for ( int i=0; i<4; i++ ) {
786 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
787 if ( shapeID != HOLE_ID )
788 theMeshDS->SetNodeInVolume( node[i], shapeID );
789 nodeAssigne[ nodeID[i] ] = shapeID;
792 if ( toMeshHoles || shapeID != HOLE_ID ) {
793 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
794 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
797 shapeIDs.insert( shapeID );
800 // Remove nodes of tetras inside holes if !toMeshHoles
801 if ( !toMeshHoles ) {
802 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
803 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
804 ID = itOnNode->first;
805 if ( nodeAssigne[ ID ] == HOLE_ID )
806 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
811 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
812 munmap(mapPtr, length);
820 delete [] nodeAssigne;
823 if ( shapeIDs.size() != nbShape ) {
824 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
825 for (int i=0; i<nbShape; i++) {
826 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
827 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
828 cout << " Solid #" << shapeID << " not found" << endl;
836 //=======================================================================
837 //function : readResultFile
839 //=======================================================================
841 static bool readResultFile(const int fileOpen,
842 SMESHDS_Mesh* theMeshDS,
844 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
854 int nbElems, nbNodes, nbInputNodes;
855 int nodeId, triangleId;
861 const SMDS_MeshNode **node;
864 coord = new double[3];
865 node = new const SMDS_MeshNode*[4];
867 SMDS_MeshNode * aNewNode;
868 map <int,const SMDS_MeshNode*>::iterator IdNode;
869 SMDS_MeshElement* aTet;
871 // Read the file state
872 fileStat = fstat(fileOpen, &status);
873 length = status.st_size;
875 // Mapping the result file into memory
876 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
879 ptr = readMapIntLine(ptr, tab);
884 nbInputNodes = tab[2];
886 theNodeByGhs3dId.resize( nbNodes );
888 // Reading the nodeId
889 for (int i=0; i < 4*nbElems; i++)
890 nodeId = strtol(ptr, &ptr, 10);
892 // Reading the nodeCoor and update the nodeMap
893 shapeID = theMeshDS->ShapeToIndex( aSolid );
894 for (int iNode=0; iNode < nbNodes; iNode++) {
895 for (int iCoor=0; iCoor < 3; iCoor++)
896 coord[ iCoor ] = strtod(ptr, &ptr);
897 if ((iNode+1) > nbInputNodes) {
898 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
899 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
900 theNodeByGhs3dId[ iNode ] = aNewNode;
904 // Reading the triangles
905 nbTriangle = strtol(ptr, &ptr, 10);
907 for (int i=0; i < 3*nbTriangle; i++)
908 triangleId = strtol(ptr, &ptr, 10);
912 // Associating the tetrahedrons to the shapes
913 for (int iElem = 0; iElem < nbElems; iElem++) {
914 for (int iNode = 0; iNode < 4; iNode++) {
915 ID = strtol(tetraPtr, &tetraPtr, 10);
916 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
918 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
919 shapeID = theMeshDS->ShapeToIndex( aSolid );
920 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
923 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
924 munmap(mapPtr, length);
934 //=============================================================================
936 *Here we are going to use the GHS3D mesher
938 //=============================================================================
940 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
941 const TopoDS_Shape& theShape)
944 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
946 _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
948 // create bounding box for every shape inside the compound
951 TopoDS_Shape* tabShape;
953 tabShape = new TopoDS_Shape[_nbShape];
954 tabBox = new double*[_nbShape];
955 for (int i=0; i<_nbShape; i++)
956 tabBox[i] = new double[6];
957 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
959 TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
960 for (; expBox.More(); expBox.Next()) {
961 tabShape[iShape] = expBox.Current();
963 BRepBndLib::Add(expBox.Current(), BoundingBox);
964 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
965 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
966 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
967 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
971 // a unique working file name
972 // to avoid access to the same files by eg different users
973 TCollection_AsciiString aGenericName
974 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
976 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
977 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
978 aFacesFileName = aGenericName + ".faces"; // in faces
979 aPointsFileName = aGenericName + ".points"; // in points
980 aResultFileName = aGenericName + ".noboite";// out points and volumes
981 aBadResFileName = aGenericName + ".boite"; // out bad result
982 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
983 aLogFileName = aGenericName + ".log"; // log
989 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
990 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
994 aFacesFile->is_open() && aPointsFile->is_open();
996 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
999 INFOS( "Can't write into " << aFacesFileName);
1000 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1002 map <int,int> aSmdsToGhs3dIdMap;
1003 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1005 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1006 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1009 aPointsFile.close();
1012 if ( !_keepFiles ) {
1013 OSD_File( aFacesFileName ).Remove();
1014 OSD_File( aPointsFileName ).Remove();
1016 return error(COMPERR_BAD_INPUT_MESH);
1018 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1020 // -----------------
1022 // -----------------
1024 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1025 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1026 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1029 cout << "Ghs3d execution..." << endl;
1030 cout << cmd << endl;
1032 system( cmd.ToCString() ); // run
1035 cout << "End of Ghs3d execution !" << endl;
1041 // Mapping the result file
1044 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1045 if ( fileOpen < 0 ) {
1047 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1048 cout << "Log: " << aLogFileName << endl;
1053 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1054 Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1058 // ---------------------
1059 // remove working files
1060 // ---------------------
1065 OSD_File( aLogFileName ).Remove();
1067 else if ( OSD_File( aLogFileName ).Size() > 0 )
1069 // get problem description from the log file
1070 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1071 storeErrorDescription( aLogFileName, conv );
1075 // the log file is empty
1076 OSD_File( aLogFileName ).Remove();
1077 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1078 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1081 if ( !_keepFiles ) {
1082 OSD_File( aFacesFileName ).Remove();
1083 OSD_File( aPointsFileName ).Remove();
1084 OSD_File( aResultFileName ).Remove();
1085 OSD_File( aBadResFileName ).Remove();
1086 OSD_File( aBbResFileName ).Remove();
1088 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1091 cout << "treated !" << endl;
1094 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1101 //=============================================================================
1103 *Here we are going to use the GHS3D mesher w/o geometry
1105 //=============================================================================
1106 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1107 SMESH_MesherHelper* aHelper)
1109 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1111 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1112 TopoDS_Shape theShape = aHelper->GetSubShape();
1114 // a unique working file name
1115 // to avoid access to the same files by eg different users
1116 TCollection_AsciiString aGenericName
1117 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1119 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1120 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1121 aFacesFileName = aGenericName + ".faces"; // in faces
1122 aPointsFileName = aGenericName + ".points"; // in points
1123 aResultFileName = aGenericName + ".noboite";// out points and volumes
1124 aBadResFileName = aGenericName + ".boite"; // out bad result
1125 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1126 aLogFileName = aGenericName + ".log"; // log
1128 // -----------------
1130 // -----------------
1132 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1133 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1136 aFacesFile->is_open() && aPointsFile->is_open();
1138 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1142 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1144 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1146 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1147 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1150 aPointsFile.close();
1153 if ( !_keepFiles ) {
1154 OSD_File( aFacesFileName ).Remove();
1155 OSD_File( aPointsFileName ).Remove();
1157 return error(COMPERR_BAD_INPUT_MESH);
1159 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1161 // -----------------
1163 // -----------------
1165 TCollection_AsciiString cmd =
1166 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1167 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1168 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1170 system( cmd.ToCString() ); // run
1176 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1177 if ( fileOpen < 0 ) {
1179 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1180 cout << "Log: " << aLogFileName << endl;
1185 Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1188 // ---------------------
1189 // remove working files
1190 // ---------------------
1195 OSD_File( aLogFileName ).Remove();
1197 else if ( OSD_File( aLogFileName ).Size() > 0 )
1199 // get problem description from the log file
1200 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1201 storeErrorDescription( aLogFileName, conv );
1204 // the log file is empty
1205 OSD_File( aLogFileName ).Remove();
1206 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1207 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1212 OSD_File( aFacesFileName ).Remove();
1213 OSD_File( aPointsFileName ).Remove();
1214 OSD_File( aResultFileName ).Remove();
1215 OSD_File( aBadResFileName ).Remove();
1216 OSD_File( aBbResFileName ).Remove();
1222 //================================================================================
1224 * \brief Provide human readable text by error code reported by ghs3d
1226 //================================================================================
1228 static string translateError(const int errNum)
1232 return "The surface mesh includes a face of type other than edge, "
1233 "triangle or quadrilateral. This face type is not supported.";
1235 return "Not enough memory for the face table.";
1237 return "Not enough memory.";
1239 return "Not enough memory.";
1241 return "Face is ignored.";
1243 return "End of file. Some data are missing in the file.";
1245 return "Read error on the file. There are wrong data in the file.";
1247 return "the metric file is inadequate (dimension other than 3).";
1249 return "the metric file is inadequate (values not per vertices).";
1251 return "the metric file contains more than one field.";
1253 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1254 "value of number of mesh vertices in the \".noboite\" file.";
1256 return "Too many sub-domains.";
1258 return "the number of vertices is negative or null.";
1260 return "the number of faces is negative or null.";
1262 return "A face has a null vertex.";
1264 return "incompatible data.";
1266 return "the number of vertices is negative or null.";
1268 return "the number of vertices is negative or null (in the \".mesh\" file).";
1270 return "the number of faces is negative or null.";
1272 return "A face appears more than once in the input surface mesh.";
1274 return "An edge appears more than once in the input surface mesh.";
1276 return "A face has a vertex negative or null.";
1278 return "NOT ENOUGH MEMORY.";
1280 return "Not enough available memory.";
1282 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1283 "in terms of quality or the input list of points is wrong.";
1285 return "Some vertices are too close to one another or coincident.";
1287 return "Some vertices are too close to one another or coincident.";
1289 return "A vertex cannot be inserted.";
1291 return "There are at least two points considered as coincident.";
1293 return "Some vertices are too close to one another or coincident.";
1295 return "The surface mesh regeneration step has failed.";
1297 return "Constrained edge cannot be enforced.";
1299 return "Constrained face cannot be enforced.";
1301 return "Missing faces.";
1303 return "No guess to start the definition of the connected component(s).";
1305 return "The surface mesh includes at least one hole. The domain is not well defined.";
1307 return "Impossible to define a component.";
1309 return "The surface edge intersects another surface edge.";
1311 return "The surface edge intersects the surface face.";
1313 return "One boundary point lies within a surface face.";
1315 return "One surface edge intersects a surface face.";
1317 return "One boundary point lies within a surface edge.";
1319 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1320 "to too many swaps.";
1322 return "Edge is unique (i.e., bounds a hole in the surface).";
1324 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1326 return "Too many components, too many sub-domain.";
1328 return "The surface mesh includes at least one hole. "
1329 "Therefore there is no domain properly defined.";
1331 return "Statistics.";
1333 return "Statistics.";
1335 return "Warning, it is dramatically tedious to enforce the boundary items.";
1337 return "Not enough memory at this time, nevertheless, the program continues. "
1338 "The expected mesh will be correct but not really as large as required.";
1340 return "see above error code, resulting quality may be poor.";
1342 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1344 return "Unknown face type.";
1347 return "End of file. Some data are missing in the file.";
1349 return "A too small volume element is detected.";
1351 return "There exists at least a null or negative volume element.";
1353 return "There exist null or negative volume elements.";
1355 return "A too small volume element is detected. A face is considered being degenerated.";
1357 return "Some element is suspected to be very bad shaped or wrong.";
1359 return "A too bad quality face is detected. This face is considered degenerated.";
1361 return "A too bad quality face is detected. This face is degenerated.";
1363 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1365 return "Abnormal error occured, contact hotline.";
1367 return "Not enough memory for the face table.";
1369 return "The algorithm cannot run further. "
1370 "The surface mesh is probably very bad in terms of quality.";
1372 return "Bad vertex number.";
1377 //================================================================================
1379 * \brief Retrieve from a string given number of integers
1381 //================================================================================
1383 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1386 ids.reserve( nbIds );
1389 while ( !isdigit( *ptr )) ++ptr;
1390 if ( ptr[-1] == '-' ) --ptr;
1391 ids.push_back( strtol( ptr, &ptr, 10 ));
1397 //================================================================================
1399 * \brief Retrieve problem description form a log file
1400 * \retval bool - always false
1402 //================================================================================
1404 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1405 const _Ghs2smdsConvertor & toSmdsConvertor )
1409 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1411 int file = ::open (logFile.ToCString(), O_RDONLY);
1414 return error( SMESH_Comment("See ") << logFile << " for problem description");
1417 // struct stat status;
1418 // fstat(file, &status);
1419 // size_t length = status.st_size;
1420 off_t length = lseek( file, 0, SEEK_END);
1421 lseek( file, 0, SEEK_SET);
1424 vector< char > buf( length );
1425 int nBytesRead = ::read (file, & buf[0], length);
1427 char* ptr = & buf[0];
1428 char* bufEnd = ptr + nBytesRead;
1430 SMESH_Comment errDescription;
1432 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1434 // look for errors "ERR #"
1436 set<string> foundErrorStr; // to avoid reporting same error several times
1437 set<int> elemErrorNums; // not to report different types of errors with bad elements
1438 while ( ++ptr < bufEnd )
1440 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1443 list<const SMDS_MeshElement*> badElems;
1444 vector<int> nodeIds;
1448 int errNum = strtol(ptr, &ptr, 10);
1449 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1451 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1452 ptr = getIds(ptr, NODE, nodeIds);
1453 ptr = getIds(ptr, TRIA, nodeIds);
1454 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1456 case 1000: // ERR 1000 : 1 3 2
1457 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1458 ptr = getIds(ptr, TRIA, nodeIds);
1459 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1462 // Edge (e1, e2) appears more than once in the input surface mesh
1463 ptr = getIds(ptr, EDGE, nodeIds);
1464 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1467 // Face (f 1, f 2, f 3) has a vertex negative or null
1468 ptr = getIds(ptr, TRIA, nodeIds);
1469 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1472 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1473 ptr = getIds(ptr, NODE, nodeIds);
1474 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1475 ptr = getIds(ptr, NODE, nodeIds);
1476 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1479 // Vertex v1 cannot be inserted (warning).
1480 ptr = getIds(ptr, NODE, nodeIds);
1481 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1484 // There are at least two points whose distance is dist, i.e., considered as coincident
1485 case 2103: // ERR 2103 : 16 WITH 3
1486 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1487 ptr = getIds(ptr, NODE, nodeIds);
1488 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1489 ptr = getIds(ptr, NODE, nodeIds);
1490 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1493 // Constrained edge (e1, e2) cannot be enforced (warning).
1494 ptr = getIds(ptr, EDGE, nodeIds);
1495 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1498 // Constrained face (f 1, f 2, f 3) cannot be enforced
1499 ptr = getIds(ptr, TRIA, nodeIds);
1500 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1502 case 3103: // ERR 3103 : 1 2 WITH 7 3
1503 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1504 ptr = getIds(ptr, EDGE, nodeIds);
1505 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1506 ptr = getIds(ptr, EDGE, nodeIds);
1507 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1509 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1510 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1511 ptr = getIds(ptr, EDGE, nodeIds);
1512 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1513 ptr = getIds(ptr, TRIA, nodeIds);
1514 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1516 case 3105: // ERR 3105 : 8 IN 2 3 5
1517 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1518 ptr = getIds(ptr, NODE, nodeIds);
1519 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1520 ptr = getIds(ptr, TRIA, nodeIds);
1521 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1524 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1525 ptr = getIds(ptr, EDGE, nodeIds);
1526 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1527 ptr = getIds(ptr, TRIA, nodeIds);
1528 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1530 case 3107: // ERR 3107 : 2 IN 4 1
1531 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1532 ptr = getIds(ptr, NODE, nodeIds);
1533 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1534 ptr = getIds(ptr, EDGE, nodeIds);
1535 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1537 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1538 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1539 ptr = getIds(ptr, EDGE, nodeIds);
1540 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1542 case 9000: // ERR 9000
1543 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1544 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1545 // A too small volume element is detected. Are reported the index of the element,
1546 // its four vertex indices, its volume and the tolerance threshold value
1547 ptr = getIds(ptr, ID, nodeIds);
1548 ptr = getIds(ptr, VOL, nodeIds);
1549 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1550 // even if all nodes found, volume it most probably invisible,
1551 // add its faces to demenstrate it anyhow
1553 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1554 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1555 faceNodes[2] = nodeIds[3]; // 013
1556 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1557 faceNodes[1] = nodeIds[2]; // 023
1558 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1559 faceNodes[0] = nodeIds[1]; // 123
1560 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1563 case 9001: // ERR 9001
1564 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1565 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1566 // %% NUMBER OF NULL VOLUME TETS : 0
1567 // There exists at least a null or negative volume element
1570 // There exist n null or negative volume elements
1573 // A too small volume element is detected
1576 // A too bad quality face is detected. This face is considered degenerated,
1577 // its index, its three vertex indices together with its quality value are reported
1578 break; // same as next
1579 case 9112: // ERR 9112
1580 // FACE 2 WITH VERTICES : 4 2 5
1581 // SMALL INRADIUS : 0.
1582 // A too bad quality face is detected. This face is degenerated,
1583 // its index, its three vertex indices together with its inradius are reported
1584 ptr = getIds(ptr, ID, nodeIds);
1585 ptr = getIds(ptr, TRIA, nodeIds);
1586 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1587 // add triangle edges as it most probably has zero area and hence invisible
1589 vector<int> edgeNodes(2);
1590 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1591 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1592 edgeNodes[1] = nodeIds[2]; // 0-2
1593 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1594 edgeNodes[0] = nodeIds[1]; // 1-2
1595 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1600 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1602 continue; // not to report same error several times
1604 // const SMDS_MeshElement* nullElem = 0;
1605 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1607 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1608 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1609 // if ( oneMoreErrorType )
1610 // continue; // not to report different types of errors with bad elements
1613 // store bad elements
1614 //if ( allElemsOk ) {
1615 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1616 for ( ; elem != badElems.end(); ++elem )
1617 addBadInputElement( *elem );
1621 string text = translateError( errNum );
1622 if ( errDescription.find( text ) == text.npos ) {
1623 if ( !errDescription.empty() )
1624 errDescription << "\n";
1625 errDescription << text;
1630 if ( errDescription.empty() ) { // no errors found
1631 char msg[] = "connection to server failed";
1632 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1633 errDescription << "Licence problems.";
1636 if ( errDescription.empty() )
1637 errDescription << "See " << logFile << " for problem description";
1639 errDescription << "\nSee " << logFile << " for more information";
1641 return error( errDescription );
1644 //================================================================================
1646 * \brief Creates _Ghs2smdsConvertor
1648 //================================================================================
1650 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1651 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1655 //================================================================================
1657 * \brief Creates _Ghs2smdsConvertor
1659 //================================================================================
1661 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1662 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1666 //================================================================================
1668 * \brief Return SMDS element by ids of GHS3D nodes
1670 //================================================================================
1672 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1674 size_t nbNodes = ghsNodes.size();
1675 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1676 for ( size_t i = 0; i < nbNodes; ++i ) {
1677 int ghsNode = ghsNodes[ i ];
1678 if ( _ghs2NodeMap ) {
1679 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1680 if ( in == _ghs2NodeMap->end() )
1682 nodes[ i ] = in->second;
1685 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1687 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1693 if ( nbNodes == 2 ) {
1694 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1696 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1699 if ( nbNodes == 3 ) {
1700 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1702 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1706 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );