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