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