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