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