]> SALOME platform Git repositories - plugins/ghs3dplugin.git/blob - src/GHS3DPlugin/GHS3DPlugin_GHS3D.cxx
Salome HOME
#19896 [CEA] GHS3DPlugin_Hypothesis' object has no attribute 'SetPrintLogInFile'
[plugins/ghs3dplugin.git] / src / GHS3DPlugin / GHS3DPlugin_GHS3D.cxx
1 // Copyright (C) 2004-2020  CEA/DEN, EDF R&D
2 //
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, or (at your option) any later version.
7 //
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.
12 //
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
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 //=============================================================================
21 // File      : GHS3DPlugin_GHS3D.cxx
22 // Created   : 
23 // Author    : Edward AGAPOV, modified by Lioka RAZAFINDRAZAKA (CEA) 09/02/2007
24 // Project   : SALOME
25 //=============================================================================
26 //
27 #include "GHS3DPlugin_GHS3D.hxx"
28 #include "GHS3DPlugin_Hypothesis.hxx"
29 #include "MG_Tetra_API.hxx"
30
31 #include <SMDS_FaceOfNodes.hxx>
32 #include <SMDS_LinearEdge.hxx>
33 #include <SMDS_MeshElement.hxx>
34 #include <SMDS_MeshNode.hxx>
35 #include <SMDS_VolumeOfNodes.hxx>
36 #include <SMESHDS_Group.hxx>
37 #include <SMESHDS_Mesh.hxx>
38 #include <SMESH_Comment.hxx>
39 #include <SMESH_File.hxx>
40 #include <SMESH_Group.hxx>
41 #include <SMESH_HypoFilter.hxx>
42 #include <SMESH_Mesh.hxx>
43 #include <SMESH_MeshAlgos.hxx>
44 #include <SMESH_MeshEditor.hxx>
45 #include <SMESH_MesherHelper.hxx>
46 #include <SMESH_OctreeNode.hxx>
47 #include <SMESH_subMeshEventListener.hxx>
48 #include <StdMeshers_QuadToTriaAdaptor.hxx>
49 #include <StdMeshers_ViscousLayers.hxx>
50
51 #include <BRepAdaptor_Surface.hxx>
52 #include <BRepBndLib.hxx>
53 #include <BRepBuilderAPI_MakeVertex.hxx>
54 #include <BRepClass3d.hxx>
55 #include <BRepClass3d_SolidClassifier.hxx>
56 #include <BRepExtrema_DistShapeShape.hxx>
57 #include <BRepGProp.hxx>
58 #include <BRepTools.hxx>
59 #include <BRep_Tool.hxx>
60 #include <Bnd_Box.hxx>
61 #include <GProp_GProps.hxx>
62 #include <GeomAPI_ProjectPointOnSurf.hxx>
63 #include <Precision.hxx>
64 #include <Standard_ErrorHandler.hxx>
65 #include <Standard_Failure.hxx>
66 #include <Standard_ProgramError.hxx>
67 #include <TopExp.hxx>
68 #include <TopExp_Explorer.hxx>
69 #include <TopTools_IndexedMapOfShape.hxx>
70 #include <TopTools_ListIteratorOfListOfShape.hxx>
71 #include <TopTools_MapOfShape.hxx>
72 #include <TopoDS.hxx>
73 #include <TopoDS_Shell.hxx>
74 #include <TopoDS_Solid.hxx>
75
76 #include <Basics_Utils.hxx>
77 #include <utilities.h>
78
79 #include <algorithm>
80 #include <errno.h>
81
82 #ifdef _DEBUG_
83 //#define _MY_DEBUG_
84 #endif
85
86 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
87
88 using namespace std;
89
90 #define HOLE_ID -1
91
92 // flags returning state of enforced entities, returned from writeGMFFile
93 enum InvalidEnforcedFlags { FLAG_BAD_ENF_VERT = 1,
94                             FLAG_BAD_ENF_NODE = 2,
95                             FLAG_BAD_ENF_EDGE = 4,
96                             FLAG_BAD_ENF_TRIA = 8
97 };
98 static std::string flagsToErrorStr( int anInvalidEnforcedFlags )
99 {
100   std::string str;
101   if ( anInvalidEnforcedFlags != 0 )
102   {
103     if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_VERT )
104       str = "There are enforced vertices incorrectly defined.\n";
105     if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_NODE )
106       str += "There are enforced nodes incorrectly defined.\n";
107     if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_EDGE )
108       str += "There are enforced edge incorrectly defined.\n";
109     if ( anInvalidEnforcedFlags & FLAG_BAD_ENF_TRIA )
110       str += "There are enforced triangles incorrectly defined.\n";
111   }
112   return str;
113 }
114
115 typedef const list<const SMDS_MeshFace*> TTriaList;
116
117 static const char theDomainGroupNamePrefix[] = "Domain_";
118
119 static void removeFile( const TCollection_AsciiString& fileName )
120 {
121   try {
122     SMESH_File( fileName.ToCString() ).remove();
123   }
124   catch ( ... ) {
125     MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
126   }
127 }
128
129 //=============================================================================
130 /*!
131  *  
132  */
133 //=============================================================================
134
135 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, SMESH_Gen* gen)
136   : SMESH_3D_Algo(hypId, gen), _isLibUsed( false )
137 {
138   _name = Name();
139   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
140   _onlyUnaryInput = false; // Compute() will be called on a compound of solids
141   _iShape=0;
142   _nbShape=0;
143   _compatibleHypothesis.push_back( GHS3DPlugin_Hypothesis::GetHypType());
144   _compatibleHypothesis.push_back( StdMeshers_ViscousLayers::GetHypType() );
145   _requireShape = false; // can work without shape
146   
147   _computeCanceled = false;
148   _progressAdvance = 1e-4;
149 }
150
151 //=============================================================================
152 /*!
153  *
154  */
155 //=============================================================================
156
157 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
158 {
159 }
160
161 //=============================================================================
162 /*!
163  *
164  */
165 //=============================================================================
166
167 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh&         aMesh,
168                                           const TopoDS_Shape& aShape,
169                                           Hypothesis_Status&  aStatus )
170 {
171   aStatus = SMESH_Hypothesis::HYP_OK;
172
173   _hyp = 0;
174   _viscousLayersHyp = 0;
175   _keepFiles = false;
176   _removeLogOnSuccess = true;
177   _logInStandardOutput = false;
178
179   const list <const SMESHDS_Hypothesis * >& hyps =
180     GetUsedHypothesis(aMesh, aShape, /*ignoreAuxiliary=*/false);
181   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
182   for ( ; h != hyps.end(); ++h )
183   {
184     if ( !_hyp )
185       _hyp = dynamic_cast< const GHS3DPlugin_Hypothesis*> ( *h );
186     if ( !_viscousLayersHyp )
187       _viscousLayersHyp = dynamic_cast< const StdMeshers_ViscousLayers*> ( *h );
188   }
189   if ( _hyp )
190   {
191     _keepFiles           = _hyp->GetKeepFiles();
192     _removeLogOnSuccess  = _hyp->GetRemoveLogOnSuccess();
193     _logInStandardOutput = _hyp->GetStandardOutputLog();
194   }
195
196   if ( _viscousLayersHyp )
197     error( _viscousLayersHyp->CheckHypothesis( aMesh, aShape, aStatus ));
198
199   return aStatus == HYP_OK;
200 }
201
202
203 //=======================================================================
204 //function : entryToShape
205 //purpose  : 
206 //=======================================================================
207
208 TopoDS_Shape GHS3DPlugin_GHS3D::entryToShape(std::string entry)
209 {
210   if ( SMESH_Gen_i::getStudyServant()->_is_nil() )
211     throw SALOME_Exception("MG-Tetra plugin can't work w/o publishing in the study");
212
213   GEOM::GEOM_Object_var aGeomObj;
214   TopoDS_Shape S = TopoDS_Shape();
215   SALOMEDS::SObject_var aSObj = SMESH_Gen_i::getStudyServant()->FindObjectID( entry.c_str() );
216   if (!aSObj->_is_nil() ) {
217     CORBA::Object_var obj = aSObj->GetObject();
218     aGeomObj = GEOM::GEOM_Object::_narrow(obj);
219     aSObj->UnRegister();
220   }
221   if ( !aGeomObj->_is_nil() )
222     S = SMESH_Gen_i::GetSMESHGen()->GeomObjectToShape( aGeomObj.in() );
223   return S;
224 }
225
226 //================================================================================
227 /*!
228  * \brief returns id of a solid if a triangle defined by the nodes is a temporary face
229  * either on a side facet of pyramid or a top of pentahedron and defines sub-domian
230  * outside the volume; else returns HOLE_ID
231  */
232 //================================================================================
233
234 static int checkTmpFace(const SMDS_MeshNode* node1,
235                         const SMDS_MeshNode* node2,
236                         const SMDS_MeshNode* node3)
237 {
238   // find a pyramid sharing the 3 nodes
239   SMDS_ElemIteratorPtr vIt1 = node1->GetInverseElementIterator(SMDSAbs_Volume);
240   while ( vIt1->more() )
241   {
242     const SMDS_MeshElement* vol = vIt1->next();
243     const int           nbNodes = vol->NbCornerNodes();
244     if ( nbNodes != 5 && nbNodes != 6 ) continue;
245     int i2, i3;
246     if ( (i2 = vol->GetNodeIndex( node2 )) >= 0 &&
247          (i3 = vol->GetNodeIndex( node3 )) >= 0 )
248     {
249       if ( nbNodes == 5 )
250       {
251         // Triangle defines sub-domian inside the pyramid if it's
252         // normal points out of the vol
253
254         // make i2 and i3 hold indices of base nodes of the vol while
255         // keeping the nodes order in the triangle
256         const int iApex = 4;
257         if ( i2 == iApex )
258           i2 = i3, i3 = vol->GetNodeIndex( node1 );
259         else if ( i3 == iApex )
260           i3 = i2, i2 = vol->GetNodeIndex( node1 );
261
262         int i3base = (i2+1) % 4; // next index after i2 within the pyramid base
263         bool isDomainInPyramid = ( i3base != i3 );
264         return isDomainInPyramid ? HOLE_ID : vol->getshapeId();
265       }
266       else
267       {
268         return vol->getshapeId(); // triangle is a prism top
269       }
270     }
271   }
272   return HOLE_ID;
273 }
274
275 //=======================================================================
276 //function : findShapeID
277 //purpose  : find the solid corresponding to MG-Tetra sub-domain following
278 //           the technique proposed in MG-Tetra manual (available within
279 //           MG-Tetra installation) in chapter "B.4 Subdomain (sub-region) assignment".
280 //           In brief: normal of the triangle defined by the given nodes
281 //           points out of the domain it is associated to
282 //=======================================================================
283
284 static int findShapeID(SMESH_Mesh&          mesh,
285                        const SMDS_MeshNode* node1,
286                        const SMDS_MeshNode* node2,
287                        const SMDS_MeshNode* node3,
288                        const bool           toMeshHoles)
289 {
290   const int invalidID = 0;
291   SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
292
293   // face the nodes belong to
294   vector<const SMDS_MeshNode *> nodes(3);
295   nodes[0] = node1;
296   nodes[1] = node2;
297   nodes[2] = node3;
298   const SMDS_MeshElement * face = meshDS->FindElement( nodes, SMDSAbs_Face, /*noMedium=*/true);
299   if ( !face )
300     return checkTmpFace(node1, node2, node3);
301 #ifdef _MY_DEBUG_
302   std::cout << "bnd face " << face->GetID() << " - ";
303 #endif
304   // geom face the face assigned to
305   SMESH_MeshEditor editor(&mesh);
306   int geomFaceID = editor.FindShape( face );
307   if ( !geomFaceID )
308     return checkTmpFace(node1, node2, node3);
309   TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
310   if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
311     return invalidID;
312   TopoDS_Face geomFace = TopoDS::Face( shape );
313
314   // solids bounded by geom face
315   TopTools_IndexedMapOfShape solids, shells;
316   TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
317   for ( ; ansIt.More(); ansIt.Next() ) {
318     switch ( ansIt.Value().ShapeType() ) {
319     case TopAbs_SOLID:
320       solids.Add( ansIt.Value() ); break;
321     case TopAbs_SHELL:
322       shells.Add( ansIt.Value() ); break;
323     default:;
324     }
325   }
326   // analyse found solids
327   if ( solids.Extent() == 0 || shells.Extent() == 0)
328     return invalidID;
329
330   const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
331   if ( solids.Extent() == 1 )
332   {
333     if ( toMeshHoles )
334       return meshDS->ShapeToIndex( solid1 );
335
336     // - Are we at a hole boundary face?
337     if ( shells(1).IsSame( BRepClass3d::OuterShell( solid1 )) )
338     { // - No, but maybe a hole is bound by two shapes? Does shells(1) touch another shell?
339       bool touch = false;
340       TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
341       // check if any edge of shells(1) belongs to another shell
342       for ( ; eExp.More() && !touch; eExp.Next() ) {
343         ansIt = mesh.GetAncestors( eExp.Current() );
344         for ( ; ansIt.More() && !touch; ansIt.Next() ) {
345           if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
346             touch = ( !ansIt.Value().IsSame( shells(1) ));
347         }
348       }
349       if (!touch)
350         return meshDS->ShapeToIndex( solid1 );
351     }
352   }
353   // find orientation of geom face within the first solid
354   TopExp_Explorer fExp( solid1, TopAbs_FACE );
355   for ( ; fExp.More(); fExp.Next() )
356     if ( geomFace.IsSame( fExp.Current() )) {
357       geomFace = TopoDS::Face( fExp.Current() );
358       break;
359     }
360   if ( !fExp.More() )
361     return invalidID; // face not found
362
363   // normale to triangle
364   gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
365   gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
366   gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
367   gp_Vec vec12( node1Pnt, node2Pnt );
368   gp_Vec vec13( node1Pnt, node3Pnt );
369   gp_Vec meshNormal = vec12 ^ vec13;
370   if ( meshNormal.SquareMagnitude() < DBL_MIN )
371     return invalidID;
372
373   // get normale to geomFace at any node
374   bool geomNormalOK = false;
375   gp_Vec geomNormal;
376   SMESH_MesherHelper helper( mesh ); helper.SetSubShape( geomFace );
377   for ( int i = 0; !geomNormalOK && i < 3; ++i )
378   {
379     // find UV of i-th node on geomFace
380     const SMDS_MeshNode* nNotOnSeamEdge = 0;
381     if ( helper.IsSeamShape( nodes[i]->getshapeId() )) {
382       if ( helper.IsSeamShape( nodes[(i+1)%3]->getshapeId() ))
383         nNotOnSeamEdge = nodes[(i+2)%3];
384       else
385         nNotOnSeamEdge = nodes[(i+1)%3];
386     }
387     bool uvOK = true;
388     gp_XY uv = helper.GetNodeUV( geomFace, nodes[i], nNotOnSeamEdge, &uvOK );
389     // check that uv is correct
390     if (uvOK) {
391       double tol = 1e-6;
392       TopoDS_Shape nodeShape = helper.GetSubShapeByNode( nodes[i], meshDS );
393       if ( !nodeShape.IsNull() )
394         switch ( nodeShape.ShapeType() )
395         {
396         case TopAbs_FACE:   tol = BRep_Tool::Tolerance( TopoDS::Face( nodeShape )); break;
397         case TopAbs_EDGE:   tol = BRep_Tool::Tolerance( TopoDS::Edge( nodeShape )); break;
398         case TopAbs_VERTEX: tol = BRep_Tool::Tolerance( TopoDS::Vertex( nodeShape )); break;
399         default:;
400         }
401       gp_Pnt nodePnt ( nodes[i]->X(), nodes[i]->Y(), nodes[i]->Z() );
402       BRepAdaptor_Surface surface( geomFace );
403       uvOK = ( nodePnt.Distance( surface.Value( uv.X(), uv.Y() )) < 2 * tol );
404       if ( uvOK ) {
405         // normale to geomFace at UV
406         gp_Vec du, dv;
407         surface.D1( uv.X(), uv.Y(), nodePnt, du, dv );
408         geomNormal = du ^ dv;
409         if ( geomFace.Orientation() == TopAbs_REVERSED )
410           geomNormal.Reverse();
411         geomNormalOK = ( geomNormal.SquareMagnitude() > DBL_MIN * 1e3 );
412       }
413     }
414   }
415   if ( !geomNormalOK)
416     return invalidID;
417
418   // compare normals
419   bool isReverse = ( meshNormal * geomNormal ) < 0;
420   if ( !isReverse )
421     return meshDS->ShapeToIndex( solid1 );
422
423   if ( solids.Extent() == 1 )
424     return HOLE_ID; // we are inside a hole
425   else
426     return meshDS->ShapeToIndex( solids(2) );
427 }
428
429 //=======================================================================
430 //function : addElemInMeshGroup
431 //purpose  : Update or create groups in mesh
432 //=======================================================================
433
434 static void addElemInMeshGroup(SMESH_Mesh*             theMesh,
435                                const SMDS_MeshElement* anElem,
436                                std::string&            groupName,
437                                std::set<std::string>&  groupsToRemove)
438 {
439   if ( !anElem ) return; // issue 0021776
440
441   bool groupDone = false;
442   SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
443   while (grIt->more()) {
444     SMESH_Group * group = grIt->next();
445     if ( !group ) continue;
446     SMESHDS_GroupBase* groupDS = group->GetGroupDS();
447     if ( !groupDS ) continue;
448     if ( groupDS->GetType()==anElem->GetType() &&groupName.compare(group->GetName())==0) {
449       SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
450       aGroupDS->SMDSGroup().Add(anElem);
451       groupDone = true;
452       break;
453     }
454   }
455   
456   if (!groupDone)
457   {
458     SMESH_Group* aGroup = theMesh->AddGroup(anElem->GetType(), groupName.c_str());
459     aGroup->SetName( groupName.c_str() );
460     SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
461     aGroupDS->SMDSGroup().Add(anElem);
462     groupDone = true;
463   }
464   if (!groupDone)
465     throw SALOME_Exception(LOCALIZED("A given element was not added to a group"));
466 }
467
468
469 //=======================================================================
470 //function : updateMeshGroups
471 //purpose  : Update or create groups in mesh
472 //=======================================================================
473
474 static void updateMeshGroups(SMESH_Mesh* theMesh, std::set<std::string> groupsToRemove)
475 {
476   SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
477   while (grIt->more()) {
478     SMESH_Group * group = grIt->next();
479     if ( !group ) continue;
480     SMESHDS_GroupBase* groupDS = group->GetGroupDS();
481     if ( !groupDS ) continue;
482     std::string currentGroupName = (string)group->GetName();
483     if (groupDS->IsEmpty() && groupsToRemove.find(currentGroupName) != groupsToRemove.end()) {
484       // Previous group created by enforced elements
485       theMesh->RemoveGroup(groupDS->GetID());
486     }
487   }
488 }
489
490 //=======================================================================
491 //function : removeEmptyGroupsOfDomains
492 //purpose  : remove empty groups named "Domain_nb" created due to
493 //           "To make groups of domains" option.
494 //=======================================================================
495
496 static void removeEmptyGroupsOfDomains(SMESH_Mesh* mesh,
497                                        bool        notEmptyAsWell = false)
498 {
499   const char* refName = theDomainGroupNamePrefix;
500   const size_t refLen = strlen( theDomainGroupNamePrefix );
501
502   std::list<int> groupIDs = mesh->GetGroupIds();
503   std::list<int>::const_iterator id = groupIDs.begin();
504   for ( ; id != groupIDs.end(); ++id )
505   {
506     SMESH_Group* group = mesh->GetGroup( *id );
507     if ( !group || ( !group->GetGroupDS()->IsEmpty() && !notEmptyAsWell ))
508       continue;
509     const char* name = group->GetName();
510     char* end;
511     // check the name
512     if ( strncmp( name, refName, refLen ) == 0 && // starts from refName;
513          isdigit( *( name + refLen )) &&          // refName is followed by a digit;
514          strtol( name + refLen, &end, 10) >= 0 && // there are only digits ...
515          *end == '\0')                            // ... till a string end.
516     {
517       mesh->RemoveGroup( *id );
518     }
519   }
520 }
521
522 //================================================================================
523 /*!
524  * \brief Create the groups corresponding to domains
525  */
526 //================================================================================
527
528 static void makeDomainGroups( std::vector< std::vector< const SMDS_MeshElement* > >& elemsOfDomain,
529                               SMESH_MesherHelper*                                    theHelper)
530 {
531   // int nbDomains = 0;
532   // for ( size_t i = 0; i < elemsOfDomain.size(); ++i )
533   //   nbDomains += ( elemsOfDomain[i].size() > 0 );
534
535   // if ( nbDomains > 1 )
536   for ( size_t iDomain = 0; iDomain < elemsOfDomain.size(); ++iDomain )
537   {
538     std::vector< const SMDS_MeshElement* > & elems = elemsOfDomain[ iDomain ];
539     if ( elems.empty() ) continue;
540
541     // find existing groups
542     std::vector< SMESH_Group* > groupOfType( SMDSAbs_NbElementTypes, (SMESH_Group*)NULL );
543     const std::string domainName = ( SMESH_Comment( theDomainGroupNamePrefix ) << iDomain );
544     SMESH_Mesh::GroupIteratorPtr groupIt = theHelper->GetMesh()->GetGroups();
545     while ( groupIt->more() )
546     {
547       SMESH_Group* group = groupIt->next();
548       if ( domainName == group->GetName() &&
549            dynamic_cast< SMESHDS_Group* >( group->GetGroupDS()) )
550         groupOfType[ group->GetGroupDS()->GetType() ] = group;
551     }
552     // create and fill the groups
553     size_t iElem = 0;
554     do
555     {
556       SMESH_Group* group = groupOfType[ elems[ iElem ]->GetType() ];
557       if ( !group )
558         group = theHelper->GetMesh()->AddGroup( elems[ iElem ]->GetType(),
559                                                 domainName.c_str() );
560       SMDS_MeshGroup& groupDS =
561         static_cast< SMESHDS_Group* >( group->GetGroupDS() )->SMDSGroup();
562
563       while ( iElem < elems.size() && groupDS.Add( elems[iElem] ))
564         ++iElem;
565
566     } while ( iElem < elems.size() );
567   }
568 }
569
570 //=======================================================================
571 //function : readGMFFile
572 //purpose  : read GMF file w/o geometry associated to mesh
573 //=======================================================================
574
575 static bool readGMFFile(MG_Tetra_API*                   MGOutput,
576                         const char*                     theFile,
577                         GHS3DPlugin_GHS3D*              theAlgo,
578                         SMESH_MesherHelper*             theHelper,
579                         std::vector <const SMDS_MeshNode*> &    theNodeByGhs3dId,
580                         std::vector <const SMDS_MeshElement*> & theFaceByGhs3dId,
581                         map<const SMDS_MeshNode*,int> & theNodeToGhs3dIdMap,
582                         std::vector<std::string> &      aNodeGroupByGhs3dId,
583                         std::vector<std::string> &      anEdgeGroupByGhs3dId,
584                         std::vector<std::string> &      aFaceGroupByGhs3dId,
585                         std::set<std::string> &         groupsToRemove,
586                         bool                            toMakeGroupsOfDomains=false,
587                         bool                            toMeshHoles=true)
588 {
589   std::string tmpStr;
590   SMESHDS_Mesh* theMeshDS = theHelper->GetMeshDS();
591   const bool hasGeom = ( theHelper->GetMesh()->HasShapeToMesh() );
592
593   int nbInitialNodes = theNodeByGhs3dId.size();
594   
595 #ifdef _MY_DEBUG_
596   const bool isQuadMesh = 
597     theHelper->GetMesh()->NbEdges( ORDER_QUADRATIC ) ||
598     theHelper->GetMesh()->NbFaces( ORDER_QUADRATIC ) ||
599     theHelper->GetMesh()->NbVolumes( ORDER_QUADRATIC );
600   std::cout << "theNodeByGhs3dId.size(): " << nbInitialNodes << std::endl;
601   std::cout << "theHelper->GetMesh()->NbNodes(): " << theMeshDS->NbNodes() << std::endl;
602   std::cout << "isQuadMesh: " << isQuadMesh << std::endl;
603 #endif
604   
605   // ---------------------------------
606   // Read generated elements and nodes
607   // ---------------------------------
608
609   int nbElem = 0, nbRef = 0;
610   int aGMFNodeID = 0;
611   std::vector< const SMDS_MeshNode*> GMFNode;
612 #ifdef _MY_DEBUG_
613   std::map<int, std::set<int> > subdomainId2tetraId;
614 #endif
615   std::map <GmfKwdCod,int> tabRef;
616   const bool force3d = !hasGeom;
617   const int  noID    = 0;
618
619   tabRef[GmfVertices]       = 3; // for new nodes and enforced nodes
620   tabRef[GmfCorners]        = 1;
621   tabRef[GmfEdges]          = 2; // for enforced edges
622   tabRef[GmfRidges]         = 1;
623   tabRef[GmfTriangles]      = 3; // for enforced faces
624   tabRef[GmfQuadrilaterals] = 4;
625   tabRef[GmfTetrahedra]     = 4; // for new tetras
626   tabRef[GmfHexahedra]      = 8;
627
628   int ver, dim;
629   int InpMsh = MGOutput->GmfOpenMesh( theFile, GmfRead, &ver, &dim);
630   if (!InpMsh)
631     return false;
632
633   // Read ids of domains
634   vector< int > solidIDByDomain;
635   if ( hasGeom )
636   {
637     int solid1; // id used in case of 1 domain or some reading failure
638     if ( theHelper->GetSubShape().ShapeType() == TopAbs_SOLID )
639       solid1 = theHelper->GetSubShapeID();
640     else
641       solid1 = theMeshDS->ShapeToIndex
642         ( TopExp_Explorer( theHelper->GetSubShape(), TopAbs_SOLID ).Current() );
643
644     int nbDomains = MGOutput->GmfStatKwd( InpMsh, GmfSubDomainFromGeom );
645     if ( nbDomains > 1 )
646     {
647       solidIDByDomain.resize( nbDomains+1, theHelper->GetSubShapeID() );
648       int faceNbNodes, faceIndex, orientation, domainNb;
649       MGOutput->GmfGotoKwd( InpMsh, GmfSubDomainFromGeom );
650       for ( int i = 0; i < nbDomains; ++i )
651       {
652         faceIndex = 0;
653         MGOutput->GmfGetLin( InpMsh, GmfSubDomainFromGeom,
654                    &faceNbNodes, &faceIndex, &orientation, &domainNb, i);
655         solidIDByDomain[ domainNb ] = 1;
656         if ( 0 < faceIndex && faceIndex-1 < (int)theFaceByGhs3dId.size() )
657         {
658           const SMDS_MeshElement* face = theFaceByGhs3dId[ faceIndex-1 ];
659           const SMDS_MeshNode* nn[3] = { face->GetNode(0),
660                                          face->GetNode(1),
661                                          face->GetNode(2) };
662           if ( orientation < 0 )
663             std::swap( nn[1], nn[2] );
664           solidIDByDomain[ domainNb ] =
665             findShapeID( *theHelper->GetMesh(), nn[0], nn[1], nn[2], toMeshHoles );
666           if ( solidIDByDomain[ domainNb ] > 0 )
667           {
668 #ifdef _MY_DEBUG_
669             std::cout << "solid " << solidIDByDomain[ domainNb ] << std::endl;
670 #endif
671             const TopoDS_Shape& foundShape = theMeshDS->IndexToShape( solidIDByDomain[ domainNb ] );
672             if ( ! theHelper->IsSubShape( foundShape, theHelper->GetSubShape() ))
673               solidIDByDomain[ domainNb ] = HOLE_ID;
674           }
675         }
676       }
677     }
678     if ( solidIDByDomain.size() < 2 )
679       solidIDByDomain.resize( 2, solid1 );
680   }
681
682   // Issue 0020682. Avoid creating nodes and tetras at place where
683   // volumic elements already exist
684   SMESH_ElementSearcher* elemSearcher = 0;
685   std::vector< const SMDS_MeshElement* > foundVolumes;
686   if ( !hasGeom && theHelper->GetMesh()->NbVolumes() > 0 )
687     elemSearcher = SMESH_MeshAlgos::GetElementSearcher( *theMeshDS );
688   auto_ptr< SMESH_ElementSearcher > elemSearcherDeleter( elemSearcher );
689
690   // IMP 0022172: [CEA 790] create the groups corresponding to domains
691   std::vector< std::vector< const SMDS_MeshElement* > > elemsOfDomain;
692
693   int nbVertices = MGOutput->GmfStatKwd( InpMsh, GmfVertices ) - nbInitialNodes;
694   if ( nbVertices < 0 )
695     return false;
696   GMFNode.resize( nbVertices + 1 );
697
698   std::map <GmfKwdCod,int>::const_iterator it = tabRef.begin();
699   for ( ; it != tabRef.end() ; ++it)
700   {
701     if(theAlgo->computeCanceled()) {
702       return false;
703     }
704     int dummy, solidID;
705     GmfKwdCod token = it->first;
706     nbRef           = it->second;
707
708     nbElem = MGOutput->GmfStatKwd( InpMsh, token);
709     if (nbElem > 0) {
710       MGOutput->GmfGotoKwd( InpMsh, token);
711       std::cout << "Read " << nbElem;
712     }
713     else
714       continue;
715
716     std::vector<int> id (nbElem*tabRef[token]); // node ids
717     std::vector<int> domainID( nbElem ); // domain
718
719     if (token == GmfVertices) {
720       (nbElem <= 1) ? tmpStr = " vertex" : tmpStr = " vertices";
721 //       std::cout << nbInitialNodes << " from input mesh " << std::endl;
722
723       // Remove orphan nodes from previous enforced mesh which was cleared
724 //       if ( nbElem < nbMeshNodes ) {
725 //         const SMDS_MeshNode* node;
726 //         SMDS_NodeIteratorPtr nodeIt = theMeshDS->nodesIterator();
727 //         while ( nodeIt->more() )
728 //         {
729 //           node = nodeIt->next();
730 //           if (theNodeToGhs3dIdMap.find(node) != theNodeToGhs3dIdMap.end())
731 //             theMeshDS->RemoveNode(node);
732 //         }
733 //       }
734
735       
736       int aGMFID;
737
738       float VerTab_f[3];
739       double x, y, z;
740       const SMDS_MeshNode * aGMFNode;
741
742       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
743         if(theAlgo->computeCanceled()) {
744           return false;
745         }
746         if (ver == GmfFloat) {
747           MGOutput->GmfGetLin( InpMsh, token, &VerTab_f[0], &VerTab_f[1], &VerTab_f[2], &dummy);
748           x = VerTab_f[0];
749           y = VerTab_f[1];
750           z = VerTab_f[2];
751         }
752         else {
753           MGOutput->GmfGetLin( InpMsh, token, &x, &y, &z, &dummy);
754         }
755         if (iElem >= nbInitialNodes) {
756           if ( elemSearcher &&
757                 elemSearcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_Volume, foundVolumes))
758             aGMFNode = 0;
759           else
760             aGMFNode = theHelper->AddNode(x, y, z);
761           
762           aGMFID = iElem -nbInitialNodes +1;
763           GMFNode[ aGMFID ] = aGMFNode;
764           if (aGMFID-1 < (int)aNodeGroupByGhs3dId.size() && !aNodeGroupByGhs3dId.at(aGMFID-1).empty())
765             addElemInMeshGroup(theHelper->GetMesh(), aGMFNode, aNodeGroupByGhs3dId.at(aGMFID-1), groupsToRemove);
766         }
767       }
768     }
769     else if (token == GmfCorners && nbElem > 0) {
770       (nbElem <= 1) ? tmpStr = " corner" : tmpStr = " corners";
771       for ( int iElem = 0; iElem < nbElem; iElem++ )
772         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]]);
773     }
774     else if (token == GmfRidges && nbElem > 0) {
775       (nbElem <= 1) ? tmpStr = " ridge" : tmpStr = " ridges";
776       for ( int iElem = 0; iElem < nbElem; iElem++ )
777         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]]);
778     }
779     else if (token == GmfEdges && nbElem > 0) {
780       (nbElem <= 1) ? tmpStr = " edge" : tmpStr = " edges";
781       for ( int iElem = 0; iElem < nbElem; iElem++ )
782         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &domainID[iElem]);
783     }
784     else if (token == GmfTriangles && nbElem > 0) {
785       (nbElem <= 1) ? tmpStr = " triangle" : tmpStr = " triangles";
786       for ( int iElem = 0; iElem < nbElem; iElem++ )
787         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &domainID[iElem]);
788     }
789     else if (token == GmfQuadrilaterals && nbElem > 0) {
790       (nbElem <= 1) ? tmpStr = " Quadrilateral" : tmpStr = " Quadrilaterals";
791       for ( int iElem = 0; iElem < nbElem; iElem++ )
792         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3], &domainID[iElem]);
793     }
794     else if (token == GmfTetrahedra && nbElem > 0) {
795       (nbElem <= 1) ? tmpStr = " Tetrahedron" : tmpStr = " Tetrahedra";
796       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
797         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3], &domainID[iElem]);
798 #ifdef _MY_DEBUG_
799         subdomainId2tetraId[dummy].insert(iElem+1);
800 #endif
801       }
802     }
803     else if (token == GmfHexahedra && nbElem > 0) {
804       (nbElem <= 1) ? tmpStr = " Hexahedron" : tmpStr = " Hexahedra";
805       for ( int iElem = 0; iElem < nbElem; iElem++ )
806         MGOutput->GmfGetLin( InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
807                   &id[iElem*tabRef[token]+4], &id[iElem*tabRef[token]+5], &id[iElem*tabRef[token]+6], &id[iElem*tabRef[token]+7], &domainID[iElem]);
808     }
809     std::cout << tmpStr << std::endl;
810     std::cout << std::endl;
811
812     switch (token) {
813     case GmfCorners:
814     case GmfRidges:
815     case GmfEdges:
816     case GmfTriangles:
817     case GmfQuadrilaterals:
818     case GmfTetrahedra:
819     case GmfHexahedra:
820     {
821       std::vector< const SMDS_MeshNode* > node( nbRef );
822       std::vector< int >          nodeID( nbRef );
823       std::vector< SMDS_MeshNode* > enfNode( nbRef );
824       const SMDS_MeshElement* aCreatedElem;
825
826       for ( int iElem = 0; iElem < nbElem; iElem++ )
827       {
828         if(theAlgo->computeCanceled()) {
829           return false;
830         }
831         // Check if elem is already in input mesh. If yes => skip
832         bool fullyCreatedElement = false; // if at least one of the nodes was created
833         for ( int iRef = 0; iRef < nbRef; iRef++ )
834         {
835           aGMFNodeID = id[iElem*tabRef[token]+iRef]; // read nbRef aGMFNodeID
836           if (aGMFNodeID <= nbInitialNodes) // input nodes
837           {
838             aGMFNodeID--;
839             node[ iRef ] = theNodeByGhs3dId[aGMFNodeID];
840           }
841           else
842           {
843             fullyCreatedElement = true;
844             aGMFNodeID -= nbInitialNodes;
845             nodeID[ iRef ] = aGMFNodeID ;
846             node  [ iRef ] = GMFNode[ aGMFNodeID ];
847           }
848         }
849         aCreatedElem = 0;
850         switch (token)
851         {
852         case GmfEdges:
853           if (fullyCreatedElement) {
854             aCreatedElem = theHelper->AddEdge( node[0], node[1], noID, force3d );
855             if (anEdgeGroupByGhs3dId.size() && !anEdgeGroupByGhs3dId[iElem].empty())
856               addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, anEdgeGroupByGhs3dId[iElem], groupsToRemove);
857           }
858           break;
859         case GmfTriangles:
860           if (fullyCreatedElement) {
861             aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], noID, force3d );
862             if (aFaceGroupByGhs3dId.size() && !aFaceGroupByGhs3dId[iElem].empty())
863               addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, aFaceGroupByGhs3dId[iElem], groupsToRemove);
864           }
865           break;
866         case GmfQuadrilaterals:
867           if (fullyCreatedElement) {
868             aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], node[3], noID, force3d );
869           }
870           break;
871         case GmfTetrahedra:
872           if ( hasGeom )
873           {
874             solidID = solidIDByDomain[ domainID[iElem]];
875             if ( solidID != HOLE_ID )
876             {
877               aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
878                                                    noID, force3d );
879               theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
880               for ( int iN = 0; iN < 4; ++iN )
881                 if ( node[iN]->getshapeId() < 1 )
882                   theMeshDS->SetNodeInVolume( node[iN], solidID );
883             }
884           }
885           else
886           {
887             if ( elemSearcher ) {
888               // Issue 0020682. Avoid creating nodes and tetras at place where
889               // volumic elements already exist
890               if ( !node[1] || !node[0] || !node[2] || !node[3] )
891                 continue;
892               if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
893                                                       SMESH_TNodeXYZ(node[1]) +
894                                                       SMESH_TNodeXYZ(node[2]) +
895                                                       SMESH_TNodeXYZ(node[3]) ) / 4.,
896                                                      SMDSAbs_Volume, foundVolumes ))
897                 break;
898             }
899             aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
900                                                  noID, force3d );
901           }
902           break;
903         case GmfHexahedra:
904           if ( hasGeom )
905           {
906             solidID = solidIDByDomain[ domainID[iElem]];
907             if ( solidID != HOLE_ID )
908             {
909               aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
910                                                    node[4], node[7], node[6], node[5],
911                                                    noID, force3d );
912               theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
913               for ( int iN = 0; iN < 8; ++iN )
914                 if ( node[iN]->getshapeId() < 1 )
915                   theMeshDS->SetNodeInVolume( node[iN], solidID );
916             }
917           }
918           else
919           {
920             if ( elemSearcher ) {
921               // Issue 0020682. Avoid creating nodes and tetras at place where
922               // volumic elements already exist
923               if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] || !node[6] || !node[7])
924                 continue;
925               if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
926                                                       SMESH_TNodeXYZ(node[1]) +
927                                                       SMESH_TNodeXYZ(node[2]) +
928                                                       SMESH_TNodeXYZ(node[3]) +
929                                                       SMESH_TNodeXYZ(node[4]) +
930                                                       SMESH_TNodeXYZ(node[5]) +
931                                                       SMESH_TNodeXYZ(node[6]) +
932                                                       SMESH_TNodeXYZ(node[7])) / 8.,
933                                                      SMDSAbs_Volume, foundVolumes ))
934                 break;
935             }
936             aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
937                                                  node[4], node[7], node[6], node[5],
938                                                  noID, force3d );
939           }
940           break;
941         default: continue;
942         } // switch (token)
943
944         // care about medium nodes
945         if ( aCreatedElem &&
946              aCreatedElem->IsQuadratic() &&
947              ( solidID = aCreatedElem->getshapeId() ) > 0 )
948         {
949           int iN = aCreatedElem->NbCornerNodes(), nbN = aCreatedElem->NbNodes();
950           for ( ; iN < nbN; ++iN )
951           {
952             const SMDS_MeshNode* n = aCreatedElem->GetNode(iN);
953             if ( n->getshapeId() < 1 )
954               theMeshDS->SetNodeInVolume( n, solidID );
955           }
956         }
957
958         if ( aCreatedElem && toMakeGroupsOfDomains )
959         {
960           if ( domainID[iElem] >= (int) elemsOfDomain.size() )
961             elemsOfDomain.resize( domainID[iElem] + 1 );
962           elemsOfDomain[ domainID[iElem] ].push_back( aCreatedElem );
963         }
964       } // loop on elements of one type
965       break;
966     } // case ...
967     default:;
968     } // switch (token)
969   } // loop on tabRef
970
971   // remove nodes in holes
972   if ( hasGeom )
973   {
974     for ( int i = 1; i <= nbVertices; ++i )
975       if ( GMFNode[i]->NbInverseElements() == 0 )
976         theMeshDS->RemoveFreeNode( GMFNode[i], /*sm=*/0, /*fromGroups=*/false );
977   }
978
979   MGOutput->GmfCloseMesh( InpMsh);
980
981   // 0022172: [CEA 790] create the groups corresponding to domains
982   if ( toMakeGroupsOfDomains )
983     makeDomainGroups( elemsOfDomain, theHelper );
984
985 #ifdef _MY_DEBUG_
986   std::map<int, std::set<int> >::const_iterator subdomainIt = subdomainId2tetraId.begin();
987   TCollection_AsciiString aSubdomainFileName = theFile;
988   aSubdomainFileName = aSubdomainFileName + ".subdomain";
989   ofstream aSubdomainFile  ( aSubdomainFileName.ToCString()  , ios::out);
990
991   aSubdomainFile << "Nb subdomains " << subdomainId2tetraId.size() << std::endl;
992   for(;subdomainIt != subdomainId2tetraId.end() ; ++subdomainIt) {
993     int subdomainId = subdomainIt->first;
994     std::set<int> tetraIds = subdomainIt->second;
995     std::set<int>::const_iterator tetraIdsIt = tetraIds.begin();
996     aSubdomainFile << subdomainId << std::endl;
997     for(;tetraIdsIt != tetraIds.end() ; ++tetraIdsIt) {
998       aSubdomainFile << (*tetraIdsIt) << " ";
999     }
1000     aSubdomainFile << std::endl;
1001   }
1002   aSubdomainFile.close();
1003 #endif  
1004   
1005   return true;
1006 }
1007
1008
1009 static bool writeGMFFile(MG_Tetra_API*                                   MGInput,
1010                          const char*                                     theMeshFileName,
1011                          const char*                                     theRequiredFileName,
1012                          const char*                                     theSolFileName,
1013                          const SMESH_ProxyMesh&                          theProxyMesh,
1014                          SMESH_MesherHelper&                             theHelper,
1015                          std::vector <const SMDS_MeshNode*> &            theNodeByGhs3dId,
1016                          std::vector <const SMDS_MeshElement*> &         theFaceByGhs3dId,
1017                          std::map<const SMDS_MeshNode*,int> &            aNodeToGhs3dIdMap,
1018                          std::vector<std::string> &                      aNodeGroupByGhs3dId,
1019                          std::vector<std::string> &                      anEdgeGroupByGhs3dId,
1020                          std::vector<std::string> &                      aFaceGroupByGhs3dId,
1021                          GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap & theEnforcedNodes,
1022                          GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedEdges,
1023                          GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedTriangles,
1024                          std::map<std::vector<double>, std::string> &    enfVerticesWithGroup,
1025                          GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues & theEnforcedVertices,
1026                          int &                                           theInvalidEnforcedFlags)
1027 {
1028   std::string tmpStr;
1029   int idx, idxRequired = 0, idxSol = 0;
1030   const int dummyint = 0;
1031   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues::const_iterator vertexIt;
1032   std::vector<double> enfVertexSizes;
1033   const SMDS_MeshElement* elem;
1034   TIDSortedElemSet anElemSet, theKeptEnforcedEdges, theKeptEnforcedTriangles;
1035   SMDS_ElemIteratorPtr nodeIt;
1036   std::vector <const SMDS_MeshNode*> theEnforcedNodeByGhs3dId;
1037   map<const SMDS_MeshNode*,int> anEnforcedNodeToGhs3dIdMap, anExistingEnforcedNodeToGhs3dIdMap;
1038   std::vector< const SMDS_MeshElement* > foundElems;
1039   map<const SMDS_MeshNode*,TopAbs_State> aNodeToTopAbs_StateMap;
1040   int nbFoundElems;
1041   GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap::iterator elemIt;
1042   TIDSortedElemSet::iterator elemSetIt;
1043   bool isOK;
1044   SMESH_Mesh* theMesh = theHelper.GetMesh();
1045   const bool hasGeom = theMesh->HasShapeToMesh();
1046   SMESHUtils::Deleter< SMESH_ElementSearcher > pntCls
1047     ( SMESH_MeshAlgos::GetElementSearcher(*theMesh->GetMeshDS()));
1048   
1049   int nbEnforcedVertices = theEnforcedVertices.size();
1050   theInvalidEnforcedFlags = 0;
1051
1052   // count faces
1053   int nbFaces = theProxyMesh.NbFaces();
1054   int nbNodes;
1055   theFaceByGhs3dId.reserve( nbFaces );
1056   
1057   // groups management
1058   int usedEnforcedNodes = 0;
1059   std::string gn = "";
1060
1061   if ( nbFaces == 0 )
1062     return false;
1063   
1064   idx = MGInput->GmfOpenMesh( theMeshFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1065   if (!idx)
1066     return false;
1067   
1068   /* ========================== FACES ========================== */
1069   /* TRIANGLES ========================== */
1070   SMDS_ElemIteratorPtr eIt =
1071     hasGeom ? theProxyMesh.GetFaces( theHelper.GetSubShape()) : theProxyMesh.GetFaces();
1072   while ( eIt->more() )
1073   {
1074     elem = eIt->next();
1075     anElemSet.insert(elem);
1076     nodeIt = elem->nodesIterator();
1077     nbNodes = elem->NbCornerNodes();
1078     while ( nodeIt->more() && nbNodes--)
1079     {
1080       // find MG-Tetra ID
1081       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1082       int newId = aNodeToGhs3dIdMap.size() + 1; // MG-Tetra ids count from 1
1083       aNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1084     }
1085   }
1086   if ( !anElemSet.empty() &&
1087        (*anElemSet.begin())->IsQuadratic() &&
1088        theProxyMesh.NbProxySubMeshes() > 0 )
1089   {
1090     // add medium nodes of proxy triangles to theHelper (#16843)
1091     for ( elemSetIt = anElemSet.begin(); elemSetIt != anElemSet.end(); ++elemSetIt )
1092       theHelper.AddTLinks( static_cast< const SMDS_MeshFace* >( *elemSetIt ));
1093   }
1094
1095   /* EDGES ========================== */
1096
1097   // Iterate over the enforced edges
1098   for(elemIt = theEnforcedEdges.begin() ; elemIt != theEnforcedEdges.end() ; ++elemIt) {
1099     elem = elemIt->first;
1100     isOK = true;
1101     nodeIt = elem->nodesIterator();
1102     nbNodes = 2;
1103     while ( nodeIt->more() && nbNodes-- ) {
1104       // find MG-Tetra ID
1105       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1106       // Test if point is inside shape to mesh
1107       gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1108       TopAbs_State result = pntCls->GetPointState( myPoint );
1109       if ( result == TopAbs_OUT ) {
1110         isOK = false;
1111         theInvalidEnforcedFlags |= FLAG_BAD_ENF_EDGE;
1112         break;
1113       }
1114       aNodeToTopAbs_StateMap.insert( make_pair( node, result ));
1115     }
1116     if (isOK) {
1117       nodeIt = elem->nodesIterator();
1118       nbNodes = 2;
1119       int newId = -1;
1120       while ( nodeIt->more() && nbNodes-- ) {
1121         // find MG-Tetra ID
1122         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1123         gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1124         nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1125 #ifdef _MY_DEBUG_
1126         std::cout << "Node at "<<node->X()<<", "<<node->Y()<<", "<<node->Z()<<std::endl;
1127         std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
1128 #endif
1129         if (nbFoundElems ==0) {
1130           if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1131             newId = aNodeToGhs3dIdMap.size() + anEnforcedNodeToGhs3dIdMap.size() + 1; // MG-Tetra ids count from 1
1132             anEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1133           }
1134         }
1135         else if (nbFoundElems ==1) {
1136           const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1137           newId = (*aNodeToGhs3dIdMap.find(existingNode)).second;
1138           anExistingEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1139         }
1140         else
1141           isOK = false;
1142 #ifdef _MY_DEBUG_
1143         std::cout << "MG-Tetra node ID: "<<newId<<std::endl;
1144 #endif
1145       }
1146       if (isOK)
1147         theKeptEnforcedEdges.insert(elem);
1148       else
1149         theInvalidEnforcedFlags |= FLAG_BAD_ENF_EDGE;
1150     }
1151   }
1152   
1153   /* ENFORCED TRIANGLES ========================== */
1154   
1155   // Iterate over the enforced triangles
1156   for(elemIt = theEnforcedTriangles.begin() ; elemIt != theEnforcedTriangles.end() ; ++elemIt) {
1157     elem = elemIt->first;
1158     isOK = true;
1159     nodeIt = elem->nodesIterator();
1160     nbNodes = 3;
1161     while ( nodeIt->more() && nbNodes--) {
1162       // find MG-Tetra ID
1163       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1164       // Test if point is inside shape to mesh
1165       gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1166       TopAbs_State result = pntCls->GetPointState( myPoint );
1167       if ( result == TopAbs_OUT ) {
1168         isOK = false;
1169         theInvalidEnforcedFlags |= FLAG_BAD_ENF_TRIA;
1170         break;
1171       }
1172       aNodeToTopAbs_StateMap.insert( make_pair( node, result ));
1173     }
1174     if (isOK) {
1175       nodeIt = elem->nodesIterator();
1176       nbNodes = 3;
1177       int newId = -1;
1178       while ( nodeIt->more() && nbNodes--) {
1179         // find MG-Tetra ID
1180         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1181         gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1182         nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1183 #ifdef _MY_DEBUG_
1184         std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
1185 #endif
1186         if (nbFoundElems ==0) {
1187           if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1188             newId = aNodeToGhs3dIdMap.size() + anEnforcedNodeToGhs3dIdMap.size() + 1; // MG-Tetra ids count from 1
1189             anEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1190           }
1191         }
1192         else if (nbFoundElems ==1) {
1193           const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1194           newId = (*aNodeToGhs3dIdMap.find(existingNode)).second;
1195           anExistingEnforcedNodeToGhs3dIdMap.insert( make_pair( node, newId ));
1196         }
1197         else
1198           isOK = false;
1199 #ifdef _MY_DEBUG_
1200         std::cout << "MG-Tetra node ID: "<<newId<<std::endl;
1201 #endif
1202       }
1203       if (isOK)
1204         theKeptEnforcedTriangles.insert(elem);
1205       else
1206         theInvalidEnforcedFlags |= FLAG_BAD_ENF_TRIA;
1207     }
1208   }
1209   
1210   // put nodes to theNodeByGhs3dId vector
1211 #ifdef _MY_DEBUG_
1212   std::cout << "aNodeToGhs3dIdMap.size(): "<<aNodeToGhs3dIdMap.size()<<std::endl;
1213 #endif
1214   theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
1215   map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
1216   for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
1217   {
1218 //     std::cout << "n2id->first: "<<n2id->first<<std::endl;
1219     theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // MG-Tetra ids count from 1
1220   }
1221
1222   // put nodes to anEnforcedNodeToGhs3dIdMap vector
1223 #ifdef _MY_DEBUG_
1224   std::cout << "anEnforcedNodeToGhs3dIdMap.size(): "<<anEnforcedNodeToGhs3dIdMap.size()<<std::endl;
1225 #endif
1226   theEnforcedNodeByGhs3dId.resize( anEnforcedNodeToGhs3dIdMap.size());
1227   n2id = anEnforcedNodeToGhs3dIdMap.begin();
1228   for ( ; n2id != anEnforcedNodeToGhs3dIdMap.end(); ++ n2id)
1229   {
1230     if (n2id->second > (int)aNodeToGhs3dIdMap.size()) {
1231       theEnforcedNodeByGhs3dId[ n2id->second - aNodeToGhs3dIdMap.size() - 1 ] = n2id->first; // MG-Tetra ids count from 1
1232     }
1233   }
1234   
1235   
1236   /* ========================== NODES ========================== */
1237   vector<const SMDS_MeshNode*> theOrderedNodes, theRequiredNodes;
1238   std::set< std::vector<double> > nodesCoords;
1239   vector<const SMDS_MeshNode*>::const_iterator ghs3dNodeIt = theNodeByGhs3dId.begin();
1240   vector<const SMDS_MeshNode*>::const_iterator after  = theNodeByGhs3dId.end();
1241   
1242   (theNodeByGhs3dId.size() <= 1) ? tmpStr = " node" : " nodes";
1243   std::cout << theNodeByGhs3dId.size() << tmpStr << " from mesh ..." << std::endl;
1244   for ( ; ghs3dNodeIt != after; ++ghs3dNodeIt )
1245   {
1246     const SMDS_MeshNode* node = *ghs3dNodeIt;
1247     std::vector<double> coords;
1248     coords.push_back(node->X());
1249     coords.push_back(node->Y());
1250     coords.push_back(node->Z());
1251     nodesCoords.insert(coords);
1252     theOrderedNodes.push_back(node);
1253   }
1254   
1255   // Iterate over the enforced nodes given by enforced elements
1256   ghs3dNodeIt = theEnforcedNodeByGhs3dId.begin();
1257   after  = theEnforcedNodeByGhs3dId.end();
1258   (theEnforcedNodeByGhs3dId.size() <= 1) ? tmpStr = " node" : " nodes";
1259   std::cout << theEnforcedNodeByGhs3dId.size() << tmpStr << " from enforced elements ..." << std::endl;
1260   for ( ; ghs3dNodeIt != after; ++ghs3dNodeIt )
1261   {
1262     const SMDS_MeshNode* node = *ghs3dNodeIt;
1263     std::vector<double> coords;
1264     coords.push_back(node->X());
1265     coords.push_back(node->Y());
1266     coords.push_back(node->Z());
1267 #ifdef _MY_DEBUG_
1268     std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1269 #endif
1270     
1271     if (nodesCoords.find(coords) != nodesCoords.end()) {
1272       // node already exists in original mesh
1273 #ifdef _MY_DEBUG_
1274       std::cout << " found" << std::endl;
1275 #endif
1276       continue;
1277     }
1278     
1279     if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1280       // node already exists in enforced vertices
1281 #ifdef _MY_DEBUG_
1282       std::cout << " found" << std::endl;
1283 #endif
1284       continue;
1285     }
1286     
1287 //     gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1288 //     nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1289 //     if (nbFoundElems ==0) {
1290 //       std::cout << " not found" << std::endl;
1291 //       if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1292 //         nodesCoords.insert(coords);
1293 //         theOrderedNodes.push_back(node);
1294 //       }
1295 //     }
1296 //     else {
1297 //       std::cout << " found in initial mesh" << std::endl;
1298 //       const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1299 //       nodesCoords.insert(coords);
1300 //       theOrderedNodes.push_back(existingNode);
1301 //     }
1302     
1303 #ifdef _MY_DEBUG_
1304     std::cout << " not found" << std::endl;
1305 #endif
1306     
1307     nodesCoords.insert(coords);
1308     theOrderedNodes.push_back(node);
1309 //     theRequiredNodes.push_back(node);
1310   }
1311   
1312   
1313   // Iterate over the enforced nodes
1314   GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap::const_iterator enfNodeIt;
1315   (theEnforcedNodes.size() <= 1) ? tmpStr = " node" : " nodes";
1316   std::cout << theEnforcedNodes.size() << tmpStr << " from enforced nodes ..." << std::endl;
1317   for(enfNodeIt = theEnforcedNodes.begin() ; enfNodeIt != theEnforcedNodes.end() ; ++enfNodeIt)
1318   {
1319     const SMDS_MeshNode* node = enfNodeIt->first;
1320     std::vector<double> coords;
1321     coords.push_back(node->X());
1322     coords.push_back(node->Y());
1323     coords.push_back(node->Z());
1324 #ifdef _MY_DEBUG_
1325     std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1326 #endif
1327     
1328     // Test if point is inside shape to mesh
1329     gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1330     TopAbs_State result = pntCls->GetPointState( myPoint );
1331     if ( result == TopAbs_OUT ) {
1332 #ifdef _MY_DEBUG_
1333       std::cout << " out of volume" << std::endl;
1334 #endif
1335       theInvalidEnforcedFlags |= FLAG_BAD_ENF_NODE;
1336       continue;
1337     }
1338     
1339     if (nodesCoords.find(coords) != nodesCoords.end()) {
1340 #ifdef _MY_DEBUG_
1341       std::cout << " found in nodesCoords" << std::endl;
1342 #endif
1343 //       theRequiredNodes.push_back(node);
1344       continue;
1345     }
1346
1347     if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1348 #ifdef _MY_DEBUG_
1349       std::cout << " found in theEnforcedVertices" << std::endl;
1350 #endif
1351       continue;
1352     }
1353     
1354 //     nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1355 //     if (nbFoundElems ==0) {
1356 //       std::cout << " not found" << std::endl;
1357 //       if (result == TopAbs_IN) {
1358 //         nodesCoords.insert(coords);
1359 //         theRequiredNodes.push_back(node);
1360 //       }
1361 //     }
1362 //     else {
1363 //       std::cout << " found in initial mesh" << std::endl;
1364 //       const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1365 // //       nodesCoords.insert(coords);
1366 //       theRequiredNodes.push_back(existingNode);
1367 //     }
1368 //     
1369 //     
1370 //     
1371 //     if (pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems) == 0)
1372 //       continue;
1373
1374 //     if ( result != TopAbs_IN )
1375 //       continue;
1376     
1377 #ifdef _MY_DEBUG_
1378     std::cout << " not found" << std::endl;
1379 #endif
1380     nodesCoords.insert(coords);
1381 //     theOrderedNodes.push_back(node);
1382     theRequiredNodes.push_back(node);
1383   }
1384   int requiredNodes = theRequiredNodes.size();
1385   
1386   int solSize = 0;
1387   std::vector<std::vector<double> > ReqVerTab;
1388   if (nbEnforcedVertices) {
1389 //    ReqVerTab.clear();
1390     (nbEnforcedVertices <= 1) ? tmpStr = " node" : " nodes";
1391     std::cout << nbEnforcedVertices << tmpStr << " from enforced vertices ..." << std::endl;
1392     // Iterate over the enforced vertices
1393     for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
1394       double x = vertexIt->first[0];
1395       double y = vertexIt->first[1];
1396       double z = vertexIt->first[2];
1397       // Test if point is inside shape to mesh
1398       gp_Pnt myPoint(x,y,z);
1399       TopAbs_State result = pntCls->GetPointState( myPoint );
1400       if ( result == TopAbs_OUT )
1401       {
1402         std::cout << "Warning: enforced vertex at ( " << x << "," << y << "," << z << " ) is out of the meshed domain!!!" << std::endl;
1403         theInvalidEnforcedFlags |= FLAG_BAD_ENF_VERT;
1404         //continue;
1405       }
1406       std::vector<double> coords;
1407       coords.push_back(x);
1408       coords.push_back(y);
1409       coords.push_back(z);
1410       ReqVerTab.push_back(coords);
1411       enfVertexSizes.push_back(vertexIt->second);
1412       solSize++;
1413     }
1414   }
1415   
1416   
1417   // GmfVertices
1418   std::cout << "Begin writting required nodes in GmfVertices" << std::endl;
1419   std::cout << "Nb vertices: " << theOrderedNodes.size() << std::endl;
1420   MGInput->GmfSetKwd( idx, GmfVertices, theOrderedNodes.size()/*+solSize*/);
1421   for (ghs3dNodeIt = theOrderedNodes.begin();ghs3dNodeIt != theOrderedNodes.end();++ghs3dNodeIt) {
1422     MGInput->GmfSetLin( idx, GmfVertices, (*ghs3dNodeIt)->X(), (*ghs3dNodeIt)->Y(), (*ghs3dNodeIt)->Z(), dummyint);
1423   }
1424
1425   std::cout << "End writting required nodes in GmfVertices" << std::endl;
1426
1427   if (requiredNodes + solSize) {
1428     std::cout << "Begin writting in req and sol file" << std::endl;
1429     aNodeGroupByGhs3dId.resize( requiredNodes + solSize );
1430     idxRequired = MGInput->GmfOpenMesh( theRequiredFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1431     if (!idxRequired) {
1432       return false;
1433     }
1434     idxSol = MGInput->GmfOpenMesh( theSolFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1435     if (!idxSol) {
1436       return false;
1437     }
1438     int TypTab[] = {GmfSca};
1439     double ValTab[] = {0.0};
1440     MGInput->GmfSetKwd( idxRequired, GmfVertices, requiredNodes + solSize);
1441     MGInput->GmfSetKwd( idxSol, GmfSolAtVertices, requiredNodes + solSize, 1, TypTab);
1442 //     int usedEnforcedNodes = 0;
1443 //     std::string gn = "";
1444     for (ghs3dNodeIt = theRequiredNodes.begin();ghs3dNodeIt != theRequiredNodes.end();++ghs3dNodeIt) {
1445       MGInput->GmfSetLin( idxRequired, GmfVertices, (*ghs3dNodeIt)->X(), (*ghs3dNodeIt)->Y(), (*ghs3dNodeIt)->Z(), dummyint);
1446       MGInput->GmfSetLin( idxSol, GmfSolAtVertices, ValTab);
1447       if (theEnforcedNodes.find((*ghs3dNodeIt)) != theEnforcedNodes.end())
1448         gn = theEnforcedNodes.find((*ghs3dNodeIt))->second;
1449       aNodeGroupByGhs3dId[usedEnforcedNodes] = gn;
1450       usedEnforcedNodes++;
1451     }
1452
1453     for (int i=0;i<solSize;i++) {
1454       std::cout << ReqVerTab[i][0] <<" "<< ReqVerTab[i][1] << " "<< ReqVerTab[i][2] << std::endl;
1455 #ifdef _MY_DEBUG_
1456       std::cout << "enfVertexSizes.at("<<i<<"): " << enfVertexSizes.at(i) << std::endl;
1457 #endif
1458       double solTab[] = {enfVertexSizes.at(i)};
1459       MGInput->GmfSetLin( idxRequired, GmfVertices, ReqVerTab[i][0], ReqVerTab[i][1], ReqVerTab[i][2], dummyint);
1460       MGInput->GmfSetLin( idxSol, GmfSolAtVertices, solTab);
1461       aNodeGroupByGhs3dId[usedEnforcedNodes] = enfVerticesWithGroup.find(ReqVerTab[i])->second;
1462 #ifdef _MY_DEBUG_
1463       std::cout << "aNodeGroupByGhs3dId["<<usedEnforcedNodes<<"] = \""<<aNodeGroupByGhs3dId[usedEnforcedNodes]<<"\""<<std::endl;
1464 #endif
1465       usedEnforcedNodes++;
1466     }
1467     std::cout << "End writting in req and sol file" << std::endl;
1468   }
1469
1470   int nedge[2], ntri[3];
1471     
1472   // GmfEdges
1473   int usedEnforcedEdges = 0;
1474   if (theKeptEnforcedEdges.size()) {
1475     anEdgeGroupByGhs3dId.resize( theKeptEnforcedEdges.size() );
1476 //    idxRequired = MGInput->GmfOpenMesh( theRequiredFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1477 //    if (!idxRequired)
1478 //      return false;
1479     MGInput->GmfSetKwd( idx, GmfEdges, theKeptEnforcedEdges.size());
1480 //    MGInput->GmfSetKwd( idxRequired, GmfEdges, theKeptEnforcedEdges.size());
1481     for(elemSetIt = theKeptEnforcedEdges.begin() ; elemSetIt != theKeptEnforcedEdges.end() ; ++elemSetIt) {
1482       elem = (*elemSetIt);
1483       nodeIt = elem->nodesIterator();
1484       int index=0;
1485       while ( nodeIt->more() ) {
1486         // find MG-Tetra ID
1487         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1488         map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToGhs3dIdMap.find(node);
1489         if (it == anEnforcedNodeToGhs3dIdMap.end()) {
1490           it = anExistingEnforcedNodeToGhs3dIdMap.find(node);
1491           if (it == anEnforcedNodeToGhs3dIdMap.end())
1492             throw "Node not found";
1493         }
1494         nedge[index] = it->second;
1495         index++;
1496       }
1497       MGInput->GmfSetLin( idx, GmfEdges, nedge[0], nedge[1], dummyint);
1498       anEdgeGroupByGhs3dId[usedEnforcedEdges] = theEnforcedEdges.find(elem)->second;
1499 //      MGInput->GmfSetLin( idxRequired, GmfEdges, nedge[0], nedge[1], dummyint);
1500       usedEnforcedEdges++;
1501     }
1502   }
1503
1504
1505   if (usedEnforcedEdges) {
1506     MGInput->GmfSetKwd( idx, GmfRequiredEdges, usedEnforcedEdges);
1507     for (int enfID=1;enfID<=usedEnforcedEdges;enfID++) {
1508       MGInput->GmfSetLin( idx, GmfRequiredEdges, enfID);
1509     }
1510   }
1511
1512   // GmfTriangles
1513   int usedEnforcedTriangles = 0;
1514   if (anElemSet.size()+theKeptEnforcedTriangles.size()) {
1515     aFaceGroupByGhs3dId.resize( anElemSet.size()+theKeptEnforcedTriangles.size() );
1516     MGInput->GmfSetKwd( idx, GmfTriangles, anElemSet.size()+theKeptEnforcedTriangles.size());
1517     int k=0;
1518     for(elemSetIt = anElemSet.begin() ; elemSetIt != anElemSet.end() ; ++elemSetIt,++k) {
1519       elem = (*elemSetIt);
1520       theFaceByGhs3dId.push_back( elem );
1521       nodeIt = elem->nodesIterator();
1522       int index=0;
1523       for ( int j = 0; j < 3; ++j ) {
1524         // find MG-Tetra ID
1525         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1526         map< const SMDS_MeshNode*,int >::iterator it = aNodeToGhs3dIdMap.find(node);
1527         if (it == aNodeToGhs3dIdMap.end())
1528           throw "Node not found";
1529         ntri[index] = it->second;
1530         index++;
1531       }
1532       MGInput->GmfSetLin( idx, GmfTriangles, ntri[0], ntri[1], ntri[2], dummyint);
1533       aFaceGroupByGhs3dId[k] = "";
1534     }
1535     if ( !theHelper.GetMesh()->HasShapeToMesh() )
1536       SMESHUtils::FreeVector( theFaceByGhs3dId );
1537     if (theKeptEnforcedTriangles.size()) {
1538       for(elemSetIt = theKeptEnforcedTriangles.begin() ; elemSetIt != theKeptEnforcedTriangles.end() ; ++elemSetIt,++k) {
1539         elem = (*elemSetIt);
1540         nodeIt = elem->nodesIterator();
1541         int index=0;
1542         for ( int j = 0; j < 3; ++j ) {
1543           // find MG-Tetra ID
1544           const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1545           map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToGhs3dIdMap.find(node);
1546           if (it == anEnforcedNodeToGhs3dIdMap.end()) {
1547             it = anExistingEnforcedNodeToGhs3dIdMap.find(node);
1548             if (it == anEnforcedNodeToGhs3dIdMap.end())
1549               throw "Node not found";
1550           }
1551           ntri[index] = it->second;
1552           index++;
1553         }
1554         MGInput->GmfSetLin( idx, GmfTriangles, ntri[0], ntri[1], ntri[2], dummyint);
1555         aFaceGroupByGhs3dId[k] = theEnforcedTriangles.find(elem)->second;
1556         usedEnforcedTriangles++;
1557       }
1558     }
1559   }
1560
1561   
1562   if (usedEnforcedTriangles) {
1563     MGInput->GmfSetKwd( idx, GmfRequiredTriangles, usedEnforcedTriangles);
1564     for (int enfID=1;enfID<=usedEnforcedTriangles;enfID++)
1565       MGInput->GmfSetLin( idx, GmfRequiredTriangles, anElemSet.size()+enfID);
1566   }
1567
1568   MGInput->GmfCloseMesh(idx);
1569   if (idxRequired)
1570     MGInput->GmfCloseMesh(idxRequired);
1571   if (idxSol)
1572     MGInput->GmfCloseMesh(idxSol);
1573
1574   return true;
1575 }
1576
1577 //=============================================================================
1578 /*!
1579  *Here we are going to use the MG-Tetra mesher with geometry
1580  */
1581 //=============================================================================
1582
1583 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
1584                                 const TopoDS_Shape& theShape)
1585 {
1586   bool Ok(false);
1587   TopExp_Explorer expBox ( theShape, TopAbs_SOLID );
1588
1589   // a unique working file name
1590   // to avoid access to the same files by eg different users
1591   _genericName = GHS3DPlugin_Hypothesis::GetFileName(_hyp);
1592   TCollection_AsciiString aGenericName((char*) _genericName.c_str() );
1593   TCollection_AsciiString aGenericNameRequired = aGenericName + "_required";
1594
1595   TCollection_AsciiString aLogFileName    = aGenericName + ".log";    // log
1596   TCollection_AsciiString aResultFileName;
1597
1598   TCollection_AsciiString aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1599   aGMFFileName              = aGenericName + ".mesh"; // GMF mesh file
1600   aResultFileName           = aGenericName + "Vol.mesh"; // GMF mesh file
1601   aResSolFileName           = aGenericName + "Vol.sol"; // GMF mesh file
1602   aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1603   aSolFileName              = aGenericNameRequired + ".sol"; // GMF solution file
1604   
1605   std::map <int,int> aNodeId2NodeIndexMap, aSmdsToGhs3dIdMap, anEnforcedNodeIdToGhs3dIdMap;
1606   std::map <int, int> nodeID2nodeIndexMap;
1607   std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1608   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues coordsSizeMap = GHS3DPlugin_Hypothesis::GetEnforcedVerticesCoordsSize(_hyp);
1609   GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = GHS3DPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1610   GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = GHS3DPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1611   GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = GHS3DPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1612   GHS3DPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = GHS3DPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1613
1614   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList enfVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1615   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1616   std::vector<double> coords;
1617
1618   for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1619   {
1620     GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertex* enfVertex = (*enfVerIt);
1621     if (enfVertex->coords.size()) {
1622       coordsSizeMap.insert(make_pair(enfVertex->coords,enfVertex->size));
1623       enfVerticesWithGroup.insert(make_pair(enfVertex->coords,enfVertex->groupName));
1624     }
1625     else {
1626       TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
1627       for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1628         coords.clear();
1629         if (it.Value().ShapeType() == TopAbs_VERTEX){
1630           gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1631           coords.push_back(aPnt.X());
1632           coords.push_back(aPnt.Y());
1633           coords.push_back(aPnt.Z());
1634           if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1635             coordsSizeMap.insert(make_pair(coords,enfVertex->size));
1636             enfVerticesWithGroup.insert(make_pair(coords,enfVertex->groupName));
1637           }
1638         }
1639       }
1640     }
1641   }
1642   int nbEnforcedVertices = coordsSizeMap.size();
1643   int nbEnforcedNodes = enforcedNodes.size();
1644
1645   std::string tmpStr;
1646   (nbEnforcedNodes <= 1) ? tmpStr = "node" : "nodes";
1647   std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1648   (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : "vertices";
1649   std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1650
1651   SMESH_MesherHelper helper( theMesh );
1652   helper.SetSubShape( theShape );
1653
1654   std::vector <const SMDS_MeshNode*> aNodeByGhs3dId, anEnforcedNodeByGhs3dId;
1655   std::vector <const SMDS_MeshElement*> aFaceByGhs3dId;
1656   std::map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
1657   std::vector<std::string> aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId;
1658
1659   MG_Tetra_API mgTetra( _computeCanceled, _progress );
1660
1661   _isLibUsed = mgTetra.IsLibrary();
1662   if ( theMesh.NbQuadrangles() > 0 )
1663     _progressAdvance /= 10;
1664   if ( _viscousLayersHyp )
1665     _progressAdvance /= 10;
1666
1667   // proxyMesh must live till readGMFFile() as a proxy face can be used by
1668   // MG-Tetra for domain indication
1669   SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1670
1671   // make prisms on quadrangles and viscous layers
1672   if ( theMesh.NbQuadrangles() > 0 || _viscousLayersHyp )
1673   {
1674     vector<SMESH_ProxyMesh::Ptr> components;
1675     for (expBox.ReInit(); expBox.More(); expBox.Next())
1676     {
1677       if ( _viscousLayersHyp )
1678       {
1679         proxyMesh = _viscousLayersHyp->Compute( theMesh, expBox.Current() );
1680         if ( !proxyMesh )
1681           return false;
1682         if ( theMesh.NbQuadrangles() == 0 )
1683           components.push_back( proxyMesh );
1684       }
1685       if ( theMesh.NbQuadrangles() > 0 )
1686       {
1687         StdMeshers_QuadToTriaAdaptor* q2t = new StdMeshers_QuadToTriaAdaptor;
1688         Ok = q2t->Compute( theMesh, expBox.Current(), proxyMesh.get() );
1689         components.push_back( SMESH_ProxyMesh::Ptr( q2t ));
1690         if ( !Ok )
1691           return false;
1692       }
1693     }
1694     proxyMesh.reset( new SMESH_ProxyMesh( components ));
1695   }
1696   // build viscous layers
1697   // else if ( _viscousLayersHyp )
1698   // {
1699   //   proxyMesh = _viscousLayersHyp->Compute( theMesh, theShape );
1700   //   if ( !proxyMesh )
1701   //     return false;
1702   // }
1703
1704   int anInvalidEnforcedFlags = 0;
1705   Ok = writeGMFFile(&mgTetra,
1706                     aGMFFileName.ToCString(),
1707                     aRequiredVerticesFileName.ToCString(),
1708                     aSolFileName.ToCString(),
1709                     *proxyMesh, helper,
1710                     aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
1711                     aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
1712                     enforcedNodes, enforcedEdges, enforcedTriangles,
1713                     enfVerticesWithGroup, coordsSizeMap, anInvalidEnforcedFlags);
1714
1715   // Write aSmdsToGhs3dIdMap to temp file
1716   TCollection_AsciiString aSmdsToGhs3dIdMapFileName;
1717   aSmdsToGhs3dIdMapFileName = aGenericName + ".ids";  // ids relation
1718   ofstream aIdsFile  ( aSmdsToGhs3dIdMapFileName.ToCString()  , ios::out);
1719   Ok = aIdsFile.rdbuf()->is_open();
1720   if (!Ok) {
1721     INFOS( "Can't write into " << aSmdsToGhs3dIdMapFileName);
1722     //return error(SMESH_Comment("Can't write into ") << aSmdsToGhs3dIdMapFileName);
1723   }
1724   else
1725   {
1726     INFOS( "Writing ids relation into " << aSmdsToGhs3dIdMapFileName);
1727     aIdsFile << "Smds MG-Tetra" << std::endl;
1728     map <int,int>::const_iterator myit;
1729     for (myit=aSmdsToGhs3dIdMap.begin() ; myit != aSmdsToGhs3dIdMap.end() ; ++myit) {
1730       aIdsFile << myit->first << " " << myit->second << std::endl;
1731     }
1732   }
1733   aIdsFile.close();
1734   if ( ! Ok ) {
1735     if ( !_keepFiles ) {
1736       removeFile( aGMFFileName );
1737       removeFile( aRequiredVerticesFileName );
1738       removeFile( aSolFileName );
1739       removeFile( aSmdsToGhs3dIdMapFileName );
1740     }
1741     return error(COMPERR_BAD_INPUT_MESH);
1742   }
1743   removeFile( aResultFileName ); // needed for boundary recovery module usage
1744
1745   // -----------------
1746   // run MG-Tetra mesher
1747   // -----------------
1748
1749   TCollection_AsciiString cmd = GHS3DPlugin_Hypothesis::CommandToRun( _hyp, true, mgTetra.IsExecutable() ).c_str();
1750
1751   if ( mgTetra.IsExecutable() )
1752   {
1753     cmd += TCollection_AsciiString(" --in ") + aGMFFileName;
1754     if ( nbEnforcedVertices + nbEnforcedNodes)
1755       cmd += TCollection_AsciiString(" --required_vertices ") + aGenericNameRequired;
1756     cmd += TCollection_AsciiString(" --out ") + aResultFileName;
1757   }
1758   if ( !_logInStandardOutput )
1759   {
1760     mgTetra.SetLogFile( aLogFileName.ToCString() );
1761     cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1762   }
1763
1764   BRIEF_INFOS("")
1765   BRIEF_INFOS("MG-Tetra execution...")
1766   BRIEF_INFOS(cmd)
1767
1768   _computeCanceled = false;
1769
1770   std::string errStr;
1771   Ok = mgTetra.Compute( cmd.ToCString(), errStr ); // run
1772
1773   if ( _logInStandardOutput && mgTetra.IsLibrary() ) {
1774     BRIEF_INFOS("");
1775     BRIEF_INFOS(mgTetra.GetLog());
1776     BRIEF_INFOS("")
1777   }
1778   if ( Ok ) {
1779     BRIEF_INFOS("");
1780     BRIEF_INFOS("End of MG-Tetra execution !");
1781     BRIEF_INFOS("")
1782   }
1783
1784   // --------------
1785   // read a result
1786   // --------------
1787
1788   GHS3DPlugin_Hypothesis::TSetStrings groupsToRemove = GHS3DPlugin_Hypothesis::GetGroupsToRemove(_hyp);
1789   bool toMeshHoles =
1790     _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1791   const bool toMakeGroupsOfDomains = GHS3DPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
1792
1793   helper.IsQuadraticSubMesh( theShape );
1794   helper.SetElementsOnShape( false );
1795
1796   Ok = readGMFFile(&mgTetra,
1797                    aResultFileName.ToCString(),
1798                    this,
1799                    &helper, aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
1800                    aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
1801                    groupsToRemove, toMakeGroupsOfDomains, toMeshHoles);
1802
1803   removeEmptyGroupsOfDomains( helper.GetMesh(), /*notEmptyAsWell =*/ !toMakeGroupsOfDomains );
1804
1805
1806
1807   // ---------------------
1808   // remove working files
1809   // ---------------------
1810
1811   if ( Ok )
1812   {
1813     if ( anInvalidEnforcedFlags )
1814       error( COMPERR_WARNING, flagsToErrorStr( anInvalidEnforcedFlags ));
1815     if ( _removeLogOnSuccess )
1816       removeFile( aLogFileName );
1817     // if ( _hyp && _hyp->GetToMakeGroupsOfDomains() )
1818     //   error( COMPERR_WARNING, "'toMakeGroupsOfDomains' is ignored since the mesh is on shape" );
1819   }
1820   else if ( mgTetra.HasLog() )
1821   {
1822     if( _computeCanceled )
1823       error( "interruption initiated by user" );
1824     else
1825     {
1826       // get problem description from the log file
1827       _Ghs2smdsConvertor conv( aNodeByGhs3dId, proxyMesh );
1828       error( getErrorDescription( _logInStandardOutput ? 0 : aLogFileName.ToCString(),
1829                                   mgTetra.GetLog(), conv ));
1830     }
1831   }
1832   else if ( !errStr.empty() )
1833   {
1834     // the log file is empty
1835     removeFile( aLogFileName );
1836     INFOS( "MG-Tetra Error, " << errStr);
1837     error(COMPERR_ALGO_FAILED, errStr);
1838   }
1839
1840   if ( !_keepFiles ) {
1841     if (! Ok && _computeCanceled )
1842       removeFile( aLogFileName );
1843     removeFile( aGMFFileName );
1844     removeFile( aRequiredVerticesFileName );
1845     removeFile( aSolFileName );
1846     removeFile( aResSolFileName );
1847     removeFile( aResultFileName );
1848     removeFile( aSmdsToGhs3dIdMapFileName );
1849   }
1850   if ( mgTetra.IsExecutable() )
1851   {
1852     std::cout << "<" << aResultFileName.ToCString() << "> MG-Tetra output file ";
1853     if ( !Ok )
1854       std::cout << "not ";
1855     std::cout << "treated !" << std::endl;
1856     std::cout << std::endl;
1857   }
1858   else
1859   {
1860     std::cout << "MG-Tetra " << ( Ok ? "succeeded" : "failed") << std::endl;
1861   }
1862   return Ok;
1863 }
1864
1865 //=============================================================================
1866 /*!
1867  *Here we are going to use the MG-Tetra mesher w/o geometry
1868  */
1869 //=============================================================================
1870 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
1871                                 SMESH_MesherHelper* theHelper)
1872 {
1873   theHelper->IsQuadraticSubMesh( theHelper->GetSubShape() );
1874
1875   // a unique working file name
1876   // to avoid access to the same files by eg different users
1877   _genericName = GHS3DPlugin_Hypothesis::GetFileName(_hyp);
1878   TCollection_AsciiString aGenericName((char*) _genericName.c_str() );
1879   TCollection_AsciiString aGenericNameRequired = aGenericName + "_required";
1880
1881   TCollection_AsciiString aLogFileName    = aGenericName + ".log";    // log
1882   TCollection_AsciiString aResultFileName;
1883   bool Ok;
1884
1885   TCollection_AsciiString aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1886   aGMFFileName              = aGenericName + ".mesh"; // GMF mesh file
1887   aResultFileName           = aGenericName + "Vol.mesh"; // GMF mesh file
1888   aResSolFileName           = aGenericName + "Vol.sol"; // GMF mesh file
1889   aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1890   aSolFileName              = aGenericNameRequired + ".sol"; // GMF solution file
1891
1892   std::map <int, int> nodeID2nodeIndexMap;
1893   std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1894   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexCoordsValues coordsSizeMap;
1895   TopoDS_Shape GeomShape;
1896   std::vector<double> coords;
1897   gp_Pnt aPnt;
1898   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertex* enfVertex;
1899
1900   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList enfVertices = GHS3DPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1901   GHS3DPlugin_Hypothesis::TGHS3DEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1902
1903   for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1904   {
1905     enfVertex = (*enfVerIt);
1906     if (enfVertex->coords.size()) {
1907       coordsSizeMap.insert(make_pair(enfVertex->coords,enfVertex->size));
1908       enfVerticesWithGroup.insert(make_pair(enfVertex->coords,enfVertex->groupName));
1909     }
1910     else {
1911       GeomShape = entryToShape(enfVertex->geomEntry);
1912       for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1913         coords.clear();
1914         if (it.Value().ShapeType() == TopAbs_VERTEX){
1915           aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1916           coords.push_back(aPnt.X());
1917           coords.push_back(aPnt.Y());
1918           coords.push_back(aPnt.Z());
1919           if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1920             coordsSizeMap.insert(make_pair(coords,enfVertex->size));
1921             enfVerticesWithGroup.insert(make_pair(coords,enfVertex->groupName));
1922           }
1923         }
1924       }
1925     }
1926   }
1927
1928   GHS3DPlugin_Hypothesis::TIDSortedNodeGroupMap     enforcedNodes = GHS3DPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1929   GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap     enforcedEdges = GHS3DPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1930   GHS3DPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = GHS3DPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1931   GHS3DPlugin_Hypothesis::TID2SizeMap             nodeIDToSizeMap = GHS3DPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1932
1933   std::string tmpStr;
1934
1935   int nbEnforcedVertices = coordsSizeMap.size();
1936   int nbEnforcedNodes = enforcedNodes.size();
1937   (nbEnforcedNodes <= 1) ? tmpStr = "node" : tmpStr = "nodes";
1938   std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1939   (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : tmpStr = "vertices";
1940   std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1941
1942   std::vector <const SMDS_MeshNode*> aNodeByGhs3dId, anEnforcedNodeByGhs3dId;
1943   std::vector <const SMDS_MeshElement*> aFaceByGhs3dId;
1944   std::map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
1945   std::vector<std::string> aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId;
1946
1947
1948   MG_Tetra_API mgTetra( _computeCanceled, _progress );
1949
1950   _isLibUsed = mgTetra.IsLibrary();
1951   if ( theMesh.NbQuadrangles() > 0 )
1952     _progressAdvance /= 10;
1953
1954   // proxyMesh must live till readGMFFile() as a proxy face can be used by
1955   // MG-Tetra for domain indication
1956   SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1957   if ( theMesh.NbQuadrangles() > 0 )
1958   {
1959     StdMeshers_QuadToTriaAdaptor* aQuad2Trias = new StdMeshers_QuadToTriaAdaptor;
1960     Ok = aQuad2Trias->Compute( theMesh );
1961     proxyMesh.reset( aQuad2Trias );
1962     if ( !Ok )
1963       return false;
1964   }
1965
1966   int anInvalidEnforcedFlags = 0;
1967   Ok = writeGMFFile(&mgTetra,
1968                     aGMFFileName.ToCString(), aRequiredVerticesFileName.ToCString(), aSolFileName.ToCString(),
1969                     *proxyMesh, *theHelper,
1970                     aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
1971                     aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
1972                     enforcedNodes, enforcedEdges, enforcedTriangles,
1973                     enfVerticesWithGroup, coordsSizeMap, anInvalidEnforcedFlags);
1974
1975   // -----------------
1976   // run MG-Tetra mesher
1977   // -----------------
1978
1979   TCollection_AsciiString cmd = GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false, mgTetra.IsExecutable() ).c_str();
1980
1981   if ( mgTetra.IsExecutable() )
1982   {
1983     cmd += TCollection_AsciiString(" --in ") + aGMFFileName;
1984     if ( nbEnforcedVertices + nbEnforcedNodes)
1985       cmd += TCollection_AsciiString(" --required_vertices ") + aGenericNameRequired;
1986     cmd += TCollection_AsciiString(" --out ") + aResultFileName;
1987   }
1988   if ( !_logInStandardOutput )
1989   {
1990     mgTetra.SetLogFile( aLogFileName.ToCString() );
1991     cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1992   }
1993
1994   BRIEF_INFOS("")
1995   BRIEF_INFOS("MG-Tetra execution...")
1996   BRIEF_INFOS(cmd)
1997
1998   _computeCanceled = false;
1999
2000   std::string errStr;
2001   Ok = mgTetra.Compute( cmd.ToCString(), errStr ); // run
2002
2003   if ( _logInStandardOutput && mgTetra.IsLibrary() ) {
2004     BRIEF_INFOS("");
2005     BRIEF_INFOS(mgTetra.GetLog());
2006     BRIEF_INFOS("")
2007   }
2008   if ( Ok ) {
2009     BRIEF_INFOS("");
2010     BRIEF_INFOS("End of MG-Tetra execution !");
2011     BRIEF_INFOS("")
2012   }
2013
2014   // --------------
2015   // read a result
2016   // --------------
2017   GHS3DPlugin_Hypothesis::TSetStrings groupsToRemove = GHS3DPlugin_Hypothesis::GetGroupsToRemove(_hyp);
2018   const bool toMakeGroupsOfDomains = GHS3DPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
2019
2020   Ok = Ok && readGMFFile(&mgTetra,
2021                          aResultFileName.ToCString(),
2022                          this,
2023                          theHelper, aNodeByGhs3dId, aFaceByGhs3dId, aNodeToGhs3dIdMap,
2024                          aNodeGroupByGhs3dId, anEdgeGroupByGhs3dId, aFaceGroupByGhs3dId,
2025                          groupsToRemove, toMakeGroupsOfDomains);
2026
2027   updateMeshGroups(theHelper->GetMesh(), groupsToRemove);
2028   removeEmptyGroupsOfDomains( theHelper->GetMesh(), /*notEmptyAsWell =*/ !toMakeGroupsOfDomains );
2029
2030   if ( Ok ) {
2031     GHS3DPlugin_Hypothesis* that = (GHS3DPlugin_Hypothesis*)this->_hyp;
2032     if (that)
2033       that->ClearGroupsToRemove();
2034   }
2035   // ---------------------
2036   // remove working files
2037   // ---------------------
2038
2039   if ( Ok )
2040   {
2041     if ( anInvalidEnforcedFlags )
2042       error( COMPERR_WARNING, flagsToErrorStr( anInvalidEnforcedFlags ));
2043     if ( _removeLogOnSuccess )
2044       removeFile( aLogFileName );
2045
2046     //if ( !toMakeGroupsOfDomains && _hyp && _hyp->GetToMakeGroupsOfDomains() )
2047     //error( COMPERR_WARNING, "'toMakeGroupsOfDomains' is ignored since 'toMeshHoles' is OFF." );
2048   }
2049   else if ( mgTetra.HasLog() )
2050   {
2051     if( _computeCanceled )
2052       error( "interruption initiated by user" );
2053     else
2054     {
2055       // get problem description from the log file
2056       _Ghs2smdsConvertor conv( aNodeByGhs3dId, proxyMesh );
2057       error( getErrorDescription( _logInStandardOutput ? 0 : aLogFileName.ToCString(),
2058                                   mgTetra.GetLog(), conv ));
2059     }
2060   }
2061   else {
2062     // the log file is empty
2063     removeFile( aLogFileName );
2064     INFOS( "MG-Tetra Error, " << errStr);
2065     error(COMPERR_ALGO_FAILED, errStr);
2066   }
2067
2068   if ( !_keepFiles )
2069   {
2070     if (! Ok && _computeCanceled)
2071       removeFile( aLogFileName );
2072     removeFile( aGMFFileName );
2073     removeFile( aResultFileName );
2074     removeFile( aRequiredVerticesFileName );
2075     removeFile( aSolFileName );
2076     removeFile( aResSolFileName );
2077   }
2078   return Ok;
2079 }
2080
2081 void GHS3DPlugin_GHS3D::CancelCompute()
2082 {
2083   _computeCanceled = true;
2084 #if !defined WIN32 && !defined __APPLE__
2085   std::string cmd = "ps xo pid,args | grep " + _genericName;
2086   //cmd += " | grep -e \"^ *[0-9]\\+ \\+" + GHS3DPlugin_Hypothesis::GetExeName() + "\"";
2087   cmd += " | awk '{print $1}' | xargs kill -9 > /dev/null 2>&1";
2088   system( cmd.c_str() );
2089 #endif
2090 }
2091
2092 //================================================================================
2093 /*!
2094  * \brief Provide human readable text by error code reported by MG-Tetra
2095  */
2096 //================================================================================
2097
2098 static const char* translateError(const int errNum)
2099 {
2100   switch ( errNum ) {
2101   case 0:
2102     return "The surface mesh includes a face of type other than edge, "
2103       "triangle or quadrilateral. This face type is not supported.";
2104   case 1:
2105     return "Not enough memory for the face table.";
2106   case 2:
2107     return "Not enough memory.";
2108   case 3:
2109     return "Not enough memory.";
2110   case 4:
2111     return "Face is ignored.";
2112   case 5:
2113     return "End of file. Some data are missing in the file.";
2114   case 6:
2115     return "Read error on the file. There are wrong data in the file.";
2116   case 7:
2117     return "the metric file is inadequate (dimension other than 3).";
2118   case 8:
2119     return "the metric file is inadequate (values not per vertices).";
2120   case 9:
2121     return "the metric file contains more than one field.";
2122   case 10:
2123     return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
2124       "value of number of mesh vertices in the \".noboite\" file.";
2125   case 12:
2126     return "Too many sub-domains.";
2127   case 13:
2128     return "the number of vertices is negative or null.";
2129   case 14:
2130     return "the number of faces is negative or null.";
2131   case 15:
2132     return "A face has a null vertex.";
2133   case 22:
2134     return "incompatible data.";
2135   case 131:
2136     return "the number of vertices is negative or null.";
2137   case 132:
2138     return "the number of vertices is negative or null (in the \".mesh\" file).";
2139   case 133:
2140     return "the number of faces is negative or null.";
2141   case 1000:
2142     return "A face appears more than once in the input surface mesh.";
2143   case 1001:
2144     return "An edge appears more than once in the input surface mesh.";
2145   case 1002:
2146     return "A face has a vertex negative or null.";
2147   case 1003:
2148     return "NOT ENOUGH MEMORY.";
2149   case 2000:
2150     return "Not enough available memory.";
2151   case 2002:
2152     return "Some initial points cannot be inserted. The surface mesh is probably very bad "
2153       "in terms of quality or the input list of points is wrong.";
2154   case 2003:
2155     return "Some vertices are too close to one another or coincident.";
2156   case 2004:
2157     return "Some vertices are too close to one another or coincident.";
2158   case 2012:
2159     return "A vertex cannot be inserted.";
2160   case 2014:
2161     return "There are at least two points considered as coincident.";
2162   case 2103:
2163     return "Some vertices are too close to one another or coincident.";
2164   case 3000:
2165     return "The surface mesh regeneration step has failed.";
2166   case 3009:
2167     return "Constrained edge cannot be enforced.";
2168   case 3019:
2169     return "Constrained face cannot be enforced.";
2170   case 3029:
2171     return "Missing faces.";
2172   case 3100:
2173     return "No guess to start the definition of the connected component(s).";
2174   case 3101:
2175     return "The surface mesh includes at least one hole. The domain is not well defined.";
2176   case 3102:
2177     return "Impossible to define a component.";
2178   case 3103:
2179     return "The surface edge intersects another surface edge.";
2180   case 3104:
2181     return "The surface edge intersects the surface face.";
2182   case 3105:
2183     return "One boundary point lies within a surface face.";
2184   case 3106:
2185     return "One surface edge intersects a surface face.";
2186   case 3107:
2187     return "One boundary point lies within a surface edge.";
2188   case 3108:
2189     return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
2190       "to too many swaps.";
2191   case 3109:
2192     return "Edge is unique (i.e., bounds a hole in the surface).";
2193   case 3122:
2194     return "Presumably, the surface mesh is not compatible with the domain being processed.";
2195   case 3123:
2196     return "Too many components, too many sub-domain.";
2197   case 3209:
2198     return "The surface mesh includes at least one hole. "
2199       "Therefore there is no domain properly defined.";
2200   case 3300:
2201     return "Statistics.";
2202   case 3400:
2203     return "Statistics.";
2204   case 3500:
2205     return "Warning, it is dramatically tedious to enforce the boundary items.";
2206   case 4000:
2207     return "Not enough memory at this time, nevertheless, the program continues. "
2208       "The expected mesh will be correct but not really as large as required.";
2209   case 4002:
2210     return "see above error code, resulting quality may be poor.";
2211   case 4003:
2212     return "Not enough memory at this time, nevertheless, the program continues (warning).";
2213   case 8000:
2214     return "Unknown face type.";
2215   case 8005:
2216   case 8006:
2217     return "End of file. Some data are missing in the file.";
2218   case 9000:
2219     return "A too small volume element is detected.";
2220   case 9001:
2221     return "There exists at least a null or negative volume element.";
2222   case 9002:
2223     return "There exist null or negative volume elements.";
2224   case 9003:
2225     return "A too small volume element is detected. A face is considered being degenerated.";
2226   case 9100:
2227     return "Some element is suspected to be very bad shaped or wrong.";
2228   case 9102:
2229     return "A too bad quality face is detected. This face is considered degenerated.";
2230   case 9112:
2231     return "A too bad quality face is detected. This face is degenerated.";
2232   case 9122:
2233     return "Presumably, the surface mesh is not compatible with the domain being processed.";
2234   case 9999:
2235     return "Abnormal error occured, contact hotline.";
2236   case 23600:
2237     return "Not enough memory for the face table.";
2238   case 23601:
2239     return "The algorithm cannot run further. "
2240       "The surface mesh is probably very bad in terms of quality.";
2241   case 23602:
2242     return "Bad vertex number.";
2243   case 1001200:
2244     return "Cannot close mesh file NomFil.";
2245   case 1002010:
2246     return "There are wrong data.";
2247   case 1002120:
2248     return "The number of faces is negative or null.";
2249   case 1002170:
2250     return "The number of vertices is negative or null in the '.sol' file.";
2251   case 1002190:
2252     return "The number of tetrahedra is negative or null.";
2253   case 1002210:
2254     return "The number of vertices is negative or null.";
2255   case 1002211:
2256     return "A face has a vertex negative or null.";
2257   case 1002270:
2258     return "The field is not a size in file NomFil.";
2259   case 1002280:
2260     return "A count is wrong in the enclosing box in the .boite.mesh input "
2261       "file (option '--read_boite').";
2262   case 1002290:
2263     return "A tetrahedron has a vertex with a negative number.";
2264   case 1002300:
2265     return "the 'MeshVersionFormatted' is not 1 or 2 in the '.mesh' file or the '.sol'.";
2266  case 1002370:
2267    return "The number of values in the '.sol' (metric file) is incompatible with "
2268      "the expected value of number of mesh vertices in the '.mesh' file.";
2269   case 1003000:
2270     return "Not enough memory.";
2271   case 1003020:
2272     return "Not enough memory for the face table.";
2273   case 1003050:
2274     return "Insufficient memory ressources detected due to a bad quality "
2275       "surface mesh leading to too many swaps.";
2276   case 1005010:
2277     return "The surface coordinates of a vertex are differing from the "
2278       "volume coordinates, probably due to a precision problem.";
2279   case 1005050:
2280     return "Invalid dimension. Dimension 3 expected.";
2281   case 1005100:
2282     return "A point has a tag 0. This point is probably outside the domain which has been meshed.";
2283   case 1005103:
2284     return "The vertices of an element are too close to one another or coincident.";
2285   case 1005104:
2286     return "There are at least two points whose distance is very small, and considered as coincident.";
2287   case 1005105:
2288     return "Two vertices are too close to one another or coincident.";
2289   case 1005106:
2290     return "A vertex cannot be inserted.";
2291   case 1005107:
2292     return "Two vertices are too close to one another or coincident. Note : When "
2293       "this error occurs during the overconstrained processing phase, this is only "
2294       "a warning which means that it is difficult to break some overconstrained facets.";
2295   case 1005110:
2296     return "Two surface edges are intersecting.";
2297   case 1005120:
2298     return "A surface edge intersects a surface face.";
2299   case 1005150:
2300     return "A boundary point lies within a surface face.";
2301   case 1005160:
2302     return "A boundary point lies within a surface edge.";
2303   case 1005200:
2304     return "A surface mesh appears more than once in the input surface mesh.";
2305   case 1005210:
2306     return "An edge appears more than once in the input surface mesh.";
2307   case 1005225:
2308     return "Surface with unvalid triangles.";
2309   case 1005270:
2310     return "The metric in the '.sol' file contains more than one field.";
2311   case 1005300:
2312     return "The surface mesh includes at least one hole. The domain is not well defined.";
2313   case 1005301:
2314     return "Presumably, the surface mesh is not compatible with the domain being processed (warning).";
2315   case 1005302:
2316     return "Probable faces overlapping somewher.";
2317   case 1005320:
2318     return "The quadratic version does not work with prescribed free edges.";
2319   case 1005321:
2320     return "The quadratic version does not work with a volume mesh.";
2321   case 1005370:
2322     return "The metric in the '.sol' file is inadequate (values not per vertices).";
2323   case 1005371:
2324     return "The number of vertices in the '.sol' is different from the one in the "
2325       "'.mesh' file for the required vertices (option '--required_vertices').";
2326   case 1005372:
2327     return "More than one type in file NomFil. The type must be equal to 1 in the '.sol'"
2328       "for the required vertices (option '--required_vertices').";
2329   case 1005515:
2330     return "Bad vertex number.";
2331   case 1005560:
2332     return "No guess to start the definition of the connected component(s).";
2333   case 1005602:
2334     return "Some initial points cannot be inserted.";
2335   case 1005620:
2336     return "A too bad quality face is detected. This face is considered degenerated.";
2337   case 1005621:
2338     return "A too bad quality face is detected. This face is degenerated.";
2339   case 1005622:
2340     return "The algorithm cannot run further.";
2341   case 1005690:
2342     return "A too small volume element is detected.";
2343   case 1005691:
2344     return "A tetrahedra is suspected to be very bad shaped or wrong.";
2345   case 1005692:
2346     return "There is at least a null or negative volume element. The resulting mesh"
2347       "may be inappropriate.";
2348   case 1005693:
2349     return "There are some null or negative volume element. The resulting mesh may"
2350       "be inappropriate.";
2351   case 1005820:
2352     return "An edge is unique (i.e., bounds a hole in the surface).";
2353   case 1007000:
2354     return "Abnormal or internal error.";
2355   case 1007010:
2356     return "Too many components with respect to too many sub-domain.";
2357   case 1007400:
2358     return "An internal error has been encountered or a signal has been received. "
2359       "Current mesh will not be saved.";
2360   case 1008491:
2361     return "Impossible to define a component.";
2362   case 1008410:
2363     return "There are some overconstrained edges.";
2364   case 1008420:
2365     return "There are some overconstrained facets.";
2366   case 1008422:
2367     return "Give the number of missing faces (information given when regeneration phase failed).";
2368   case 1008423:
2369     return "A constrained face cannot be enforced (information given when regeneration phase failed).";
2370   case 1008441:
2371     return "A constrained edge cannot be enforced.";
2372   case 1008460:
2373     return "It is dramatically tedious to enforce the boundary items.";
2374   case 1008480:
2375     return "The surface mesh regeneration step has failed. A .boite.mesh and .boite.map files are created.";
2376   case 1008490:
2377     return "Invalid resulting mesh.";
2378   case 1008495:
2379     return "P2 correction not successful.";
2380   case 1009000:
2381     return "Program has received an interruption or a termination signal sent by the "
2382       "user or the system administrator. Current mesh will not be saved.";
2383   }
2384   return "";
2385 }
2386
2387 //================================================================================
2388 /*!
2389  * \brief Retrieve from a string given number of integers
2390  */
2391 //================================================================================
2392
2393 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
2394 {
2395   ids.clear();
2396   ids.reserve( nbIds );
2397   while ( nbIds )
2398   {
2399     while ( !isdigit( *ptr )) ++ptr;
2400     if ( ptr[-1] == '-' ) --ptr;
2401     ids.push_back( strtol( ptr, &ptr, 10 ));
2402     --nbIds;
2403   }
2404   return ptr;
2405 }
2406
2407 //================================================================================
2408 /*!
2409  * \brief Retrieve problem description form a log file
2410  *  \retval bool - always false
2411  */
2412 //================================================================================
2413
2414 SMESH_ComputeErrorPtr
2415 GHS3DPlugin_GHS3D::getErrorDescription(const char*                logFile,
2416                                        const std::string&         log,
2417                                        const _Ghs2smdsConvertor & toSmdsConvertor,
2418                                        const bool                 isOk/* = false*/ )
2419 {
2420   SMESH_BadInputElements* badElemsErr =
2421     new SMESH_BadInputElements( toSmdsConvertor.getMesh(), COMPERR_ALGO_FAILED );
2422   SMESH_ComputeErrorPtr err( badElemsErr );
2423
2424   char* ptr = const_cast<char*>( log.c_str() );
2425   char* buf = ptr, * bufEnd = ptr + log.size();
2426
2427
2428   SMESH_Comment errDescription;
2429
2430   enum { NODE = 1, EDGE, TRIA, VOL, SKIP_ID = 1 };
2431
2432   // look for MeshGems version
2433   // Since "MG-TETRA -- MeshGems 1.1-3 (January, 2013)" error codes change.
2434   // To discriminate old codes from new ones we add 1000000 to the new codes.
2435   // This way value of the new codes is same as absolute value of codes printed
2436   // in the log after "MGMESSAGE" string.
2437   int versionAddition = 0;
2438   {
2439     char* verPtr = ptr;
2440     while ( ++verPtr < bufEnd )
2441     {
2442       if ( strncmp( verPtr, "MG-TETRA -- MeshGems ", 21 ) != 0 )
2443         continue;
2444       if ( strcmp( verPtr, "MG-TETRA -- MeshGems 1.1-3 " ) >= 0 )
2445         versionAddition = 1000000;
2446       ptr = verPtr;
2447       break;
2448     }
2449   }
2450
2451   // look for errors "ERR #"
2452
2453   set<string> foundErrorStr; // to avoid reporting same error several times
2454   set<int>    elemErrorNums; // not to report different types of errors with bad elements
2455   while ( ++ptr < bufEnd )
2456   {
2457     if ( strncmp( ptr, "ERR ", 4 ) != 0 )
2458       continue;
2459
2460     list<const SMDS_MeshElement*>& badElems = badElemsErr->myBadElements;
2461     vector<int> nodeIds;
2462
2463     ptr += 4;
2464     char* errBeg = ptr;
2465     int   errNum = strtol(ptr, &ptr, 10) + versionAddition;
2466     // we treat errors enumerated in [SALOME platform 0019316] issue
2467     // and all errors from a new (Release 1.1) MeshGems User Manual
2468     switch ( errNum ) {
2469     case 0015: // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
2470     case 1005620 : // a too bad quality face is detected. This face is considered degenerated.
2471       ptr = getIds(ptr, SKIP_ID, nodeIds);
2472       ptr = getIds(ptr, TRIA, nodeIds);
2473       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2474       break;
2475     case 1005621 : // a too bad quality face is detected. This face is degenerated.
2476       // hence the is degenerated it is invisible, add its edges in addition
2477       ptr = getIds(ptr, SKIP_ID, nodeIds);
2478       ptr = getIds(ptr, TRIA, nodeIds);
2479       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2480       {
2481         vector<int> edgeNodes( nodeIds.begin(), --nodeIds.end() ); // 01
2482         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2483         edgeNodes[1] = nodeIds[2]; // 02
2484         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2485         edgeNodes[0] = nodeIds[1]; // 12
2486       }      
2487       break;
2488     case 1000: // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
2489       // ERR  1000 :  1 3 2
2490     case 1002: // Face (f 1, f 2, f 3) has a vertex negative or null
2491     case 3019: // Constrained face (f 1, f 2, f 3) cannot be enforced
2492     case 1002211: // a face has a vertex negative or null.
2493     case 1005200 : // a surface mesh appears more than once in the input surface mesh.
2494     case 1008423 : // a constrained face cannot be enforced (regeneration phase failed).
2495       ptr = getIds(ptr, TRIA, nodeIds);
2496       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2497       break;
2498     case 1001: // Edge (e1, e2) appears more than once in the input surface mesh
2499     case 3009: // Constrained edge (e1, e2) cannot be enforced (warning).
2500       // ERR  3109 :  EDGE  5 6 UNIQUE
2501     case 3109: // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
2502     case 1005210 : // an edge appears more than once in the input surface mesh.
2503     case 1005820 : // an edge is unique (i.e., bounds a hole in the surface).
2504     case 1008441 : // a constrained edge cannot be enforced.
2505       ptr = getIds(ptr, EDGE, nodeIds);
2506       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2507       break;
2508     case 2004: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
2509     case 2014: // at least two points whose distance is dist, i.e., considered as coincident
2510     case 2103: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
2511       // ERR  2103 :  16 WITH  3
2512     case 1005105 : // two vertices are too close to one another or coincident.
2513     case 1005107: // Two vertices are too close to one another or coincident.
2514       ptr = getIds(ptr, NODE, nodeIds);
2515       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2516       ptr = getIds(ptr, NODE, nodeIds);
2517       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2518       break;
2519     case 2012: // Vertex v1 cannot be inserted (warning).
2520     case 1005106 : // a vertex cannot be inserted.
2521       ptr = getIds(ptr, NODE, nodeIds);
2522       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2523       break;
2524     case 3103: // The surface edge (e1, e2) intersects another surface edge (e3, e4)
2525     case 1005110 : // two surface edges are intersecting.
2526       // ERR  3103 :  1 2 WITH  7 3
2527       ptr = getIds(ptr, EDGE, nodeIds);
2528       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2529       ptr = getIds(ptr, EDGE, nodeIds);
2530       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2531       break;
2532     case 3104: // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
2533       // ERR  3104 :  9 10 WITH  1 2 3
2534     case 3106: // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
2535     case 1005120 : // a surface edge intersects a surface face.
2536       ptr = getIds(ptr, EDGE, nodeIds);
2537       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2538       ptr = getIds(ptr, TRIA, nodeIds);
2539       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2540       break;
2541     case 3105: // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
2542       // ERR  3105 :  8 IN  2 3 5
2543     case 1005150 : // a boundary point lies within a surface face.
2544       ptr = getIds(ptr, NODE, nodeIds);
2545       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2546       ptr = getIds(ptr, TRIA, nodeIds);
2547       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2548       break;
2549     case 3107: // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
2550       // ERR  3107 :  2 IN  4 1
2551     case 1005160 : // a boundary point lies within a surface edge.
2552       ptr = getIds(ptr, NODE, nodeIds);
2553       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2554       ptr = getIds(ptr, EDGE, nodeIds);
2555       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2556       break;
2557     case 9000: // ERR  9000
2558       //  ELEMENT  261 WITH VERTICES :  7 396 -8 242
2559       //  VOLUME   : -1.11325045E+11 W.R.T. EPSILON   0.
2560       // A too small volume element is detected. Are reported the index of the element,
2561       // its four vertex indices, its volume and the tolerance threshold value
2562       ptr = getIds(ptr, SKIP_ID, nodeIds);
2563       ptr = getIds(ptr, VOL, nodeIds);
2564       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2565       // even if all nodes found, volume it most probably invisible,
2566       // add its faces to demonstrate it anyhow
2567       {
2568         vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
2569         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2570         faceNodes[2] = nodeIds[3]; // 013
2571         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2572         faceNodes[1] = nodeIds[2]; // 023
2573         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2574         faceNodes[0] = nodeIds[1]; // 123
2575         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2576       }
2577       break;
2578     case 9001: // ERR  9001
2579       //  %% NUMBER OF NEGATIVE VOLUME TETS  :  1
2580       //  %% THE LARGEST NEGATIVE TET        :   1.75376581E+11
2581       //  %%  NUMBER OF NULL VOLUME TETS     :  0
2582       // There exists at least a null or negative volume element
2583       break;
2584     case 9002:
2585       // There exist n null or negative volume elements
2586       break;
2587     case 9003:
2588       // A too small volume element is detected
2589       break;
2590     case 9102:
2591       // A too bad quality face is detected. This face is considered degenerated,
2592       // its index, its three vertex indices together with its quality value are reported
2593       break; // same as next
2594     case 9112: // ERR  9112
2595       //  FACE   2 WITH VERTICES :  4 2 5
2596       //  SMALL INRADIUS :   0.
2597       // A too bad quality face is detected. This face is degenerated,
2598       // its index, its three vertex indices together with its inradius are reported
2599       ptr = getIds(ptr, SKIP_ID, nodeIds);
2600       ptr = getIds(ptr, TRIA, nodeIds);
2601       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2602       // add triangle edges as it most probably has zero area and hence invisible
2603       {
2604         vector<int> edgeNodes(2);
2605         edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
2606         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2607         edgeNodes[1] = nodeIds[2]; // 0-2
2608         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2609         edgeNodes[0] = nodeIds[1]; // 1-2
2610         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2611       }
2612       break;
2613     case 1005103 : // the vertices of an element are too close to one another or coincident.
2614       ptr = getIds(ptr, TRIA, nodeIds);
2615       if ( nodeIds.back() == 0 ) // index of the third vertex of the element (0 for an edge)
2616         nodeIds.resize( EDGE );
2617       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2618       break;
2619     }
2620
2621     bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
2622     if ( !isNewError )
2623       continue; // not to report same error several times
2624
2625 //     const SMDS_MeshElement* nullElem = 0;
2626 //     bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
2627
2628 //     if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
2629 //       bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
2630 //       if ( oneMoreErrorType )
2631 //         continue; // not to report different types of errors with bad elements
2632 //     }
2633
2634     // make error text
2635     string text = translateError( errNum );
2636     if ( errDescription.find( text ) == text.npos ) {
2637       if ( !errDescription.empty() )
2638         errDescription << "\n";
2639       errDescription << text;
2640     }
2641
2642   } // end while
2643
2644   if ( errDescription.empty() ) { // no errors found
2645     char msgLic1[] = "connection to server failed";
2646     char msgLic2[] = " Dlim ";
2647     if ( search( &buf[0], bufEnd, msgLic1, msgLic1 + strlen(msgLic1)) != bufEnd ||
2648          search( &buf[0], bufEnd, msgLic2, msgLic2 + strlen(msgLic2)) != bufEnd )
2649       errDescription << "Licence problems.";
2650     else
2651     {
2652       char msg2[] = "SEGMENTATION FAULT";
2653       if ( search( &buf[0], bufEnd, msg2, msg2 + strlen(msg2)) != bufEnd )
2654         errDescription << "MG-Tetra: SEGMENTATION FAULT. ";
2655     }
2656   }
2657
2658   if ( !isOk && logFile && logFile[0] )
2659   {
2660     if ( errDescription.empty() )
2661       errDescription << "See " << logFile << " for problem description";
2662     else
2663       errDescription << "\nSee " << logFile << " for more information";
2664   }
2665
2666   err->myComment = errDescription;
2667
2668   if ( err->myComment.empty() && !err->HasBadElems() )
2669     err = SMESH_ComputeError::New(); // OK
2670
2671   return err;
2672 }
2673
2674 //================================================================================
2675 /*!
2676  * \brief Creates _Ghs2smdsConvertor
2677  */
2678 //================================================================================
2679
2680 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap,
2681                                         SMESH_ProxyMesh::Ptr                   mesh)
2682   :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 ), _mesh( mesh )
2683 {
2684 }
2685
2686 //================================================================================
2687 /*!
2688  * \brief Creates _Ghs2smdsConvertor
2689  */
2690 //================================================================================
2691
2692 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> &  nodeByGhsId,
2693                                         SMESH_ProxyMesh::Ptr                   mesh)
2694   : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId ), _mesh( mesh )
2695 {
2696 }
2697
2698 //================================================================================
2699 /*!
2700  * \brief Return SMDS element by ids of MG-Tetra nodes
2701  */
2702 //================================================================================
2703
2704 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
2705 {
2706   size_t nbNodes = ghsNodes.size();
2707   vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
2708   for ( size_t i = 0; i < nbNodes; ++i ) {
2709     int ghsNode = ghsNodes[ i ];
2710     if ( _ghs2NodeMap ) {
2711       map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
2712       if ( in == _ghs2NodeMap->end() )
2713         return 0;
2714       nodes[ i ] = in->second;
2715     }
2716     else {
2717       if ( ghsNode < 1 || ghsNode > (int)_nodeByGhsId->size() )
2718         return 0;
2719       nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
2720     }
2721   }
2722   if ( nbNodes == 1 )
2723     return nodes[0];
2724
2725   if ( nbNodes == 2 ) {
2726     const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
2727     if ( !edge || edge->GetID() < 1 || _mesh->IsTemporary( edge ))
2728       edge = new SMDS_LinearEdge( nodes[0], nodes[1] );
2729     return edge;
2730   }
2731   if ( nbNodes == 3 ) {
2732     const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
2733     if ( !face || face->GetID() < 1 || _mesh->IsTemporary( face ))
2734       face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
2735     return face;
2736   }
2737   if ( nbNodes == 4 )
2738     return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
2739
2740   return 0;
2741 }
2742
2743 //================================================================================
2744 /*!
2745  * \brief Return a mesh
2746  */
2747 //================================================================================
2748
2749 const SMDS_Mesh* _Ghs2smdsConvertor::getMesh() const
2750 {
2751   return _mesh->GetMeshDS();
2752 }
2753
2754 //=============================================================================
2755 /*!
2756  *
2757  */
2758 //=============================================================================
2759 bool GHS3DPlugin_GHS3D::Evaluate(SMESH_Mesh& aMesh,
2760                                  const TopoDS_Shape& aShape,
2761                                  MapShapeNbElems& aResMap)
2762 {
2763   int nbtri = 0, nbqua = 0;
2764   double fullArea = 0.0;
2765   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2766     TopoDS_Face F = TopoDS::Face( exp.Current() );
2767     SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2768     MapShapeNbElemsItr anIt = aResMap.find(sm);
2769     if( anIt==aResMap.end() ) {
2770       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2771       smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
2772                                             "Submesh can not be evaluated",this));
2773       return false;
2774     }
2775     std::vector<int> aVec = (*anIt).second;
2776     nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
2777     nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
2778     GProp_GProps G;
2779     BRepGProp::SurfaceProperties(F,G);
2780     double anArea = G.Mass();
2781     fullArea += anArea;
2782   }
2783
2784   // collect info from edges
2785   int nb0d_e = 0, nb1d_e = 0;
2786   bool IsQuadratic = false;
2787   bool IsFirst = true;
2788   TopTools_MapOfShape tmpMap;
2789   for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2790     TopoDS_Edge E = TopoDS::Edge(exp.Current());
2791     if( tmpMap.Contains(E) )
2792       continue;
2793     tmpMap.Add(E);
2794     SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
2795     MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
2796     std::vector<int> aVec = (*anIt).second;
2797     nb0d_e += aVec[SMDSEntity_Node];
2798     nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
2799     if(IsFirst) {
2800       IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
2801       IsFirst = false;
2802     }
2803   }
2804   tmpMap.Clear();
2805
2806   double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
2807
2808   GProp_GProps G;
2809   BRepGProp::VolumeProperties(aShape,G);
2810   double aVolume = G.Mass();
2811   double tetrVol = 0.1179*ELen*ELen*ELen;
2812   double CoeffQuality = 0.9;
2813   int nbVols = int(aVolume/tetrVol/CoeffQuality);
2814   int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2815   int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2816   std::vector<int> aVec(SMDSEntity_Last);
2817   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2818   if( IsQuadratic ) {
2819     aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2820     aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2821     aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2822   }
2823   else {
2824     aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2825     aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2826     aVec[SMDSEntity_Pyramid] = nbqua;
2827   }
2828   SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2829   aResMap.insert(std::make_pair(sm,aVec));
2830
2831   return true;
2832 }
2833
2834 bool GHS3DPlugin_GHS3D::importGMFMesh(const char* theGMFFileName, SMESH_Mesh& theMesh)
2835 {
2836   SMESH_ComputeErrorPtr err = theMesh.GMFToMesh( theGMFFileName, /*makeRequiredGroups =*/ true );
2837
2838   theMesh.GetMeshDS()->Modified();
2839
2840   return ( !err || err->IsOK());
2841 }
2842
2843 namespace
2844 {
2845   //================================================================================
2846   /*!
2847    * \brief Sub-mesh event listener setting enforced elements as soon as an enforced
2848    *        mesh is loaded
2849    */
2850   struct _EnforcedMeshRestorer : public SMESH_subMeshEventListener
2851   {
2852     _EnforcedMeshRestorer():
2853       SMESH_subMeshEventListener( /*isDeletable = */true, Name() )
2854     {}
2855
2856     //================================================================================
2857     /*!
2858      * \brief Returns an ID of listener
2859      */
2860     static const char* Name() { return "GHS3DPlugin_GHS3D::_EnforcedMeshRestorer"; }
2861
2862     //================================================================================
2863     /*!
2864      * \brief Treat events of the subMesh
2865      */
2866     void ProcessEvent(const int                       event,
2867                       const int                       eventType,
2868                       SMESH_subMesh*                  subMesh,
2869                       SMESH_subMeshEventListenerData* data,
2870                       const SMESH_Hypothesis*         hyp)
2871     {
2872       if ( SMESH_subMesh::SUBMESH_LOADED == event &&
2873            SMESH_subMesh::COMPUTE_EVENT  == eventType &&
2874            data &&
2875            !data->mySubMeshes.empty() )
2876       {
2877         // An enforced mesh (subMesh->_father) has been loaded from hdf file
2878         if ( GHS3DPlugin_Hypothesis* hyp = GetGHSHypothesis( data->mySubMeshes.front() ))
2879           hyp->RestoreEnfElemsByMeshes();
2880       }
2881     }
2882     //================================================================================
2883     /*!
2884      * \brief Returns GHS3DPlugin_Hypothesis used to compute a subMesh
2885      */
2886     static GHS3DPlugin_Hypothesis* GetGHSHypothesis( SMESH_subMesh* subMesh )
2887     {
2888       SMESH_HypoFilter ghsHypFilter
2889         ( SMESH_HypoFilter::HasName( GHS3DPlugin_Hypothesis::GetHypType() ));
2890       return (GHS3DPlugin_Hypothesis* )
2891         subMesh->GetFather()->GetHypothesis( subMesh->GetSubShape(),
2892                                              ghsHypFilter,
2893                                              /*visitAncestors=*/true);
2894     }
2895   };
2896
2897   //================================================================================
2898   /*!
2899    * \brief Sub-mesh event listener removing empty groups created due to "To make
2900    *        groups of domains".
2901    */
2902   struct _GroupsOfDomainsRemover : public SMESH_subMeshEventListener
2903   {
2904     _GroupsOfDomainsRemover():
2905       SMESH_subMeshEventListener( /*isDeletable = */true,
2906                                   "GHS3DPlugin_GHS3D::_GroupsOfDomainsRemover" ) {}
2907     /*!
2908      * \brief Treat events of the subMesh
2909      */
2910     void ProcessEvent(const int                       event,
2911                       const int                       eventType,
2912                       SMESH_subMesh*                  subMesh,
2913                       SMESH_subMeshEventListenerData* data,
2914                       const SMESH_Hypothesis*         hyp)
2915     {
2916       if (SMESH_subMesh::ALGO_EVENT == eventType &&
2917           !subMesh->GetAlgo() )
2918       {
2919         removeEmptyGroupsOfDomains( subMesh->GetFather(), /*notEmptyAsWell=*/true );
2920       }
2921     }
2922   };
2923 }
2924
2925 //================================================================================
2926 /*!
2927  * \brief Set an event listener to set enforced elements as soon as an enforced
2928  *        mesh is loaded
2929  */
2930 //================================================================================
2931
2932 void GHS3DPlugin_GHS3D::SubmeshRestored(SMESH_subMesh* subMesh)
2933 {
2934   if ( GHS3DPlugin_Hypothesis* hyp = _EnforcedMeshRestorer::GetGHSHypothesis( subMesh ))
2935   {
2936     GHS3DPlugin_Hypothesis::TGHS3DEnforcedMeshList enfMeshes = hyp->_GetEnforcedMeshes();
2937     GHS3DPlugin_Hypothesis::TGHS3DEnforcedMeshList::iterator it = enfMeshes.begin();
2938     for(;it != enfMeshes.end();++it) {
2939       GHS3DPlugin_Hypothesis::TGHS3DEnforcedMesh* enfMesh = *it;
2940       if ( SMESH_Mesh* mesh = GetMeshByPersistentID( enfMesh->persistID ))
2941       {
2942         SMESH_subMesh* smToListen = mesh->GetSubMesh( mesh->GetShapeToMesh() );
2943         // a listener set to smToListen will care of hypothesis stored in SMESH_EventListenerData
2944         subMesh->SetEventListener( new _EnforcedMeshRestorer(),
2945                                    SMESH_subMeshEventListenerData::MakeData( subMesh ),
2946                                    smToListen);
2947       }
2948     }
2949   }
2950 }
2951
2952 //================================================================================
2953 /*!
2954  * \brief Sets an event listener removing empty groups created due to "To make
2955  *        groups of domains".
2956  * \param subMesh - submesh where algo is set
2957  *
2958  * This method is called when a submesh gets HYP_OK algo_state.
2959  * After being set, event listener is notified on each event of a submesh.
2960  */
2961 //================================================================================
2962
2963 void GHS3DPlugin_GHS3D::SetEventListener(SMESH_subMesh* subMesh)
2964 {
2965   subMesh->SetEventListener( new _GroupsOfDomainsRemover(), 0, subMesh );
2966 }
2967
2968 //================================================================================
2969 /*!
2970  * \brief If possible, returns progress of computation [0.,1.]
2971  */
2972 //================================================================================
2973
2974 double GHS3DPlugin_GHS3D::GetProgress() const
2975 {
2976   if ( _isLibUsed )
2977   {
2978     // this->_progress is advanced by MG_Tetra_API according to messages from MG library
2979     // but sharply. Advance it a bit to get smoother advancement.
2980     GHS3DPlugin_GHS3D* me = const_cast<GHS3DPlugin_GHS3D*>( this );
2981     if ( _progress < 0.1 ) // the first message is at 10%
2982       me->_progress = GetProgressByTic();
2983     else if ( _progress < 0.98 )
2984       me->_progress += _progressAdvance;
2985     return _progress;
2986   }
2987
2988   return -1;
2989 }