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 <TopTools_MapOfShape.hxx>
62 #include <BRepGProp.hxx>
63 #include <GProp_GProps.hxx>
65 #include "utilities.h"
68 #include <sys/sysinfo.h>
71 //#include <Standard_Stream.hxx>
74 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
95 //=============================================================================
99 //=============================================================================
101 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
102 : SMESH_3D_Algo(hypId, studyId, gen)
104 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
106 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
107 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
110 _compatibleHypothesis.push_back("GHS3D_Parameters");
111 _requireShape = false; // can work without shape
114 //=============================================================================
118 //=============================================================================
120 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
122 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
125 //=============================================================================
129 //=============================================================================
131 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
132 const TopoDS_Shape& aShape,
133 Hypothesis_Status& aStatus )
135 aStatus = SMESH_Hypothesis::HYP_OK;
137 // there is only one compatible Hypothesis so far
140 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
142 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
144 _keepFiles = _hyp->GetKeepFiles();
149 //=======================================================================
150 //function : findShape
152 //=======================================================================
154 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
156 const TopoDS_Shape shape[],
159 TopAbs_State * state = 0)
162 int j, iShape, nbNode = 4;
164 for ( j=0; j<nbNode; j++ ) {
165 gp_XYZ p ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
166 if ( aNode[j]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE ) {
173 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
174 if (state) *state = SC.State();
175 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
176 for (iShape = 0; iShape < nShape; iShape++) {
177 aShape = shape[iShape];
178 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
179 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
180 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
181 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
182 if (state) *state = SC.State();
183 if (SC.State() == TopAbs_IN)
191 //=======================================================================
192 //function : readMapIntLine
194 //=======================================================================
196 static char* readMapIntLine(char* ptr, int tab[]) {
200 for ( int i=0; i<17; i++ ) {
201 intVal = strtol(ptr, &ptr, 10);
208 //=======================================================================
209 //function : countShape
211 //=======================================================================
213 template < class Mesh, class Shape >
214 static int countShape( Mesh* mesh, Shape shape ) {
215 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
217 for ( ; expShape.More(); expShape.Next() ) {
223 //=======================================================================
224 //function : writeFaces
226 //=======================================================================
228 static bool writeFaces (ofstream & theFile,
229 SMESHDS_Mesh * theMesh,
230 const map <int,int> & theSmdsToGhs3dIdMap)
234 // NB_ELEMS DUMMY_INT
235 // Loop from 1 to NB_ELEMS
236 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
238 int nbShape = countShape( theMesh, TopAbs_FACE );
240 int *tabID; tabID = new int[nbShape];
241 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
243 SMESHDS_SubMesh* theSubMesh;
244 const SMDS_MeshElement* aFace;
245 const char* space = " ";
246 const int dummyint = 0;
247 map<int,int>::const_iterator itOnMap;
248 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
249 int shapeID, nbNodes, aSmdsID;
252 cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
255 theFile << space << theMesh->NbFaces() << space << dummyint << endl;
257 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
258 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
260 aShape = expface.Current();
261 shapeID = theMesh->ShapeToIndex( aShape );
263 for ( int j=0; j<=i; j++) {
264 if ( shapeID == tabID[j] ) {
271 tabShape[i] = aShape;
274 for ( int i =0; i < nbShape; i++ ) {
275 if ( tabID[i] != 0 ) {
276 aShape = tabShape[i];
278 theSubMesh = theMesh->MeshElements( aShape );
279 if ( !theSubMesh ) continue;
280 itOnSubMesh = theSubMesh->GetElements();
281 while ( itOnSubMesh->more() ) {
282 aFace = itOnSubMesh->next();
283 nbNodes = aFace->NbNodes();
285 theFile << space << nbNodes;
287 itOnSubFace = aFace->nodesIterator();
288 while ( itOnSubFace->more() ) {
290 aSmdsID = itOnSubFace->next()->GetID();
291 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
292 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
294 theFile << space << (*itOnMap).second;
297 // (NB_NODES + 1) times: DUMMY_INT
298 for ( int j=0; j<=nbNodes; j++)
299 theFile << space << dummyint;
312 //=======================================================================
313 //function : writeFaces
314 //purpose : Write Faces in case if generate 3D mesh w/o geometry
315 //=======================================================================
317 static bool writeFaces (ofstream & theFile,
318 SMESHDS_Mesh * theMesh,
319 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
323 // NB_ELEMS DUMMY_INT
324 // Loop from 1 to NB_ELEMS
325 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
329 list< const SMDS_MeshElement* > faces;
330 list< const SMDS_MeshElement* >::iterator f;
331 map< const SMDS_MeshNode*,int >::iterator it;
332 SMDS_ElemIteratorPtr nodeIt;
333 const SMDS_MeshElement* elem;
336 const char* space = " ";
337 const int dummyint = 0;
339 //get all faces from mesh
340 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
341 while ( eIt->more() ) {
342 const SMDS_MeshElement* elem = eIt->next();
345 faces.push_back( elem );
352 cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
354 // NB_ELEMS DUMMY_INT
355 theFile << space << nbFaces << space << dummyint << endl;
357 // Loop from 1 to NB_ELEMS
359 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
361 for ( ; f != faces.end(); ++f )
365 nbNodes = elem->NbNodes();
366 theFile << space << nbNodes;
368 // NODE_NB_1 NODE_NB_2 ...
369 nodeIt = elem->nodesIterator();
370 while ( nodeIt->more() )
373 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
374 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
375 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
376 theFile << space << it->second;
379 // (NB_NODES + 1) times: DUMMY_INT
380 for ( int i=0; i<=nbNodes; i++)
381 theFile << space << dummyint;
386 // put nodes to theNodeByGhs3dId vector
387 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
388 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
389 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
391 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
397 //=======================================================================
398 //function : writePoints
400 //=======================================================================
402 static bool writePoints (ofstream & theFile,
403 SMESHDS_Mesh * theMesh,
404 map <int,int> & theSmdsToGhs3dIdMap,
405 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
410 // Loop from 1 to NB_NODES
413 int nbNodes = theMesh->NbNodes();
417 const char* space = " ";
418 const int dummyint = 0;
421 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
422 const SMDS_MeshNode* node;
425 theFile << space << nbNodes << endl;
427 cout << "The initial 2D mesh contains :" << endl;
428 cout << " " << nbNodes << " vertices" << endl;
430 // Loop from 1 to NB_NODES
435 theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
436 theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
441 << space << node->X()
442 << space << node->Y()
443 << space << node->Z()
444 << space << dummyint;
452 //=======================================================================
453 //function : writePoints
455 //=======================================================================
457 static bool writePoints (ofstream & theFile,
458 SMESHDS_Mesh * theMesh,
459 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
464 // Loop from 1 to NB_NODES
467 //int nbNodes = theMesh->NbNodes();
468 int nbNodes = theNodeByGhs3dId.size();
472 const char* space = " ";
473 const int dummyint = 0;
475 const SMDS_MeshNode* node;
478 theFile << space << nbNodes << endl;
479 cout << nbNodes << " nodes" << endl;
481 // Loop from 1 to NB_NODES
483 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
484 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
485 for ( ; nodeIt != after; ++nodeIt )
491 << space << node->X()
492 << space << node->Y()
493 << space << node->Z()
494 << space << dummyint;
502 //=======================================================================
503 //function : findShapeID
504 //purpose : find the solid corresponding to GHS3D sub-domain following
505 // the technique proposed in GHS3D manual in chapter
506 // "B.4 Subdomain (sub-region) assignment"
507 //=======================================================================
509 static int findShapeID(SMESH_Mesh& mesh,
510 const SMDS_MeshNode* node1,
511 const SMDS_MeshNode* node2,
512 const SMDS_MeshNode* node3,
513 const bool toMeshHoles)
515 const int invalidID = 0;
516 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
518 // face the nodes belong to
519 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
523 // geom face the face assigned to
524 SMESH_MeshEditor editor(&mesh);
525 int geomFaceID = editor.FindShape( face );
528 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
529 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
531 TopoDS_Face geomFace = TopoDS::Face( shape );
533 // solids bounded by geom face
534 TopTools_IndexedMapOfShape solids, shells;
535 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
536 for ( ; ansIt.More(); ansIt.Next() ) {
537 switch ( ansIt.Value().ShapeType() ) {
539 solids.Add( ansIt.Value() ); break;
541 shells.Add( ansIt.Value() ); break;
545 // analyse found solids
546 if ( solids.Extent() == 0 || shells.Extent() == 0)
549 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
550 if ( solids.Extent() == 1 )
553 return meshDS->ShapeToIndex( solid1 );
555 // - Are we at a hole boundary face?
556 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
557 { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
559 TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
560 // check if any edge of shells(1) belongs to another shell
561 for ( ; eExp.More() && !touch; eExp.Next() ) {
562 ansIt = mesh.GetAncestors( eExp.Current() );
563 for ( ; ansIt.More() && !touch; ansIt.Next() ) {
564 if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
565 touch = ( !ansIt.Value().IsSame( shells(1) ));
569 return meshDS->ShapeToIndex( solid1 );
572 // find orientation of geom face within the first solid
573 TopExp_Explorer fExp( solid1, TopAbs_FACE );
574 for ( ; fExp.More(); fExp.Next() )
575 if ( geomFace.IsSame( fExp.Current() )) {
576 geomFace = TopoDS::Face( fExp.Current() );
580 return invalidID; // face not found
582 // normale to triangle
583 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
584 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
585 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
586 gp_Vec vec12( node1Pnt, node2Pnt );
587 gp_Vec vec13( node1Pnt, node3Pnt );
588 gp_Vec meshNormal = vec12 ^ vec13;
589 if ( meshNormal.SquareMagnitude() < DBL_MIN )
592 // find UV of node1 on geomFace
593 SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
594 const SMDS_MeshNode* nNotOnSeamEdge = 0;
595 if ( helper.IsSeamShape( node1->GetPosition()->GetShapeId() ))
596 if ( helper.IsSeamShape( node2->GetPosition()->GetShapeId() ))
597 nNotOnSeamEdge = node3;
599 nNotOnSeamEdge = node2;
601 gp_XY uv = helper.GetNodeUV( geomFace, node1, nNotOnSeamEdge, &checkUV );
602 // check that uv is correct
603 BRepAdaptor_Surface surface( geomFace );
604 double tol = BRep_Tool::Tolerance( geomFace );
605 if ( checkUV && node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
606 // GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface(), 2 * tol );
607 // if ( !projector.IsDone() || projector.NbPoints() < 1 || projector.LowerDistance() > 2 * tol)
609 // Quantity_Parameter U,V;
610 // projector.LowerDistanceParameters(U,V);
611 // if ( node1Pnt.Distance( surface.Value( U, V )) > 2 * tol )
613 // uv.SetCoord( U,V );
615 // normale to geomFace at UV
617 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
618 gp_Vec geomNormal = du ^ dv;
619 if ( geomNormal.SquareMagnitude() < DBL_MIN )
620 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
621 if ( geomFace.Orientation() == TopAbs_REVERSED )
622 geomNormal.Reverse();
625 bool isReverse = ( meshNormal * geomNormal ) < 0;
627 return meshDS->ShapeToIndex( solid1 );
629 if ( solids.Extent() == 1 )
630 return HOLE_ID; // we are inside a hole
632 return meshDS->ShapeToIndex( solids(2) );
635 //=======================================================================
636 //function : readResultFile
638 //=======================================================================
640 static bool readResultFile(const int fileOpen,
642 TopoDS_Shape tabShape[],
645 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
655 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
658 int nbElems, nbNodes, nbInputNodes;
659 int nodeId/*, triangleId*/;
661 int ID, shapeID, ghs3dShapeID;
664 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
666 int *tab, *tabID, *nodeID, *nodeAssigne;
668 const SMDS_MeshNode **node;
671 //tabID = new int[nbShape];
673 coord = new double[3];
674 node = new const SMDS_MeshNode*[4];
677 SMDS_MeshNode * aNewNode;
678 map <int,const SMDS_MeshNode*>::iterator itOnNode;
679 SMDS_MeshElement* aTet;
684 // Read the file state
685 fileStat = fstat(fileOpen, &status);
686 length = status.st_size;
688 // Mapping the result file into memory
689 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
692 ptr = readMapIntLine(ptr, tab);
697 nbInputNodes = tab[2];
699 nodeAssigne = new int[ nbNodes+1 ];
702 aSolid = tabShape[0];
704 // Reading the nodeId
705 for (int i=0; i < 4*nbElems; i++)
706 nodeId = strtol(ptr, &ptr, 10);
708 // Reading the nodeCoor and update the nodeMap
709 for (int iNode=1; iNode <= nbNodes; iNode++) {
710 for (int iCoor=0; iCoor < 3; iCoor++)
711 coord[ iCoor ] = strtod(ptr, &ptr);
712 nodeAssigne[ iNode ] = 1;
713 if ( iNode > nbInputNodes ) {
714 nodeAssigne[ iNode ] = 0;
715 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
716 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
720 // Reading the number of triangles which corresponds to the number of sub-domains
721 nbTriangle = strtol(ptr, &ptr, 10);
723 tabID = new int[nbTriangle];
724 for (int i=0; i < nbTriangle; i++) {
726 // find the solid corresponding to GHS3D sub-domain following
727 // the technique proposed in GHS3D manual in chapter
728 // "B.4 Subdomain (sub-region) assignment"
729 int nodeId1 = strtol(ptr, &ptr, 10);
730 int nodeId2 = strtol(ptr, &ptr, 10);
731 int nodeId3 = strtol(ptr, &ptr, 10);
732 if ( nbTriangle > 1 ) {
733 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
734 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
735 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
738 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
739 // -- 0020330: Pb with ghs3d as a submesh
740 // check that found shape is to be meshed
741 if ( tabID[i] > 0 ) {
742 const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( tabID[i] );
743 bool isToBeMeshed = false;
744 for ( int iS = 0; !isToBeMeshed && iS < nbShape; ++iS )
745 isToBeMeshed = foundShape.IsSame( tabShape[ iS ]);
749 // END -- 0020330: Pb with ghs3d as a submesh
751 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
753 } catch ( Standard_Failure ) {
760 if ( nbTriangle <= nbShape ) // no holes
761 toMeshHoles = true; // not avoid creating tetras in holes
763 // Associating the tetrahedrons to the shapes
764 shapeID = compoundID;
765 for (int iElem = 0; iElem < nbElems; iElem++) {
766 for (int iNode = 0; iNode < 4; iNode++) {
767 ID = strtol(tetraPtr, &tetraPtr, 10);
768 itOnNode = theGhs3dIdToNodeMap.find(ID);
769 node[ iNode ] = itOnNode->second;
770 nodeID[ iNode ] = ID;
772 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
773 // tetras within holes depending on hypo option,
774 // so we first check if aTet is inside a hole and then create it
775 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
776 if ( nbTriangle > 1 ) {
777 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
778 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
779 if ( tabID[ ghs3dShapeID ] == 0 ) {
781 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
782 if ( toMeshHoles || state == TopAbs_IN )
783 shapeID = theMeshDS->ShapeToIndex( aSolid );
784 tabID[ ghs3dShapeID ] = shapeID;
787 shapeID = tabID[ ghs3dShapeID ];
789 else if ( nbShape > 1 ) {
790 // Case where nbTriangle == 1 while nbShape == 2 encountered
791 // with compound of 2 boxes and "To mesh holes"==False,
792 // so there are no subdomains specified for each tetrahedron.
793 // Try to guess a solid by a node already bound to shape
795 for ( int i=0; i<4 && shapeID==0; i++ ) {
796 if ( nodeAssigne[ nodeID[i] ] == 1 &&
797 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
798 node[i]->GetPosition()->GetShapeId() > 1 )
800 shapeID = node[i]->GetPosition()->GetShapeId();
804 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
805 shapeID = theMeshDS->ShapeToIndex( aSolid );
808 // set new nodes and tetrahedron onto the shape
809 for ( int i=0; i<4; i++ ) {
810 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
811 if ( shapeID != HOLE_ID )
812 theMeshDS->SetNodeInVolume( node[i], shapeID );
813 nodeAssigne[ nodeID[i] ] = shapeID;
816 if ( toMeshHoles || shapeID != HOLE_ID ) {
817 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
818 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
821 shapeIDs.insert( shapeID );
824 // Remove nodes of tetras inside holes if !toMeshHoles
825 if ( !toMeshHoles ) {
826 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
827 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
828 ID = itOnNode->first;
829 if ( nodeAssigne[ ID ] == HOLE_ID )
830 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
835 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
836 munmap(mapPtr, length);
844 delete [] nodeAssigne;
847 if ( shapeIDs.size() != nbShape ) {
848 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
849 for (int i=0; i<nbShape; i++) {
850 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
851 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
852 cout << " Solid #" << shapeID << " not found" << endl;
860 //=======================================================================
861 //function : readResultFile
863 //=======================================================================
865 static bool readResultFile(const int fileOpen,
866 SMESHDS_Mesh* theMeshDS,
868 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
878 int nbElems, nbNodes, nbInputNodes;
879 int nodeId, triangleId;
885 const SMDS_MeshNode **node;
888 coord = new double[3];
889 node = new const SMDS_MeshNode*[4];
891 SMDS_MeshNode * aNewNode;
892 map <int,const SMDS_MeshNode*>::iterator IdNode;
893 SMDS_MeshElement* aTet;
895 // Read the file state
896 fileStat = fstat(fileOpen, &status);
897 length = status.st_size;
899 // Mapping the result file into memory
900 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
903 ptr = readMapIntLine(ptr, tab);
908 nbInputNodes = tab[2];
910 theNodeByGhs3dId.resize( nbNodes );
912 // Reading the nodeId
913 for (int i=0; i < 4*nbElems; i++)
914 nodeId = strtol(ptr, &ptr, 10);
916 // Reading the nodeCoor and update the nodeMap
917 shapeID = theMeshDS->ShapeToIndex( aSolid );
918 for (int iNode=0; iNode < nbNodes; iNode++) {
919 for (int iCoor=0; iCoor < 3; iCoor++)
920 coord[ iCoor ] = strtod(ptr, &ptr);
921 if ((iNode+1) > nbInputNodes) {
922 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
923 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
924 theNodeByGhs3dId[ iNode ] = aNewNode;
928 // Reading the triangles
929 nbTriangle = strtol(ptr, &ptr, 10);
931 for (int i=0; i < 3*nbTriangle; i++)
932 triangleId = strtol(ptr, &ptr, 10);
936 // Associating the tetrahedrons to the shapes
937 for (int iElem = 0; iElem < nbElems; iElem++) {
938 for (int iNode = 0; iNode < 4; iNode++) {
939 ID = strtol(tetraPtr, &tetraPtr, 10);
940 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
942 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
943 shapeID = theMeshDS->ShapeToIndex( aSolid );
944 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
947 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
948 munmap(mapPtr, length);
958 //=============================================================================
960 *Here we are going to use the GHS3D mesher
962 //=============================================================================
964 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
965 const TopoDS_Shape& theShape)
968 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
970 // we count the number of shapes
971 // _nbShape = countShape( meshDS, TopAbs_SOLID ); -- 0020330: Pb with ghs3d as a submesh
973 TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
974 for ( ; expBox.More(); expBox.Next() )
977 // create bounding box for every shape inside the compound
980 TopoDS_Shape* tabShape;
982 tabShape = new TopoDS_Shape[_nbShape];
983 tabBox = new double*[_nbShape];
984 for (int i=0; i<_nbShape; i++)
985 tabBox[i] = new double[6];
986 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
988 // TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID); -- 0020330:...ghs3d as a submesh
989 for (expBox.ReInit(); expBox.More(); expBox.Next()) {
990 tabShape[iShape] = expBox.Current();
992 BRepBndLib::Add(expBox.Current(), BoundingBox);
993 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
994 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
995 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
996 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1000 // a unique working file name
1001 // to avoid access to the same files by eg different users
1002 TCollection_AsciiString aGenericName
1003 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1005 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1006 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1007 aFacesFileName = aGenericName + ".faces"; // in faces
1008 aPointsFileName = aGenericName + ".points"; // in points
1009 aResultFileName = aGenericName + ".noboite";// out points and volumes
1010 aBadResFileName = aGenericName + ".boite"; // out bad result
1011 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1012 aLogFileName = aGenericName + ".log"; // log
1014 // -----------------
1016 // -----------------
1018 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1019 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1023 aFacesFile->is_open() && aPointsFile->is_open();
1025 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1028 INFOS( "Can't write into " << aFacesFileName);
1029 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1031 map <int,int> aSmdsToGhs3dIdMap;
1032 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1034 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1035 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1038 aPointsFile.close();
1041 if ( !_keepFiles ) {
1042 OSD_File( aFacesFileName ).Remove();
1043 OSD_File( aPointsFileName ).Remove();
1045 return error(COMPERR_BAD_INPUT_MESH);
1047 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1049 // -----------------
1051 // -----------------
1053 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1054 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1055 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1058 cout << "Ghs3d execution..." << endl;
1059 cout << cmd << endl;
1061 system( cmd.ToCString() ); // run
1064 cout << "End of Ghs3d execution !" << endl;
1070 // Mapping the result file
1073 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1074 if ( fileOpen < 0 ) {
1076 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1077 cout << "Log: " << aLogFileName << endl;
1082 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1083 Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1087 // ---------------------
1088 // remove working files
1089 // ---------------------
1094 OSD_File( aLogFileName ).Remove();
1096 else if ( OSD_File( aLogFileName ).Size() > 0 )
1098 // get problem description from the log file
1099 _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1100 storeErrorDescription( aLogFileName, conv );
1104 // the log file is empty
1105 OSD_File( aLogFileName ).Remove();
1106 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1107 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1110 if ( !_keepFiles ) {
1111 OSD_File( aFacesFileName ).Remove();
1112 OSD_File( aPointsFileName ).Remove();
1113 OSD_File( aResultFileName ).Remove();
1114 OSD_File( aBadResFileName ).Remove();
1115 OSD_File( aBbResFileName ).Remove();
1117 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1120 cout << "treated !" << endl;
1123 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1130 //=============================================================================
1132 *Here we are going to use the GHS3D mesher w/o geometry
1134 //=============================================================================
1135 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1136 SMESH_MesherHelper* aHelper)
1138 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1140 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1141 TopoDS_Shape theShape = aHelper->GetSubShape();
1143 // a unique working file name
1144 // to avoid access to the same files by eg different users
1145 TCollection_AsciiString aGenericName
1146 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1148 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1149 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1150 aFacesFileName = aGenericName + ".faces"; // in faces
1151 aPointsFileName = aGenericName + ".points"; // in points
1152 aResultFileName = aGenericName + ".noboite";// out points and volumes
1153 aBadResFileName = aGenericName + ".boite"; // out bad result
1154 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1155 aLogFileName = aGenericName + ".log"; // log
1157 // -----------------
1159 // -----------------
1161 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1162 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1165 aFacesFile->is_open() && aPointsFile->is_open();
1167 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1171 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1173 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1175 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1176 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1179 aPointsFile.close();
1182 if ( !_keepFiles ) {
1183 OSD_File( aFacesFileName ).Remove();
1184 OSD_File( aPointsFileName ).Remove();
1186 return error(COMPERR_BAD_INPUT_MESH);
1188 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1190 // -----------------
1192 // -----------------
1194 TCollection_AsciiString cmd =
1195 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1196 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1197 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1199 system( cmd.ToCString() ); // run
1205 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1206 if ( fileOpen < 0 ) {
1208 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1209 cout << "Log: " << aLogFileName << endl;
1214 Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1217 // ---------------------
1218 // remove working files
1219 // ---------------------
1224 OSD_File( aLogFileName ).Remove();
1226 else if ( OSD_File( aLogFileName ).Size() > 0 )
1228 // get problem description from the log file
1229 _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1230 storeErrorDescription( aLogFileName, conv );
1233 // the log file is empty
1234 OSD_File( aLogFileName ).Remove();
1235 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1236 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1241 OSD_File( aFacesFileName ).Remove();
1242 OSD_File( aPointsFileName ).Remove();
1243 OSD_File( aResultFileName ).Remove();
1244 OSD_File( aBadResFileName ).Remove();
1245 OSD_File( aBbResFileName ).Remove();
1251 //================================================================================
1253 * \brief Provide human readable text by error code reported by ghs3d
1255 //================================================================================
1257 static string translateError(const int errNum)
1261 return "The surface mesh includes a face of type other than edge, "
1262 "triangle or quadrilateral. This face type is not supported.";
1264 return "Not enough memory for the face table.";
1266 return "Not enough memory.";
1268 return "Not enough memory.";
1270 return "Face is ignored.";
1272 return "End of file. Some data are missing in the file.";
1274 return "Read error on the file. There are wrong data in the file.";
1276 return "the metric file is inadequate (dimension other than 3).";
1278 return "the metric file is inadequate (values not per vertices).";
1280 return "the metric file contains more than one field.";
1282 return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1283 "value of number of mesh vertices in the \".noboite\" file.";
1285 return "Too many sub-domains.";
1287 return "the number of vertices is negative or null.";
1289 return "the number of faces is negative or null.";
1291 return "A face has a null vertex.";
1293 return "incompatible data.";
1295 return "the number of vertices is negative or null.";
1297 return "the number of vertices is negative or null (in the \".mesh\" file).";
1299 return "the number of faces is negative or null.";
1301 return "A face appears more than once in the input surface mesh.";
1303 return "An edge appears more than once in the input surface mesh.";
1305 return "A face has a vertex negative or null.";
1307 return "NOT ENOUGH MEMORY.";
1309 return "Not enough available memory.";
1311 return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1312 "in terms of quality or the input list of points is wrong.";
1314 return "Some vertices are too close to one another or coincident.";
1316 return "Some vertices are too close to one another or coincident.";
1318 return "A vertex cannot be inserted.";
1320 return "There are at least two points considered as coincident.";
1322 return "Some vertices are too close to one another or coincident.";
1324 return "The surface mesh regeneration step has failed.";
1326 return "Constrained edge cannot be enforced.";
1328 return "Constrained face cannot be enforced.";
1330 return "Missing faces.";
1332 return "No guess to start the definition of the connected component(s).";
1334 return "The surface mesh includes at least one hole. The domain is not well defined.";
1336 return "Impossible to define a component.";
1338 return "The surface edge intersects another surface edge.";
1340 return "The surface edge intersects the surface face.";
1342 return "One boundary point lies within a surface face.";
1344 return "One surface edge intersects a surface face.";
1346 return "One boundary point lies within a surface edge.";
1348 return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1349 "to too many swaps.";
1351 return "Edge is unique (i.e., bounds a hole in the surface).";
1353 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1355 return "Too many components, too many sub-domain.";
1357 return "The surface mesh includes at least one hole. "
1358 "Therefore there is no domain properly defined.";
1360 return "Statistics.";
1362 return "Statistics.";
1364 return "Warning, it is dramatically tedious to enforce the boundary items.";
1366 return "Not enough memory at this time, nevertheless, the program continues. "
1367 "The expected mesh will be correct but not really as large as required.";
1369 return "see above error code, resulting quality may be poor.";
1371 return "Not enough memory at this time, nevertheless, the program continues (warning).";
1373 return "Unknown face type.";
1376 return "End of file. Some data are missing in the file.";
1378 return "A too small volume element is detected.";
1380 return "There exists at least a null or negative volume element.";
1382 return "There exist null or negative volume elements.";
1384 return "A too small volume element is detected. A face is considered being degenerated.";
1386 return "Some element is suspected to be very bad shaped or wrong.";
1388 return "A too bad quality face is detected. This face is considered degenerated.";
1390 return "A too bad quality face is detected. This face is degenerated.";
1392 return "Presumably, the surface mesh is not compatible with the domain being processed.";
1394 return "Abnormal error occured, contact hotline.";
1396 return "Not enough memory for the face table.";
1398 return "The algorithm cannot run further. "
1399 "The surface mesh is probably very bad in terms of quality.";
1401 return "Bad vertex number.";
1406 //================================================================================
1408 * \brief Retrieve from a string given number of integers
1410 //================================================================================
1412 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1415 ids.reserve( nbIds );
1418 while ( !isdigit( *ptr )) ++ptr;
1419 if ( ptr[-1] == '-' ) --ptr;
1420 ids.push_back( strtol( ptr, &ptr, 10 ));
1426 //================================================================================
1428 * \brief Retrieve problem description form a log file
1429 * \retval bool - always false
1431 //================================================================================
1433 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1434 const _Ghs2smdsConvertor & toSmdsConvertor )
1438 int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1440 int file = ::open (logFile.ToCString(), O_RDONLY);
1443 return error( SMESH_Comment("See ") << logFile << " for problem description");
1446 // struct stat status;
1447 // fstat(file, &status);
1448 // size_t length = status.st_size;
1449 off_t length = lseek( file, 0, SEEK_END);
1450 lseek( file, 0, SEEK_SET);
1453 vector< char > buf( length );
1454 int nBytesRead = ::read (file, & buf[0], length);
1456 char* ptr = & buf[0];
1457 char* bufEnd = ptr + nBytesRead;
1459 SMESH_Comment errDescription;
1461 enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1463 // look for errors "ERR #"
1465 set<string> foundErrorStr; // to avoid reporting same error several times
1466 set<int> elemErrorNums; // not to report different types of errors with bad elements
1467 while ( ++ptr < bufEnd )
1469 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1472 list<const SMDS_MeshElement*> badElems;
1473 vector<int> nodeIds;
1477 int errNum = strtol(ptr, &ptr, 10);
1478 switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1480 // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1481 ptr = getIds(ptr, NODE, nodeIds);
1482 ptr = getIds(ptr, TRIA, nodeIds);
1483 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1485 case 1000: // ERR 1000 : 1 3 2
1486 // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1487 ptr = getIds(ptr, TRIA, nodeIds);
1488 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1491 // Edge (e1, e2) appears more than once in the input surface mesh
1492 ptr = getIds(ptr, EDGE, nodeIds);
1493 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1496 // Face (f 1, f 2, f 3) has a vertex negative or null
1497 ptr = getIds(ptr, TRIA, nodeIds);
1498 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1501 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1502 ptr = getIds(ptr, NODE, nodeIds);
1503 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1504 ptr = getIds(ptr, NODE, nodeIds);
1505 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1508 // Vertex v1 cannot be inserted (warning).
1509 ptr = getIds(ptr, NODE, nodeIds);
1510 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1513 // There are at least two points whose distance is dist, i.e., considered as coincident
1514 case 2103: // ERR 2103 : 16 WITH 3
1515 // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1516 ptr = getIds(ptr, NODE, nodeIds);
1517 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1518 ptr = getIds(ptr, NODE, nodeIds);
1519 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1522 // Constrained edge (e1, e2) cannot be enforced (warning).
1523 ptr = getIds(ptr, EDGE, nodeIds);
1524 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1527 // Constrained face (f 1, f 2, f 3) cannot be enforced
1528 ptr = getIds(ptr, TRIA, nodeIds);
1529 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1531 case 3103: // ERR 3103 : 1 2 WITH 7 3
1532 // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1533 ptr = getIds(ptr, EDGE, nodeIds);
1534 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1535 ptr = getIds(ptr, EDGE, nodeIds);
1536 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1538 case 3104: // ERR 3104 : 9 10 WITH 1 2 3
1539 // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1540 ptr = getIds(ptr, EDGE, nodeIds);
1541 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1542 ptr = getIds(ptr, TRIA, nodeIds);
1543 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1545 case 3105: // ERR 3105 : 8 IN 2 3 5
1546 // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1547 ptr = getIds(ptr, NODE, nodeIds);
1548 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1549 ptr = getIds(ptr, TRIA, nodeIds);
1550 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1553 // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1554 ptr = getIds(ptr, EDGE, nodeIds);
1555 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1556 ptr = getIds(ptr, TRIA, nodeIds);
1557 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1559 case 3107: // ERR 3107 : 2 IN 4 1
1560 // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1561 ptr = getIds(ptr, NODE, nodeIds);
1562 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1563 ptr = getIds(ptr, EDGE, nodeIds);
1564 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1566 case 3109: // ERR 3109 : EDGE 5 6 UNIQUE
1567 // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1568 ptr = getIds(ptr, EDGE, nodeIds);
1569 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1571 case 9000: // ERR 9000
1572 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
1573 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
1574 // A too small volume element is detected. Are reported the index of the element,
1575 // its four vertex indices, its volume and the tolerance threshold value
1576 ptr = getIds(ptr, ID, nodeIds);
1577 ptr = getIds(ptr, VOL, nodeIds);
1578 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1579 // even if all nodes found, volume it most probably invisible,
1580 // add its faces to demenstrate it anyhow
1582 vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1583 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1584 faceNodes[2] = nodeIds[3]; // 013
1585 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1586 faceNodes[1] = nodeIds[2]; // 023
1587 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1588 faceNodes[0] = nodeIds[1]; // 123
1589 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1592 case 9001: // ERR 9001
1593 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
1594 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
1595 // %% NUMBER OF NULL VOLUME TETS : 0
1596 // There exists at least a null or negative volume element
1599 // There exist n null or negative volume elements
1602 // A too small volume element is detected
1605 // A too bad quality face is detected. This face is considered degenerated,
1606 // its index, its three vertex indices together with its quality value are reported
1607 break; // same as next
1608 case 9112: // ERR 9112
1609 // FACE 2 WITH VERTICES : 4 2 5
1610 // SMALL INRADIUS : 0.
1611 // A too bad quality face is detected. This face is degenerated,
1612 // its index, its three vertex indices together with its inradius are reported
1613 ptr = getIds(ptr, ID, nodeIds);
1614 ptr = getIds(ptr, TRIA, nodeIds);
1615 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1616 // add triangle edges as it most probably has zero area and hence invisible
1618 vector<int> edgeNodes(2);
1619 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1620 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1621 edgeNodes[1] = nodeIds[2]; // 0-2
1622 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1623 edgeNodes[0] = nodeIds[1]; // 1-2
1624 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1629 bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1631 continue; // not to report same error several times
1633 // const SMDS_MeshElement* nullElem = 0;
1634 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1636 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1637 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1638 // if ( oneMoreErrorType )
1639 // continue; // not to report different types of errors with bad elements
1642 // store bad elements
1643 //if ( allElemsOk ) {
1644 list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1645 for ( ; elem != badElems.end(); ++elem )
1646 addBadInputElement( *elem );
1650 string text = translateError( errNum );
1651 if ( errDescription.find( text ) == text.npos ) {
1652 if ( !errDescription.empty() )
1653 errDescription << "\n";
1654 errDescription << text;
1659 if ( errDescription.empty() ) { // no errors found
1660 char msg[] = "connection to server failed";
1661 if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1662 errDescription << "Licence problems.";
1665 if ( errDescription.empty() )
1666 errDescription << "See " << logFile << " for problem description";
1668 errDescription << "\nSee " << logFile << " for more information";
1670 return error( errDescription );
1673 //================================================================================
1675 * \brief Creates _Ghs2smdsConvertor
1677 //================================================================================
1679 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1680 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1684 //================================================================================
1686 * \brief Creates _Ghs2smdsConvertor
1688 //================================================================================
1690 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> & nodeByGhsId)
1691 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1695 //================================================================================
1697 * \brief Return SMDS element by ids of GHS3D nodes
1699 //================================================================================
1701 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1703 size_t nbNodes = ghsNodes.size();
1704 vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1705 for ( size_t i = 0; i < nbNodes; ++i ) {
1706 int ghsNode = ghsNodes[ i ];
1707 if ( _ghs2NodeMap ) {
1708 map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1709 if ( in == _ghs2NodeMap->end() )
1711 nodes[ i ] = in->second;
1714 if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1716 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1722 if ( nbNodes == 2 ) {
1723 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1725 edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1728 if ( nbNodes == 3 ) {
1729 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1731 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1735 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1741 //=============================================================================
1745 //=============================================================================
1746 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
1747 const TopoDS_Shape& aShape,
1748 MapShapeNbElems& aResMap)
1750 int nbtri = 0, nbqua = 0;
1751 double fullArea = 0.0;
1752 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1753 TopoDS_Face F = TopoDS::Face( exp.Current() );
1754 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1755 MapShapeNbElemsItr anIt = aResMap.find(sm);
1756 if( anIt==aResMap.end() ) {
1757 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1758 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1759 "Submesh can not be evaluated",this));
1762 std::vector<int> aVec = (*anIt).second;
1763 nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
1764 nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1766 BRepGProp::SurfaceProperties(F,G);
1767 double anArea = G.Mass();
1771 // collect info from edges
1772 int nb0d_e = 0, nb1d_e = 0;
1773 bool IsQuadratic = false;
1774 bool IsFirst = true;
1775 TopTools_MapOfShape tmpMap;
1776 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1777 TopoDS_Edge E = TopoDS::Edge(exp.Current());
1778 if( tmpMap.Contains(E) )
1781 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
1782 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
1783 std::vector<int> aVec = (*anIt).second;
1784 nb0d_e += aVec[SMDSEntity_Node];
1785 nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
1787 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1793 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
1796 BRepGProp::VolumeProperties(aShape,G);
1797 double aVolume = G.Mass();
1798 double tetrVol = 0.1179*ELen*ELen*ELen;
1799 double CoeffQuality = 0.9;
1800 int nbVols = int(aVolume/tetrVol/CoeffQuality);
1801 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
1802 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
1803 std::vector<int> aVec(SMDSEntity_Last);
1804 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1806 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
1807 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
1808 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
1811 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
1812 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
1813 aVec[SMDSEntity_Pyramid] = nbqua;
1815 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1816 aResMap.insert(std::make_pair(sm,aVec));