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