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