Salome HOME
Merge remote-tracking branch 'origin/V8_3_BR' into gdd/python3_dev
[plugins/hybridplugin.git] / src / HYBRIDPlugin / HYBRIDPlugin_HYBRID.cxx
1 // Copyright (C) 2007-2016  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // ---
21 // File   : HYBRIDPlugin_HYBRID.cxx
22 // Author : Christian VAN WAMBEKE (CEA) (from GHS3D plugin V730)
23 // ---
24 //
25 #include "HYBRIDPlugin_HYBRID.hxx"
26 #include "HYBRIDPlugin_Hypothesis.hxx"
27 #include "MG_HYBRID_API.hxx"
28
29 #include <SMDS_FaceOfNodes.hxx>
30 #include <SMDS_LinearEdge.hxx>
31 #include <SMDS_VolumeOfNodes.hxx>
32 #include <SMESHDS_Group.hxx>
33 #include <SMESHDS_Mesh.hxx>
34 #include <SMESH_Comment.hxx>
35 #include <SMESH_File.hxx>
36 #include <SMESH_Group.hxx>
37 #include <SMESH_HypoFilter.hxx>
38 #include <SMESH_Mesh.hxx>
39 #include <SMESH_MeshAlgos.hxx>
40 #include <SMESH_MeshEditor.hxx>
41 #include <SMESH_MesherHelper.hxx>
42 #include <SMESH_ProxyMesh.hxx>
43 #include <SMESH_subMeshEventListener.hxx>
44
45 #include <BRepAdaptor_Surface.hxx>
46 #include <BRepBndLib.hxx>
47 #include <BRepBuilderAPI_MakeVertex.hxx>
48 #include <BRepClass3d.hxx>
49 #include <BRepClass3d_SolidClassifier.hxx>
50 #include <BRepExtrema_DistShapeShape.hxx>
51 #include <BRepGProp.hxx>
52 #include <BRepTools.hxx>
53 #include <BRep_Tool.hxx>
54 #include <Bnd_Box.hxx>
55 #include <GProp_GProps.hxx>
56 #include <GeomAPI_ProjectPointOnSurf.hxx>
57 #include <Precision.hxx>
58 #include <Standard_ErrorHandler.hxx>
59 #include <Standard_Failure.hxx>
60 #include <Standard_ProgramError.hxx>
61 #include <TopExp.hxx>
62 #include <TopExp_Explorer.hxx>
63 #include <TopTools_IndexedMapOfShape.hxx>
64 #include <TopTools_ListIteratorOfListOfShape.hxx>
65 #include <TopTools_MapOfShape.hxx>
66 #include <TopoDS.hxx>
67 #include <TopoDS_Shell.hxx>
68 #include <TopoDS_Solid.hxx>
69
70 #include <Basics_Utils.hxx>
71 #include <utilities.h>
72
73 #include <algorithm>
74
75 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
76
77 #ifndef GMFVERSION
78 #define GMFVERSION GmfDouble
79 #endif
80 #define GMFDIMENSION 3
81
82 #define HOLE_ID -1
83
84 typedef const std::list<const SMDS_MeshFace*> TTriaList;
85
86 static const char theDomainGroupNamePrefix[] = "Domain_";
87
88 static void removeFile( const std::string& fileName )
89 {
90   try {
91     SMESH_File( fileName ).remove();
92   }
93   catch ( ... ) {
94     MESSAGE("Can't remove file: " << fileName << " ; file does not exist or permission denied");
95   }
96 }
97
98 //=============================================================================
99 /*!
100  *  
101  */
102 //=============================================================================
103
104 HYBRIDPlugin_HYBRID::HYBRIDPlugin_HYBRID(int hypId, int studyId, SMESH_Gen* gen)
105   : SMESH_3D_Algo(hypId, studyId, gen)
106 {
107   _name = Name();
108   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
109   _onlyUnaryInput = true; // Compute() will be called on each solid
110   _iShape=0;
111   _nbShape=0;
112   _compatibleHypothesis.push_back( HYBRIDPlugin_Hypothesis::GetHypType());
113   _requireShape = false; // can work without shape_studyId
114
115   smeshGen_i = SMESH_Gen_i::GetSMESHGen();
116   CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
117   SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
118
119   myStudy = NULL;
120   myStudy = aStudyMgr->GetStudyByID(_studyId);
121   
122   _computeCanceled = false;
123 }
124
125 //=============================================================================
126 /*!
127  *  
128  */
129 //=============================================================================
130
131 HYBRIDPlugin_HYBRID::~HYBRIDPlugin_HYBRID()
132 {
133 }
134
135 //=============================================================================
136 /*!
137  *  
138  */
139 //=============================================================================
140
141 bool HYBRIDPlugin_HYBRID::CheckHypothesis ( SMESH_Mesh&         aMesh,
142                                           const TopoDS_Shape& aShape,
143                                           Hypothesis_Status&  aStatus )
144 {
145   aStatus = SMESH_Hypothesis::HYP_OK;
146
147   _hyp = 0;
148   _keepFiles = false;
149   _removeLogOnSuccess = true;
150   _logInStandardOutput = false;
151
152   const std::list <const SMESHDS_Hypothesis * >& hyps =
153     GetUsedHypothesis(aMesh, aShape, /*ignoreAuxiliary=*/false);
154   std::list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
155   for ( ; h != hyps.end(); ++h )
156   {
157     if ( !_hyp )
158       _hyp = dynamic_cast< const HYBRIDPlugin_Hypothesis*> ( *h );
159   }
160   if ( _hyp )
161   {
162     _keepFiles = _hyp->GetKeepFiles();
163     _removeLogOnSuccess = _hyp->GetRemoveLogOnSuccess();
164     _logInStandardOutput = _hyp->GetStandardOutputLog();
165   }
166
167   return true;
168 }
169
170
171 //=======================================================================
172 //function : entryToShape
173 //purpose  : 
174 //=======================================================================
175
176 TopoDS_Shape HYBRIDPlugin_HYBRID::entryToShape(std::string entry)
177 {
178   if ( myStudy->_is_nil() )
179     throw SALOME_Exception("MG-HYBRID plugin can't work w/o publishing in the study");
180   GEOM::GEOM_Object_var aGeomObj;
181   TopoDS_Shape S = TopoDS_Shape();
182   SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
183   if (!aSObj->_is_nil() ) {
184     CORBA::Object_var obj = aSObj->GetObject();
185     aGeomObj = GEOM::GEOM_Object::_narrow(obj);
186     aSObj->UnRegister();
187   }
188   if ( !aGeomObj->_is_nil() )
189     S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
190   return S;
191 }
192
193 //=======================================================================
194 //function : addElemInMeshGroup
195 //purpose  : Update or create groups in mesh
196 //=======================================================================
197
198 static void addElemInMeshGroup(SMESH_Mesh*             theMesh,
199                                const SMDS_MeshElement* anElem,
200                                std::string&            groupName,
201                                std::set<std::string>&  groupsToRemove)
202 {
203   if ( !anElem ) return; // issue 0021776
204
205   bool groupDone = false;
206   SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
207   while (grIt->more()) {
208     SMESH_Group * group = grIt->next();
209     if ( !group ) continue;
210     SMESHDS_GroupBase* groupDS = group->GetGroupDS();
211     if ( !groupDS ) continue;
212     if ( groupDS->GetType()==anElem->GetType() &&groupName.compare(group->GetName())==0) {
213       SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
214       aGroupDS->SMDSGroup().Add(anElem);
215       groupDone = true;
216       break;
217     }
218   }
219   
220   if (!groupDone)
221   {
222     int groupId;
223     SMESH_Group* aGroup = theMesh->AddGroup(anElem->GetType(), groupName.c_str(), groupId);
224     aGroup->SetName( groupName.c_str() );
225     SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
226     aGroupDS->SMDSGroup().Add(anElem);
227     groupDone = true;
228   }
229   if (!groupDone)
230     throw SALOME_Exception(LOCALIZED("A given element was not added to a group"));
231 }
232
233
234 //=======================================================================
235 //function : updateMeshGroups
236 //purpose  : Update or create groups in mesh
237 //=======================================================================
238
239 static void updateMeshGroups(SMESH_Mesh* theMesh, std::set<std::string> groupsToRemove)
240 {
241   SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
242   while (grIt->more()) {
243     SMESH_Group * group = grIt->next();
244     if ( !group ) continue;
245     SMESHDS_GroupBase* groupDS = group->GetGroupDS();
246     if ( !groupDS ) continue;
247     std::string currentGroupName = (std::string)group->GetName();
248     if (groupDS->IsEmpty() && groupsToRemove.find(currentGroupName) != groupsToRemove.end()) {
249       // Previous group created by enforced elements
250       theMesh->RemoveGroup(groupDS->GetID());
251     }
252   }
253 }
254
255 //=======================================================================
256 //function : removeEmptyGroupsOfDomains
257 //purpose  : remove empty groups named "Domain_nb" created due to
258 //           "To make groups of domains" option.
259 //=======================================================================
260
261 static void removeEmptyGroupsOfDomains(SMESH_Mesh* mesh,
262                                        bool        notEmptyAsWell = false)
263 {
264   const char* refName = theDomainGroupNamePrefix;
265   const size_t refLen = strlen( theDomainGroupNamePrefix );
266
267   std::list<int> groupIDs = mesh->GetGroupIds();
268   std::list<int>::const_iterator id = groupIDs.begin();
269   for ( ; id != groupIDs.end(); ++id )
270   {
271     SMESH_Group* group = mesh->GetGroup( *id );
272     if ( !group || ( !group->GetGroupDS()->IsEmpty() && !notEmptyAsWell ))
273       continue;
274     const char* name = group->GetName();
275     char* end;
276     // check the name
277     if ( strncmp( name, refName, refLen ) == 0 && // starts from refName;
278          isdigit( *( name + refLen )) &&          // refName is followed by a digit;
279          strtol( name + refLen, &end, 10) >= 0 && // there are only digits ...
280          *end == '\0')                            // ... till a string end.
281     {
282       mesh->RemoveGroup( *id );
283     }
284   }
285 }
286
287 //================================================================================
288 /*!
289  * \brief Create the groups corresponding to domains
290  */
291 //================================================================================
292
293 static void makeDomainGroups( std::vector< std::vector< const SMDS_MeshElement* > >& elemsOfDomain,
294                               SMESH_MesherHelper*                                    theHelper)
295 {
296   for ( size_t iDomain = 0; iDomain < elemsOfDomain.size(); ++iDomain )
297   {
298     std::vector< const SMDS_MeshElement* > & elems = elemsOfDomain[ iDomain ];
299     if ( elems.empty() ) continue;
300
301     // find existing groups
302     std::vector< SMESH_Group* > groupOfType( SMDSAbs_NbElementTypes, (SMESH_Group*)NULL );
303     const std::string domainName = ( SMESH_Comment( theDomainGroupNamePrefix ) << iDomain );
304     SMESH_Mesh::GroupIteratorPtr groupIt = theHelper->GetMesh()->GetGroups();
305     while ( groupIt->more() )
306     {
307       SMESH_Group* group = groupIt->next();
308       if ( domainName == group->GetName() &&
309            dynamic_cast< SMESHDS_Group* >( group->GetGroupDS()) )
310         groupOfType[ group->GetGroupDS()->GetType() ] = group;
311     }
312     // create and fill the groups
313     size_t iElem = 0;
314     int groupID;
315     do
316     {
317       SMESH_Group* group = groupOfType[ elems[ iElem ]->GetType() ];
318       if ( !group )
319         group = theHelper->GetMesh()->AddGroup( elems[ iElem ]->GetType(),
320                                                 domainName.c_str(), groupID );
321       SMDS_MeshGroup& groupDS =
322         static_cast< SMESHDS_Group* >( group->GetGroupDS() )->SMDSGroup();
323
324       while ( iElem < elems.size() && groupDS.Add( elems[iElem] ))
325         ++iElem;
326
327     } while ( iElem < elems.size() );
328   }
329 }
330
331 //=======================================================================
332 //function : readGMFFile
333 //purpose  : read GMF file w/o geometry associated to mesh
334 //=======================================================================
335
336 static bool readGMFFile(MG_HYBRID_API*                          MGOutput,
337                         const char*                             theFile,
338                         HYBRIDPlugin_HYBRID*                    theAlgo,
339                         SMESH_MesherHelper*                     theHelper,
340                         std::vector <const SMDS_MeshNode*> &    theNodeByHybridId,
341                         std::vector <const SMDS_MeshElement*> & theFaceByHybridId,
342                         std::map<const SMDS_MeshNode*,int> &    theNodeToHybridIdMap,
343                         std::vector<std::string> &              aNodeGroupByHybridId,
344                         std::vector<std::string> &              anEdgeGroupByHybridId,
345                         std::vector<std::string> &              aFaceGroupByHybridId,
346                         std::set<std::string> &                 groupsToRemove,
347                         bool                                    toMakeGroupsOfDomains=false,
348                         bool                                    toMeshHoles=true)
349 {
350   std::string tmpStr;
351   SMESHDS_Mesh* theMeshDS = theHelper->GetMeshDS();
352   const bool hasGeom = ( theHelper->GetMesh()->HasShapeToMesh() );
353
354   // if imprinting, the original mesh faces are modified
355   // => we clear all the faces to retrieve them from Hybrid output mesh.
356   std::vector<int> facesWithImprinting = theAlgo->getHyp()->GetFacesWithImprinting();
357
358   if ( ! facesWithImprinting.empty() ) {
359 #ifdef _DEBUG_
360       std::cout << "Imprinting => Clear original mesh" << std::endl;
361 #endif
362       SMDS_ElemIteratorPtr eIt = theMeshDS->elementsIterator();
363       while( eIt->more() )
364         theMeshDS->RemoveFreeElement( eIt->next(), /*sm=*/0 );
365       SMDS_NodeIteratorPtr nIt = theMeshDS->nodesIterator();
366       while ( nIt->more() )
367         theMeshDS->RemoveFreeNode( nIt->next(), /*sm=*/0 );
368
369       theNodeByHybridId.clear();
370       theFaceByHybridId.clear();
371   }
372
373   int nbMeshNodes = theMeshDS->NbNodes();
374   int nbInitialNodes = theNodeByHybridId.size();
375
376   const bool isQuadMesh = 
377     theHelper->GetMesh()->NbEdges( ORDER_QUADRATIC ) ||
378     theHelper->GetMesh()->NbFaces( ORDER_QUADRATIC ) ||
379     theHelper->GetMesh()->NbVolumes( ORDER_QUADRATIC );
380     
381 #ifdef _DEBUG_
382   std::cout << "theNodeByHybridId.size(): " << nbInitialNodes << std::endl;
383   std::cout << "theHelper->GetMesh()->NbNodes(): " << nbMeshNodes << std::endl;
384   std::cout << "isQuadMesh: " << isQuadMesh << std::endl;
385 #endif
386   
387   // ---------------------------------
388   // Read generated elements and nodes
389   // ---------------------------------
390
391   int nbElem = 0, nbRef = 0;
392   int aGMFNodeID = 0;
393   std::vector< const SMDS_MeshNode* > GMFNode;
394 #ifdef _DEBUG_
395   std::map<int, std::set<int> > subdomainId2tetraId;
396 #endif
397   std::map <GmfKwdCod,int> tabRef;
398   const bool force3d = !hasGeom;
399   const int  noID    = 0;
400
401   tabRef[GmfVertices]       = 3; // for new nodes and enforced nodes
402   tabRef[GmfCorners]        = 1;
403   tabRef[GmfEdges]          = 2; // for enforced edges
404   tabRef[GmfRidges]         = 1;
405   tabRef[GmfTriangles]      = 3; // for enforced faces
406   tabRef[GmfQuadrilaterals] = 4;
407   tabRef[GmfTetrahedra]     = 4; // for new tetras
408   tabRef[GmfPyramids]       = 5; // for new pyramids
409   tabRef[GmfPrisms]         = 6; // for new prisms
410   tabRef[GmfHexahedra]      = 8;
411
412   int ver, dim;
413   int InpMsh = MGOutput->GmfOpenMesh(theFile, GmfRead, &ver, &dim);
414   if (!InpMsh)
415     return false;
416
417   // Hybrid is not multi-domain => We can't (and don't need to) read ids of domains in ouput file like in GHS3DPlugin
418   // We just need to get the id of the one and only solid
419   int solidID = 1;
420   if ( hasGeom )
421   {
422     if ( theHelper->GetSubShape().ShapeType() == TopAbs_SOLID )
423       solidID = theHelper->GetSubShapeID();
424     else
425       solidID = theMeshDS->ShapeToIndex
426         ( TopExp_Explorer( theHelper->GetSubShape(), TopAbs_SOLID ).Current() );
427   }
428
429   // Issue 0020682. Avoid creating nodes and tetras at place where
430   // volumic elements already exist
431   SMESH_ElementSearcher* elemSearcher = 0;
432   std::vector< const SMDS_MeshElement* > foundVolumes;
433   if ( !hasGeom && theHelper->GetMesh()->NbVolumes() > 0 )
434     elemSearcher = SMESH_MeshAlgos::GetElementSearcher( *theMeshDS );
435   SMESHUtils::Deleter< SMESH_ElementSearcher > elemSearcherDeleter( elemSearcher );
436
437   // IMP 0022172: [CEA 790] create the groups corresponding to domains
438   std::vector< std::vector< const SMDS_MeshElement* > > elemsOfDomain;
439
440   int nbVertices = MGOutput->GmfStatKwd(InpMsh, GmfVertices) - nbInitialNodes;
441   if ( nbVertices < 0 )
442     return false;
443   GMFNode.resize( nbVertices + 1 );
444
445   std::map <GmfKwdCod,int>::const_iterator it = tabRef.begin();
446   for ( ; it != tabRef.end() ; ++it)
447   {
448     if(theAlgo->computeCanceled()) {
449       MGOutput->GmfCloseMesh(InpMsh);
450       return false;
451     }
452     int dummy;
453     GmfKwdCod token = it->first;
454     nbRef           = it->second;
455
456     nbElem = MGOutput->GmfStatKwd(InpMsh, token);
457     if (nbElem > 0) {
458       MGOutput->GmfGotoKwd(InpMsh, token);
459       std::cout << "Read " << nbElem;
460     }
461     else
462       continue;
463
464     std::vector<int> id (nbElem*tabRef[token]); // node ids
465     std::vector<int> domainID( nbElem ); // domain
466
467     if (token == GmfVertices) {
468       (nbElem <= 1) ? tmpStr = " vertex" : tmpStr = " vertices";
469       
470       int aGMFID;
471       float VerTab_f[3];
472       double x, y, z;
473       const SMDS_MeshNode * aGMFNode;
474
475       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
476         if(theAlgo->computeCanceled()) {
477           MGOutput->GmfCloseMesh(InpMsh);
478           return false;
479         }
480         if (ver == GmfFloat) {
481           MGOutput->GmfGetLin(InpMsh, token, &VerTab_f[0], &VerTab_f[1], &VerTab_f[2], &dummy);
482           x = VerTab_f[0];
483           y = VerTab_f[1];
484           z = VerTab_f[2];
485         }
486         else {
487           MGOutput->GmfGetLin(InpMsh, token, &x, &y, &z, &dummy);
488         }
489         if (iElem >= nbInitialNodes) {
490           if ( elemSearcher &&
491                 elemSearcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_Volume, foundVolumes))
492             aGMFNode = 0;
493           else
494             aGMFNode = theHelper->AddNode(x, y, z);
495           
496           aGMFID = iElem -nbInitialNodes +1;
497           GMFNode[ aGMFID ] = aGMFNode;
498           if (aGMFID-1 < (int)aNodeGroupByHybridId.size() && !aNodeGroupByHybridId.at(aGMFID-1).empty())
499             addElemInMeshGroup(theHelper->GetMesh(), aGMFNode, aNodeGroupByHybridId.at(aGMFID-1), groupsToRemove);
500         }
501       }
502     }
503     else if (token == GmfCorners && nbElem > 0) {
504       (nbElem <= 1) ? tmpStr = " corner" : tmpStr = " corners";
505       for ( int iElem = 0; iElem < nbElem; iElem++ )
506         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]]);
507     }
508     else if (token == GmfRidges && nbElem > 0) {
509       (nbElem <= 1) ? tmpStr = " ridge" : tmpStr = " ridges";
510       for ( int iElem = 0; iElem < nbElem; iElem++ )
511         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]]);
512     }
513     else if (token == GmfEdges && nbElem > 0) {
514       (nbElem <= 1) ? tmpStr = " edge" : tmpStr = " edges";
515       for ( int iElem = 0; iElem < nbElem; iElem++ )
516         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &domainID[iElem]);
517     }
518     else if (token == GmfTriangles && nbElem > 0) {
519       (nbElem <= 1) ? tmpStr = " triangle" : tmpStr = " triangles";
520       for ( int iElem = 0; iElem < nbElem; iElem++ )
521         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &domainID[iElem]);
522     }
523     else if (token == GmfQuadrilaterals && nbElem > 0) {
524       (nbElem <= 1) ? tmpStr = " Quadrilateral" : tmpStr = " Quadrilaterals";
525       for ( int iElem = 0; iElem < nbElem; iElem++ )
526         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]);
527     }
528     else if (token == GmfTetrahedra && nbElem > 0) {
529       (nbElem <= 1) ? tmpStr = " Tetrahedron" : tmpStr = " Tetrahedra";
530       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
531         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]);
532 #ifdef _DEBUG_
533         subdomainId2tetraId[dummy].insert(iElem+1);
534 #endif
535       }
536     }
537     else if (token == GmfPyramids && nbElem > 0) {
538       (nbElem <= 1) ? tmpStr = " Pyramid" : tmpStr = " Pyramids";
539       for ( int iElem = 0; iElem < nbElem; iElem++ )
540         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
541                   &id[iElem*tabRef[token]+4], &domainID[iElem]);
542     }
543     else if (token == GmfPrisms && nbElem > 0) {
544       (nbElem <= 1) ? tmpStr = " Prism" : tmpStr = " Prisms";
545       for ( int iElem = 0; iElem < nbElem; iElem++ )
546         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
547                   &id[iElem*tabRef[token]+4], &id[iElem*tabRef[token]+5], &domainID[iElem]);
548     }
549     else if (token == GmfHexahedra && nbElem > 0) {
550       (nbElem <= 1) ? tmpStr = " Hexahedron" : tmpStr = " Hexahedra";
551       for ( int iElem = 0; iElem < nbElem; iElem++ )
552         MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
553                   &id[iElem*tabRef[token]+4], &id[iElem*tabRef[token]+5], &id[iElem*tabRef[token]+6], &id[iElem*tabRef[token]+7], &domainID[iElem]);
554     }
555     std::cout << tmpStr << std::endl;
556     //std::cout << std::endl;
557
558     switch (token) {
559     case GmfCorners:
560     case GmfRidges:
561     case GmfEdges:
562     case GmfTriangles:
563     case GmfQuadrilaterals:
564     case GmfTetrahedra:
565     case GmfPyramids:
566     case GmfPrisms:
567     case GmfHexahedra:
568     {
569       std::vector< const SMDS_MeshNode* > node( nbRef );
570       std::vector< int >          nodeID( nbRef );
571       std::vector< SMDS_MeshNode* > enfNode( nbRef );
572       const SMDS_MeshElement* aCreatedElem;
573
574       for ( int iElem = 0; iElem < nbElem; iElem++ )
575       {
576         if(theAlgo->computeCanceled()) {
577           MGOutput->GmfCloseMesh(InpMsh);
578           return false;
579         }
580         // Check if elem is already in input mesh. If yes => skip
581         bool fullyCreatedElement = false; // if at least one of the nodes was created
582         for ( int iRef = 0; iRef < nbRef; iRef++ )
583         {
584           aGMFNodeID = id[iElem*tabRef[token]+iRef]; // read nbRef aGMFNodeID
585           if (aGMFNodeID <= nbInitialNodes) // input nodes
586           {
587             aGMFNodeID--;
588             node[ iRef ] = theNodeByHybridId[aGMFNodeID];
589           }
590           else
591           {
592             fullyCreatedElement = true;
593             aGMFNodeID -= nbInitialNodes;
594             nodeID[ iRef ] = aGMFNodeID ;
595             node  [ iRef ] = GMFNode[ aGMFNodeID ];
596           }
597         }
598         aCreatedElem = 0;
599         switch (token)
600         {
601         case GmfEdges:
602           if (fullyCreatedElement) {
603             aCreatedElem = theHelper->AddEdge( node[0], node[1], noID, force3d );
604             if ( !anEdgeGroupByHybridId.empty() && !anEdgeGroupByHybridId[iElem].empty())
605               addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, anEdgeGroupByHybridId[iElem], groupsToRemove);
606           }
607           break;
608         case GmfTriangles:
609           if (fullyCreatedElement) {
610             aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], noID, force3d );
611             // add iElem < aFaceGroupByHybridId.size() to avoid crash if imprinting with hexa core with MeshGems <= 2.4-5
612             if ( !aFaceGroupByHybridId.empty() && iElem < aFaceGroupByHybridId.size() && !aFaceGroupByHybridId[iElem].empty() ) {
613               addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, aFaceGroupByHybridId[iElem], groupsToRemove);
614             }
615             // add element in shape for groups on geom to work
616             theMeshDS->SetMeshElementOnShape( aCreatedElem, domainID[iElem] );
617             for ( int iN = 0; iN < 3; ++iN )
618               if ( node[iN]->getshapeId() < 1 )
619                 theMeshDS->SetNodeOnFace( node[iN], domainID[iElem] );
620           }
621           break;
622         case GmfQuadrilaterals:
623           if (fullyCreatedElement) {
624             aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], node[3], noID, force3d );
625             // add element in shape for groups on geom to work
626             theMeshDS->SetMeshElementOnShape( aCreatedElem, domainID[iElem] );
627             for ( int iN = 0; iN < 3; ++iN )
628               if ( node[iN]->getshapeId() < 1 )
629                 theMeshDS->SetNodeOnFace( node[iN], domainID[iElem] );
630           }
631           break;
632         case GmfTetrahedra:
633           if ( hasGeom )
634           {
635             if ( solidID != HOLE_ID )
636             {
637               aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
638                                                    noID, force3d );
639               theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
640               for ( int iN = 0; iN < 4; ++iN )
641                 if ( node[iN]->getshapeId() < 1 )
642                   theMeshDS->SetNodeInVolume( node[iN], solidID );
643             }
644           }
645           else
646           {
647             if ( elemSearcher ) {
648               // Issue 0020682. Avoid creating nodes and tetras at place where
649               // volumic elements already exist
650               if ( !node[1] || !node[0] || !node[2] || !node[3] )
651                 continue;
652               if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
653                                                       SMESH_TNodeXYZ(node[1]) +
654                                                       SMESH_TNodeXYZ(node[2]) +
655                                                       SMESH_TNodeXYZ(node[3]) ) / 4.,
656                                                      SMDSAbs_Volume, foundVolumes ))
657                 break;
658             }
659             aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
660                                                  noID, force3d );
661           }
662           break;
663         case GmfPyramids:
664           if ( hasGeom )
665           {
666             if ( solidID != HOLE_ID )
667             {
668               aCreatedElem = theHelper->AddVolume( node[3], node[2], node[1],
669                                                    node[0], node[4],
670                                                    noID, force3d );
671               theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
672               for ( int iN = 0; iN < 5; ++iN )
673                 if ( node[iN]->getshapeId() < 1 )
674                   theMeshDS->SetNodeInVolume( node[iN], solidID );
675             }
676           }
677           else
678           {
679             if ( elemSearcher ) {
680               // Issue 0020682. Avoid creating nodes and tetras at place where
681               // volumic elements already exist
682               if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] )
683                 continue;
684               if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
685                                                       SMESH_TNodeXYZ(node[1]) +
686                                                       SMESH_TNodeXYZ(node[2]) +
687                                                       SMESH_TNodeXYZ(node[3]) +
688                                                       SMESH_TNodeXYZ(node[4])) / 5.,
689                                                      SMDSAbs_Volume, foundVolumes ))
690                 break;
691             }
692             aCreatedElem = theHelper->AddVolume( node[3], node[2], node[1],
693                                                    node[0], node[4],
694                                                    noID, force3d );
695           }
696           break;
697         case GmfPrisms:
698           if ( hasGeom )
699           {
700             if ( solidID != HOLE_ID )
701             {
702               aCreatedElem = theHelper->AddVolume( node[0], node[2], node[1],
703                                                    node[3], node[5], node[4],
704                                                    noID, force3d );
705               theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
706               for ( int iN = 0; iN < 6; ++iN )
707                 if ( node[iN]->getshapeId() < 1 )
708                   theMeshDS->SetNodeInVolume( node[iN], solidID );
709             }
710           }
711           else
712           {
713             if ( elemSearcher ) {
714               // Issue 0020682. Avoid creating nodes and tetras at place where
715               // volumic elements already exist
716               if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] )
717                 continue;
718               if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
719                                                       SMESH_TNodeXYZ(node[1]) +
720                                                       SMESH_TNodeXYZ(node[2]) +
721                                                       SMESH_TNodeXYZ(node[3]) +
722                                                       SMESH_TNodeXYZ(node[4]) +
723                                                       SMESH_TNodeXYZ(node[5])) / 6.,
724                                                      SMDSAbs_Volume, foundVolumes ))
725                 break;
726             }
727             aCreatedElem = theHelper->AddVolume( node[0], node[2], node[1],
728                                                  node[3], node[5], node[4],
729                                                  noID, force3d );
730           }
731           break;
732         case GmfHexahedra:
733           if ( hasGeom )
734           {
735             if ( solidID != HOLE_ID )
736             {
737               aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
738                                                    node[4], node[7], node[6], node[5],
739                                                    noID, force3d );
740               theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
741               for ( int iN = 0; iN < 8; ++iN )
742                 if ( node[iN]->getshapeId() < 1 )
743                   theMeshDS->SetNodeInVolume( node[iN], solidID );
744             }
745           }
746           else
747           {
748             if ( elemSearcher ) {
749               // Issue 0020682. Avoid creating nodes and tetras at place where
750               // volumic elements already exist
751               if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] || !node[6] || !node[7])
752                 continue;
753               if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
754                                                       SMESH_TNodeXYZ(node[1]) +
755                                                       SMESH_TNodeXYZ(node[2]) +
756                                                       SMESH_TNodeXYZ(node[3]) +
757                                                       SMESH_TNodeXYZ(node[4]) +
758                                                       SMESH_TNodeXYZ(node[5]) +
759                                                       SMESH_TNodeXYZ(node[6]) +
760                                                       SMESH_TNodeXYZ(node[7])) / 8.,
761                                                      SMDSAbs_Volume, foundVolumes ))
762                 break;
763             }
764             aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
765                                                  node[4], node[7], node[6], node[5],
766                                                  noID, force3d );
767           }
768           break;
769         default: continue;
770         } // switch (token)
771
772         if ( aCreatedElem && toMakeGroupsOfDomains )
773         {
774           if ( domainID[iElem] >= (int) elemsOfDomain.size() )
775             elemsOfDomain.resize( domainID[iElem] + 1 );
776           elemsOfDomain[ domainID[iElem] ].push_back( aCreatedElem );
777         }
778       } // loop on elements of one type
779       break;
780     } // case ...
781
782     default:;
783     } // switch (token)
784   } // loop on tabRef
785
786   // remove nodes in holes
787   if ( hasGeom )
788   {
789     for ( int i = 1; i <= nbVertices; ++i )
790       if ( GMFNode[i]->NbInverseElements() == 0 )
791         theMeshDS->RemoveFreeNode( GMFNode[i], /*sm=*/0, /*fromGroups=*/false );
792   }
793
794   MGOutput->GmfCloseMesh(InpMsh);
795
796   // 0022172: [CEA 790] create the groups corresponding to domains
797   if ( toMakeGroupsOfDomains )
798     makeDomainGroups( elemsOfDomain, theHelper );
799
800 #ifdef _DEBUG_
801   std::map<int, std::set<int> >::const_iterator subdomainIt = subdomainId2tetraId.begin();
802   std::string aSubdomainFileName = theFile;
803   aSubdomainFileName = aSubdomainFileName + ".subdomain";
804   ofstream aSubdomainFile  ( aSubdomainFileName  , ios::out);
805
806   aSubdomainFile << "Nb subdomains " << subdomainId2tetraId.size() << std::endl;
807   for(;subdomainIt != subdomainId2tetraId.end() ; ++subdomainIt) {
808     int subdomainId = subdomainIt->first;
809     std::set<int> tetraIds = subdomainIt->second;
810     std::set<int>::const_iterator tetraIdsIt = tetraIds.begin();
811     aSubdomainFile << subdomainId << std::endl;
812     for(;tetraIdsIt != tetraIds.end() ; ++tetraIdsIt) {
813       aSubdomainFile << (*tetraIdsIt) << " ";
814     }
815     aSubdomainFile << std::endl;
816   }
817   aSubdomainFile.close();
818 #endif  
819   
820   return true;
821 }
822
823
824 static bool writeGMFFile(MG_HYBRID_API*                                  MGInput,
825                          const char*                                     theMeshFileName,
826                          const char*                                     theRequiredFileName,
827                          const char*                                     theSolFileName,
828                          const SMESH_ProxyMesh&                          theProxyMesh,
829                          SMESH_MesherHelper&                             theHelper,
830                          std::vector <const SMDS_MeshNode*> &            theNodeByHybridId,
831                          std::vector <const SMDS_MeshElement*> &         theFaceByHybridId,
832                          std::map<const SMDS_MeshNode*,int> &            aNodeToHybridIdMap,
833                          std::vector<std::string> &                      aNodeGroupByHybridId,
834                          std::vector<std::string> &                      anEdgeGroupByHybridId,
835                          std::vector<std::string> &                      aFaceGroupByHybridId,
836                          HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap & theEnforcedNodes,
837                          HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedEdges,
838                          HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedTriangles,
839                          std::map<std::vector<double>, std::string> &    enfVerticesWithGroup,
840                          HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues & theEnforcedVertices)
841 {
842   std::string tmpStr;
843   int idx, idxRequired = 0, idxSol = 0;
844   //tabg each dummyint
845   //const int dummyint = 0;
846   const int dummyint1 = 1;
847   const int dummyint2 = 2;
848   const int dummyint3 = 3;
849   const int dummyint4 = 4;
850   const int enforcedTag = HYBRIDPlugin_Hypothesis::EnforcedTag();
851   //const int dummyint6 = 6; //are interesting for layers
852   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues::const_iterator vertexIt;
853   std::vector<double> enfVertexSizes;
854   const SMDS_MeshElement* elem;
855   TIDSortedElemSet anElemSetTri, anElemSetQuad, theKeptEnforcedEdges, theKeptEnforcedTriangles;
856   SMDS_ElemIteratorPtr nodeIt;
857   std::vector <const SMDS_MeshNode*> theEnforcedNodeByHybridId;
858   std::map<const SMDS_MeshNode*,int> anEnforcedNodeToHybridIdMap, anExistingEnforcedNodeToHybridIdMap;
859   std::vector< const SMDS_MeshElement* > foundElems;
860   std::map<const SMDS_MeshNode*,TopAbs_State> aNodeToTopAbs_StateMap;
861   int nbFoundElems;
862   HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap::iterator elemIt;
863   TIDSortedElemSet::iterator elemSetIt;
864   bool isOK;
865   SMESH_Mesh* theMesh = theHelper.GetMesh();
866   const bool hasGeom = theMesh->HasShapeToMesh();
867   SMESHUtils::Deleter< SMESH_ElementSearcher > pntCls
868     ( SMESH_MeshAlgos::GetElementSearcher(*theMesh->GetMeshDS()));
869   
870   int nbEnforcedVertices = theEnforcedVertices.size();
871   
872   // count faces
873   int nbFaces = theProxyMesh.NbFaces();
874   int nbNodes;
875   theFaceByHybridId.reserve( nbFaces );
876   
877   // groups management
878   int usedEnforcedNodes = 0;
879   std::string gn = "";
880
881   if ( nbFaces == 0 )
882     return false;
883   
884   idx = MGInput->GmfOpenMesh(theMeshFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
885   if (!idx)
886     return false;
887   
888   // ========================== FACES ==========================
889   // TRIANGLES ==========================
890   SMDS_ElemIteratorPtr eIt =
891     hasGeom ? theProxyMesh.GetFaces( theHelper.GetSubShape()) : theProxyMesh.GetFaces();
892   while ( eIt->more() )
893   {
894     elem = eIt->next();
895     nodeIt = elem->nodesIterator();
896     nbNodes = elem->NbCornerNodes();
897     if (nbNodes == 3)
898       anElemSetTri.insert(elem);
899     else if (nbNodes == 4)
900       anElemSetQuad.insert(elem);
901     else
902     {
903       std::cout << "Unexpected number of nodes: " << nbNodes << std::endl;
904       throw ("Unexpected number of nodes" );
905     }
906     while ( nodeIt->more() && nbNodes--)
907     {
908       // find HYBRID ID
909       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
910       int newId = aNodeToHybridIdMap.size() + 1; // hybrid ids count from 1
911       aNodeToHybridIdMap.insert( std::make_pair( node, newId ));
912     }
913   }
914   
915   //EDGES ==========================
916   
917   // Iterate over the enforced edges
918   for(elemIt = theEnforcedEdges.begin() ; elemIt != theEnforcedEdges.end() ; ++elemIt) {
919     elem = elemIt->first;
920     isOK = true;
921     nodeIt = elem->nodesIterator();
922     nbNodes = 2;
923     while ( nodeIt->more() && nbNodes-- ) {
924       // find HYBRID ID
925       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
926       // Test if point is inside shape to mesh
927       gp_Pnt myPoint(node->X(),node->Y(),node->Z());
928       TopAbs_State result = pntCls->GetPointState( myPoint );
929       if ( result == TopAbs_OUT ) {
930         isOK = false;
931         break;
932       }
933       aNodeToTopAbs_StateMap.insert( std::make_pair( node, result ));
934     }
935     if (isOK) {
936       nodeIt = elem->nodesIterator();
937       nbNodes = 2;
938       int newId = -1;
939       while ( nodeIt->more() && nbNodes-- ) {
940         // find HYBRID ID
941         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
942         gp_Pnt myPoint(node->X(),node->Y(),node->Z());
943         nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
944 #ifdef _DEBUG_
945         std::cout << "Node at "<<node->X()<<", "<<node->Y()<<", "<<node->Z()<<std::endl;
946         std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
947 #endif
948         if (nbFoundElems ==0) {
949           if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
950             newId = aNodeToHybridIdMap.size() + anEnforcedNodeToHybridIdMap.size() + 1; // hybrid ids count from 1
951             anEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
952           }
953         }
954         else if (nbFoundElems ==1) {
955           const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
956           newId = (*aNodeToHybridIdMap.find(existingNode)).second;
957           anExistingEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
958         }
959         else
960           isOK = false;
961 #ifdef _DEBUG_
962         std::cout << "HYBRID node ID: "<<newId<<std::endl;
963 #endif
964       }
965       if (isOK)
966         theKeptEnforcedEdges.insert(elem);
967     }
968   }
969   
970   //ENFORCED TRIANGLES ==========================
971   
972   // Iterate over the enforced triangles
973   for(elemIt = theEnforcedTriangles.begin() ; elemIt != theEnforcedTriangles.end() ; ++elemIt) {
974     elem = elemIt->first;
975     isOK = true;
976     nodeIt = elem->nodesIterator();
977     nbNodes = 3;
978     while ( nodeIt->more() && nbNodes--) {
979       // find HYBRID ID
980       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
981       // Test if point is inside shape to mesh
982       gp_Pnt myPoint(node->X(),node->Y(),node->Z());
983       TopAbs_State result = pntCls->GetPointState( myPoint );
984       if ( result == TopAbs_OUT ) {
985         isOK = false;
986         break;
987       }
988       aNodeToTopAbs_StateMap.insert( std::make_pair( node, result ));
989     }
990     if (isOK) {
991       nodeIt = elem->nodesIterator();
992       nbNodes = 3;
993       int newId = -1;
994       while ( nodeIt->more() && nbNodes--) {
995         // find HYBRID ID
996         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
997         gp_Pnt myPoint(node->X(),node->Y(),node->Z());
998         nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
999 #ifdef _DEBUG_
1000         std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
1001 #endif
1002         if (nbFoundElems ==0) {
1003           if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1004             newId = aNodeToHybridIdMap.size() + anEnforcedNodeToHybridIdMap.size() + 1; // hybrid ids count from 1
1005             anEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
1006           }
1007         }
1008         else if (nbFoundElems ==1) {
1009           const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1010           newId = (*aNodeToHybridIdMap.find(existingNode)).second;
1011           anExistingEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
1012         }
1013         else
1014           isOK = false;
1015 #ifdef _DEBUG_
1016         std::cout << "HYBRID node ID: "<<newId<<std::endl;
1017 #endif
1018       }
1019       if (isOK)
1020         theKeptEnforcedTriangles.insert(elem);
1021     }
1022   }
1023   
1024   // put nodes to theNodeByHybridId vector
1025 #ifdef _DEBUG_
1026   std::cout << "aNodeToHybridIdMap.size(): "<<aNodeToHybridIdMap.size()<<std::endl;
1027 #endif
1028   theNodeByHybridId.resize( aNodeToHybridIdMap.size() );
1029   std::map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToHybridIdMap.begin();
1030   for ( ; n2id != aNodeToHybridIdMap.end(); ++ n2id)
1031   {
1032 //     std::cout << "n2id->first: "<<n2id->first<<std::endl;
1033     theNodeByHybridId[ n2id->second - 1 ] = n2id->first; // hybrid ids count from 1
1034   }
1035
1036   // put nodes to anEnforcedNodeToHybridIdMap vector
1037 #ifdef _DEBUG_
1038   std::cout << "anEnforcedNodeToHybridIdMap.size(): "<<anEnforcedNodeToHybridIdMap.size()<<std::endl;
1039 #endif
1040   theEnforcedNodeByHybridId.resize( anEnforcedNodeToHybridIdMap.size());
1041   n2id = anEnforcedNodeToHybridIdMap.begin();
1042   for ( ; n2id != anEnforcedNodeToHybridIdMap.end(); ++ n2id)
1043   {
1044     if (n2id->second > (int)aNodeToHybridIdMap.size()) {
1045       theEnforcedNodeByHybridId[ n2id->second - aNodeToHybridIdMap.size() - 1 ] = n2id->first; // hybrid ids count from 1
1046     }
1047   }
1048   
1049   
1050   //========================== NODES ==========================
1051   std::vector<const SMDS_MeshNode*> theOrderedNodes, theRequiredNodes;
1052   std::set< std::vector<double> > nodesCoords;
1053   std::vector<const SMDS_MeshNode*>::const_iterator hybridNodeIt = theNodeByHybridId.begin();
1054   std::vector<const SMDS_MeshNode*>::const_iterator after  = theNodeByHybridId.end();
1055   
1056   (theNodeByHybridId.size() <= 1) ? tmpStr = " node" : " nodes";
1057   std::cout << theNodeByHybridId.size() << tmpStr << " from mesh ..." << std::endl;
1058   for ( ; hybridNodeIt != after; ++hybridNodeIt )
1059   {
1060     const SMDS_MeshNode* node = *hybridNodeIt;
1061     std::vector<double> coords;
1062     coords.push_back(node->X());
1063     coords.push_back(node->Y());
1064     coords.push_back(node->Z());
1065     nodesCoords.insert(coords);
1066     theOrderedNodes.push_back(node);
1067   }
1068   
1069   // Iterate over the enforced nodes given by enforced elements
1070   hybridNodeIt = theEnforcedNodeByHybridId.begin();
1071   after  = theEnforcedNodeByHybridId.end();
1072   (theEnforcedNodeByHybridId.size() <= 1) ? tmpStr = " node" : " nodes";
1073   std::cout << theEnforcedNodeByHybridId.size() << tmpStr << " from enforced elements ..." << std::endl;
1074   for ( ; hybridNodeIt != after; ++hybridNodeIt )
1075   {
1076     const SMDS_MeshNode* node = *hybridNodeIt;
1077     std::vector<double> coords;
1078     coords.push_back(node->X());
1079     coords.push_back(node->Y());
1080     coords.push_back(node->Z());
1081 #ifdef _DEBUG_
1082     std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1083 #endif
1084     
1085     if (nodesCoords.find(coords) != nodesCoords.end()) {
1086       // node already exists in original mesh
1087 #ifdef _DEBUG_
1088       std::cout << " found" << std::endl;
1089 #endif
1090       continue;
1091     }
1092     
1093     if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1094       // node already exists in enforced vertices
1095 #ifdef _DEBUG_
1096       std::cout << " found" << std::endl;
1097 #endif
1098       continue;
1099     }
1100     
1101 #ifdef _DEBUG_
1102     std::cout << " not found" << std::endl;
1103 #endif
1104     
1105     nodesCoords.insert(coords);
1106     theOrderedNodes.push_back(node);
1107 //     theRequiredNodes.push_back(node);
1108   }
1109   
1110   
1111   // Iterate over the enforced nodes
1112   HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap::const_iterator enfNodeIt;
1113   (theEnforcedNodes.size() <= 1) ? tmpStr = " node" : " nodes";
1114   std::cout << theEnforcedNodes.size() << tmpStr << " from enforced nodes ..." << std::endl;
1115   for(enfNodeIt = theEnforcedNodes.begin() ; enfNodeIt != theEnforcedNodes.end() ; ++enfNodeIt)
1116   {
1117     const SMDS_MeshNode* node = enfNodeIt->first;
1118     std::vector<double> coords;
1119     coords.push_back(node->X());
1120     coords.push_back(node->Y());
1121     coords.push_back(node->Z());
1122 #ifdef _DEBUG_
1123     std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1124 #endif
1125     
1126     // Test if point is inside shape to mesh
1127     gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1128     TopAbs_State result = pntCls->GetPointState( myPoint );
1129     if ( result == TopAbs_OUT ) {
1130 #ifdef _DEBUG_
1131       std::cout << " out of volume" << std::endl;
1132 #endif
1133       continue;
1134     }
1135     
1136     if (nodesCoords.find(coords) != nodesCoords.end()) {
1137 #ifdef _DEBUG_
1138       std::cout << " found in nodesCoords" << std::endl;
1139 #endif
1140 //       theRequiredNodes.push_back(node);
1141       continue;
1142     }
1143
1144     if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1145 #ifdef _DEBUG_
1146       std::cout << " found in theEnforcedVertices" << std::endl;
1147 #endif
1148       continue;
1149     }
1150     
1151 #ifdef _DEBUG_
1152     std::cout << " not found" << std::endl;
1153 #endif
1154     nodesCoords.insert(coords);
1155 //     theOrderedNodes.push_back(node);
1156     theRequiredNodes.push_back(node);
1157   }
1158   int requiredNodes = theRequiredNodes.size();
1159   
1160   int solSize = 0;
1161   std::vector<std::vector<double> > ReqVerTab;
1162   if (nbEnforcedVertices) {
1163     (nbEnforcedVertices <= 1) ? tmpStr = " node" : " nodes";
1164     std::cout << nbEnforcedVertices << tmpStr << " from enforced vertices ..." << std::endl;
1165     // Iterate over the enforced vertices
1166     for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
1167       double x = vertexIt->first[0];
1168       double y = vertexIt->first[1];
1169       double z = vertexIt->first[2];
1170       // Test if point is inside shape to mesh
1171       gp_Pnt myPoint(x,y,z);
1172       TopAbs_State result = pntCls->GetPointState( myPoint );
1173       if ( result == TopAbs_OUT )
1174         continue;
1175       std::vector<double> coords;
1176       coords.push_back(x);
1177       coords.push_back(y);
1178       coords.push_back(z);
1179       ReqVerTab.push_back(coords);
1180       enfVertexSizes.push_back(vertexIt->second);
1181       solSize++;
1182     }
1183   }
1184   
1185   
1186   // GmfVertices
1187   std::cout << "Begin writing required nodes in GmfVertices" << std::endl;
1188   std::cout << "Nb vertices: " << theOrderedNodes.size() << std::endl;
1189   MGInput->GmfSetKwd(idx, GmfVertices, theOrderedNodes.size());
1190   for (hybridNodeIt = theOrderedNodes.begin();hybridNodeIt != theOrderedNodes.end();++hybridNodeIt) {
1191     MGInput->GmfSetLin(idx, GmfVertices, (*hybridNodeIt)->X(), (*hybridNodeIt)->Y(), (*hybridNodeIt)->Z(), dummyint1);
1192   }
1193
1194   std::cout << "End writing required nodes in GmfVertices" << std::endl;
1195
1196   if (requiredNodes + solSize) {
1197     std::cout << "Begin writing in req and sol file" << std::endl;
1198     aNodeGroupByHybridId.resize( requiredNodes + solSize );
1199     idxRequired = MGInput->GmfOpenMesh(theRequiredFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1200     if (!idxRequired) {
1201       MGInput->GmfCloseMesh(idx);
1202       return false;
1203     }
1204     idxSol = MGInput->GmfOpenMesh(theSolFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1205     if (!idxSol) {
1206       MGInput->GmfCloseMesh(idx);
1207       if (idxRequired)
1208         MGInput->GmfCloseMesh(idxRequired);
1209       return false;
1210     }
1211     int TypTab[] = {GmfSca};
1212     double ValTab[] = {0.0};
1213     MGInput->GmfSetKwd(idxRequired, GmfVertices, requiredNodes + solSize);
1214     MGInput->GmfSetKwd(idxSol, GmfSolAtVertices, requiredNodes + solSize, 1, TypTab);
1215     for (hybridNodeIt = theRequiredNodes.begin();hybridNodeIt != theRequiredNodes.end();++hybridNodeIt) {
1216       MGInput->GmfSetLin(idxRequired, GmfVertices, (*hybridNodeIt)->X(), (*hybridNodeIt)->Y(), (*hybridNodeIt)->Z(), dummyint2);
1217       MGInput->GmfSetLin(idxSol, GmfSolAtVertices, ValTab);
1218       if (theEnforcedNodes.find((*hybridNodeIt)) != theEnforcedNodes.end())
1219         gn = theEnforcedNodes.find((*hybridNodeIt))->second;
1220       aNodeGroupByHybridId[usedEnforcedNodes] = gn;
1221       usedEnforcedNodes++;
1222     }
1223
1224     for (int i=0;i<solSize;i++) {
1225       std::cout << ReqVerTab[i][0] <<" "<< ReqVerTab[i][1] << " "<< ReqVerTab[i][2] << std::endl;
1226 #ifdef _DEBUG_
1227       std::cout << "enfVertexSizes.at("<<i<<"): " << enfVertexSizes.at(i) << std::endl;
1228 #endif
1229       double solTab[] = {enfVertexSizes.at(i)};
1230       MGInput->GmfSetLin(idxRequired, GmfVertices, ReqVerTab[i][0], ReqVerTab[i][1], ReqVerTab[i][2], dummyint3);
1231       MGInput->GmfSetLin(idxSol, GmfSolAtVertices, solTab);
1232       aNodeGroupByHybridId[usedEnforcedNodes] = enfVerticesWithGroup.find(ReqVerTab[i])->second;
1233 #ifdef _DEBUG_
1234       std::cout << "aNodeGroupByHybridId["<<usedEnforcedNodes<<"] = \""<<aNodeGroupByHybridId[usedEnforcedNodes]<<"\""<<std::endl;
1235 #endif
1236       usedEnforcedNodes++;
1237     }
1238     std::cout << "End writing in req and sol file" << std::endl;
1239   }
1240
1241   int nedge[2], ntri[3], nquad[4];
1242     
1243   // GmfEdges
1244   int usedEnforcedEdges = 0;
1245   if (theKeptEnforcedEdges.size()) {
1246     anEdgeGroupByHybridId.resize( theKeptEnforcedEdges.size() );
1247     MGInput->GmfSetKwd(idx, GmfEdges, theKeptEnforcedEdges.size());
1248     for(elemSetIt = theKeptEnforcedEdges.begin() ; elemSetIt != theKeptEnforcedEdges.end() ; ++elemSetIt) {
1249       elem = (*elemSetIt);
1250       nodeIt = elem->nodesIterator();
1251       int index=0;
1252       while ( nodeIt->more() ) {
1253         // find HYBRID ID
1254         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1255         std::map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToHybridIdMap.find(node);
1256         if (it == anEnforcedNodeToHybridIdMap.end()) {
1257           it = anExistingEnforcedNodeToHybridIdMap.find(node);
1258           if (it == anEnforcedNodeToHybridIdMap.end())
1259             throw "Node not found";
1260         }
1261         nedge[index] = it->second;
1262         index++;
1263       }
1264       MGInput->GmfSetLin(idx, GmfEdges, nedge[0], nedge[1], dummyint4);
1265       anEdgeGroupByHybridId[usedEnforcedEdges] = theEnforcedEdges.find(elem)->second;
1266       usedEnforcedEdges++;
1267     }
1268   }
1269
1270
1271   if (usedEnforcedEdges) {
1272     MGInput->GmfSetKwd(idx, GmfRequiredEdges, usedEnforcedEdges);
1273     for (int enfID=1;enfID<=usedEnforcedEdges;enfID++) {
1274       MGInput->GmfSetLin(idx, GmfRequiredEdges, enfID);
1275     }
1276   }
1277
1278   // GmfTriangles
1279   int usedEnforcedTriangles = 0;
1280   if (anElemSetTri.size()+theKeptEnforcedTriangles.size())
1281   {
1282     aFaceGroupByHybridId.resize( anElemSetTri.size()+theKeptEnforcedTriangles.size() );
1283     MGInput->GmfSetKwd(idx, GmfTriangles, anElemSetTri.size()+theKeptEnforcedTriangles.size());
1284     int k=0;
1285     for(elemSetIt = anElemSetTri.begin() ; elemSetIt != anElemSetTri.end() ; ++elemSetIt,++k)
1286     {
1287       elem = (*elemSetIt);
1288       theFaceByHybridId.push_back( elem );
1289       nodeIt = elem->nodesIterator();
1290       int index=0;
1291       for ( int j = 0; j < 3; ++j )
1292       {
1293         // find HYBRID ID
1294         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1295         std::map< const SMDS_MeshNode*,int >::iterator it = aNodeToHybridIdMap.find(node);
1296         if (it == aNodeToHybridIdMap.end())
1297           throw "Node not found";
1298         ntri[index] = it->second;
1299         index++;
1300       }
1301       MGInput->GmfSetLin(idx, GmfTriangles, ntri[0], ntri[1], ntri[2], /*tag=*/elem->getshapeId() );
1302       aFaceGroupByHybridId[k] = "";
1303     }
1304     
1305     if ( !theHelper.GetMesh()->HasShapeToMesh() ) SMESHUtils::FreeVector( theFaceByHybridId ); 
1306     std::cout << "Enforced triangles size " << theKeptEnforcedTriangles.size() << std::endl;
1307     if (theKeptEnforcedTriangles.size())
1308     {
1309       for(elemSetIt = theKeptEnforcedTriangles.begin() ; elemSetIt != theKeptEnforcedTriangles.end() ; ++elemSetIt,++k)
1310       {
1311         elem = (*elemSetIt);
1312         nodeIt = elem->nodesIterator();
1313         int index=0;
1314         for ( int j = 0; j < 3; ++j )
1315         {
1316           // find HYBRID ID
1317           const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1318           std::map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToHybridIdMap.find(node);
1319           if (it == anEnforcedNodeToHybridIdMap.end())
1320           {
1321             it = anExistingEnforcedNodeToHybridIdMap.find(node);
1322             if (it == anEnforcedNodeToHybridIdMap.end())
1323               throw "Node not found";
1324           }
1325           ntri[index] = it->second;
1326           index++;
1327         }
1328         MGInput->GmfSetLin(idx, GmfTriangles, ntri[0], ntri[1], ntri[2], enforcedTag);
1329         aFaceGroupByHybridId[k] = theEnforcedTriangles.find(elem)->second;
1330         usedEnforcedTriangles++;
1331       }
1332     }
1333   }
1334
1335   
1336   if (usedEnforcedTriangles)
1337   {
1338     MGInput->GmfSetKwd(idx, GmfRequiredTriangles, usedEnforcedTriangles);
1339     for (int enfID=1;enfID<=usedEnforcedTriangles;enfID++)
1340       MGInput->GmfSetLin(idx, GmfRequiredTriangles, anElemSetTri.size()+enfID);
1341   }
1342   
1343   if (anElemSetQuad.size())
1344   {
1345     MGInput->GmfSetKwd(idx, GmfQuadrilaterals, anElemSetQuad.size());
1346     int k=0;
1347     for(elemSetIt = anElemSetQuad.begin() ; elemSetIt != anElemSetQuad.end() ; ++elemSetIt,++k)
1348     {
1349       elem = (*elemSetIt);
1350       theFaceByHybridId.push_back( elem );
1351       nodeIt = elem->nodesIterator();
1352       int index=0;
1353       for ( int j = 0; j < 4; ++j )
1354       {
1355         // find HYBRID ID
1356         const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1357         std::map< const SMDS_MeshNode*,int >::iterator it = aNodeToHybridIdMap.find(node);
1358         if (it == aNodeToHybridIdMap.end())
1359           throw "Node not found";
1360         nquad[index] = it->second;
1361         index++;
1362       }
1363       MGInput->GmfSetLin(idx, GmfQuadrilaterals, nquad[0], nquad[1], nquad[2], nquad[3],
1364                          /*tag=*/elem->getshapeId() );
1365       // _CEA_cbo what is it for???
1366       //aFaceGroupByHybridId[k] = "";
1367     }
1368   }
1369
1370   MGInput->GmfCloseMesh(idx);
1371   if (idxRequired)
1372     MGInput->GmfCloseMesh(idxRequired);
1373   if (idxSol)
1374     MGInput->GmfCloseMesh(idxSol);
1375   
1376   return true;
1377   
1378 }
1379
1380 //=============================================================================
1381 /*!
1382  *Here we are going to use the HYBRID mesher with geometry
1383  */
1384 //=============================================================================
1385
1386 bool HYBRIDPlugin_HYBRID::Compute(SMESH_Mesh&         theMesh,
1387                                   const TopoDS_Shape& theShape)
1388 {
1389   bool Ok = false;
1390
1391   // a unique working file name
1392   // to avoid access to the same files by eg different users
1393   _genericName = HYBRIDPlugin_Hypothesis::GetFileName(_hyp);
1394   std::string aGenericName         = _genericName;
1395   std::string aGenericNameRequired = aGenericName + "_required";
1396
1397   std::string aLogFileName    = aGenericName + ".log";    // log
1398   std::string aResultFileName;
1399
1400   std::string aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1401   aGMFFileName              = aGenericName + ".mesh"; // GMF mesh file
1402   aResultFileName           = aGenericName + "Vol.mesh"; // GMF mesh file
1403   aResSolFileName           = aGenericName + "Vol.sol"; // GMF mesh file
1404   aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1405   aSolFileName              = aGenericNameRequired + ".sol"; // GMF solution file
1406   
1407   std::map <int,int> aNodeId2NodeIndexMap, aSmdsToHybridIdMap, anEnforcedNodeIdToHybridIdMap;
1408   std::map <int, int> nodeID2nodeIndexMap;
1409   std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1410   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues coordsSizeMap = HYBRIDPlugin_Hypothesis::GetEnforcedVerticesCoordsSize(_hyp);
1411   HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = HYBRIDPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1412   HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = HYBRIDPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1413   HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = HYBRIDPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1414   HYBRIDPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = HYBRIDPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1415
1416   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList enfVertices = HYBRIDPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1417   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1418   std::vector<double> coords;
1419
1420   for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1421   {
1422     HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertex* enfVertex = (*enfVerIt);
1423     if (enfVertex->coords.size()) {
1424       coordsSizeMap.insert(std::make_pair(enfVertex->coords,enfVertex->size));
1425       enfVerticesWithGroup.insert(std::make_pair(enfVertex->coords,enfVertex->groupName));
1426     }
1427     else {
1428       TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
1429       for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1430         coords.clear();
1431         if (it.Value().ShapeType() == TopAbs_VERTEX){
1432           gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1433           coords.push_back(aPnt.X());
1434           coords.push_back(aPnt.Y());
1435           coords.push_back(aPnt.Z());
1436           if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1437             coordsSizeMap.insert(std::make_pair(coords,enfVertex->size));
1438             enfVerticesWithGroup.insert(std::make_pair(coords,enfVertex->groupName));
1439           }
1440         }
1441       }
1442     }
1443   }
1444   int nbEnforcedVertices = coordsSizeMap.size();
1445   int nbEnforcedNodes = enforcedNodes.size();
1446
1447   std::string tmpStr;
1448   (nbEnforcedNodes <= 1) ? tmpStr = "node" : "nodes";
1449   std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1450   (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : "vertices";
1451   std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1452
1453   SMESH_MesherHelper helper( theMesh );
1454   helper.SetSubShape( theShape );
1455
1456   std::vector <const SMDS_MeshNode*> aNodeByHybridId, anEnforcedNodeByHybridId;
1457   std::vector <const SMDS_MeshElement*> aFaceByHybridId;
1458   std::map<const SMDS_MeshNode*,int> aNodeToHybridIdMap;
1459   std::vector<std::string> aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId;
1460
1461   SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1462
1463   MG_HYBRID_API mgHybrid( _computeCanceled, _progress );
1464
1465   Ok = writeGMFFile(&mgHybrid,
1466                     aGMFFileName.c_str(),
1467                     aRequiredVerticesFileName.c_str(),
1468                     aSolFileName.c_str(),
1469                     *proxyMesh, helper,
1470                     aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1471                     aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1472                     enforcedNodes, enforcedEdges, enforcedTriangles, /*enforcedQuadrangles,*/
1473                     enfVerticesWithGroup, coordsSizeMap);
1474
1475   // Write aSmdsToHybridIdMap to temp file
1476   std::string aSmdsToHybridIdMapFileName;
1477   aSmdsToHybridIdMapFileName = aGenericName + ".ids";  // ids relation
1478   ofstream aIdsFile  ( aSmdsToHybridIdMapFileName  , ios::out);
1479   Ok = aIdsFile.rdbuf()->is_open();
1480   if (!Ok) {
1481     INFOS( "Can't write into " << aSmdsToHybridIdMapFileName);
1482     return error(SMESH_Comment("Can't write into ") << aSmdsToHybridIdMapFileName);
1483   }
1484   INFOS( "Writing ids relation into " << aSmdsToHybridIdMapFileName);
1485   aIdsFile << "Smds Hybrid" << std::endl;
1486   std::map <int,int>::const_iterator myit;
1487   for (myit=aSmdsToHybridIdMap.begin() ; myit != aSmdsToHybridIdMap.end() ; ++myit) {
1488     aIdsFile << myit->first << " " << myit->second << std::endl;
1489   }
1490
1491   aIdsFile.close();
1492
1493   if ( ! Ok ) {
1494     if ( !_keepFiles ) {
1495       removeFile( aGMFFileName );
1496       removeFile( aRequiredVerticesFileName );
1497       removeFile( aSolFileName );
1498       removeFile( aSmdsToHybridIdMapFileName );
1499     }
1500     return error(COMPERR_BAD_INPUT_MESH);
1501   }
1502   removeFile( aResultFileName ); // needed for boundary recovery module usage
1503
1504   // -----------------
1505   // run hybrid mesher
1506   // -----------------
1507
1508   std::string cmd = HYBRIDPlugin_Hypothesis::CommandToRun( _hyp, theMesh );
1509
1510   if ( mgHybrid.IsExecutable() )
1511   {
1512     cmd += " --in "  + aGMFFileName;
1513     cmd += " --out " + aResultFileName;
1514   }
1515   std::cout << std::endl;
1516   std::cout << "Hybrid execution with geometry..." << std::endl;
1517   std::cout << cmd;
1518   if ( !_logInStandardOutput )
1519   {
1520     mgHybrid.SetLogFile( aLogFileName );
1521     if ( mgHybrid.IsExecutable() )
1522       cmd += " 1>" + aLogFileName;  // dump into file
1523     std::cout << " 1> " << aLogFileName;
1524   }
1525   std::cout << std::endl;
1526
1527   _computeCanceled = false;
1528
1529   std::string errStr;
1530   Ok = mgHybrid.Compute( cmd, errStr ); // run
1531
1532   if ( _logInStandardOutput && mgHybrid.IsLibrary() )
1533     std::cout << std::endl << mgHybrid.GetLog() << std::endl;
1534   if ( Ok )
1535     std::cout << "End of Hybrid execution !" << std::endl;
1536
1537   // --------------
1538   // read a result
1539   // --------------
1540
1541   // Mapping the result file
1542
1543   HYBRIDPlugin_Hypothesis::TSetStrings groupsToRemove = HYBRIDPlugin_Hypothesis::GetGroupsToRemove(_hyp);
1544   bool toMeshHoles =
1545     _hyp ? _hyp->GetToMeshHoles(true) : HYBRIDPlugin_Hypothesis::DefaultMeshHoles();
1546   const bool toMakeGroupsOfDomains = HYBRIDPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
1547
1548   helper.IsQuadraticSubMesh( theShape );
1549   helper.SetElementsOnShape( false );
1550
1551   Ok = readGMFFile(&mgHybrid, aResultFileName.c_str(),
1552                    this,
1553                    &helper, aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1554                    aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1555                    groupsToRemove, toMakeGroupsOfDomains, toMeshHoles);
1556
1557   removeEmptyGroupsOfDomains( helper.GetMesh(), !toMakeGroupsOfDomains );
1558
1559
1560
1561   // ---------------------
1562   // remove working files
1563   // ---------------------
1564
1565   if ( Ok )
1566   {
1567     if ( _removeLogOnSuccess )
1568       removeFile( aLogFileName );
1569   }
1570   else if ( mgHybrid.HasLog() )
1571   {
1572     // get problem description from the log file
1573     _Ghs2smdsConvertor conv( aNodeByHybridId );
1574     storeErrorDescription( _logInStandardOutput ? 0 : aLogFileName.c_str(),
1575                            mgHybrid.GetLog(), conv );
1576   }
1577   else if ( !errStr.empty() )
1578   {
1579     // the log file is empty
1580     removeFile( aLogFileName );
1581     INFOS( "HYBRID Error, " << errStr );
1582     error(COMPERR_ALGO_FAILED, errStr );
1583   }
1584
1585   if ( !_keepFiles ) {
1586     if (! Ok && _computeCanceled)
1587       removeFile( aLogFileName );
1588     removeFile( aGMFFileName );
1589     removeFile( aRequiredVerticesFileName );
1590     removeFile( aSolFileName );
1591     removeFile( aResSolFileName );
1592     removeFile( aResultFileName );
1593     removeFile( aSmdsToHybridIdMapFileName );
1594   }
1595   if ( mgHybrid.IsExecutable() )
1596   {
1597     std::cout << "<" << aResultFileName << "> HYBRID output file ";
1598     if ( !Ok )
1599       std::cout << "not ";
1600     std::cout << "treated !" << std::endl;
1601     std::cout << std::endl;
1602   }
1603   else
1604   {
1605     std::cout << "MG-HYBRID " << ( Ok ? "succeeded" : "failed") << std::endl;
1606   }
1607
1608   return Ok;
1609 }
1610
1611 //=============================================================================
1612 /*!
1613  *Here we are going to use the HYBRID mesher w/o geometry
1614  */
1615 //=============================================================================
1616 bool HYBRIDPlugin_HYBRID::Compute(SMESH_Mesh&         theMesh,
1617                                   SMESH_MesherHelper* theHelper)
1618 {
1619   theHelper->IsQuadraticSubMesh( theHelper->GetSubShape() );
1620
1621   // a unique working file name
1622   // to avoid access to the same files by eg different users
1623   _genericName = HYBRIDPlugin_Hypothesis::GetFileName(_hyp);
1624   std::string aGenericName((char*) _genericName.c_str() );
1625   std::string aGenericNameRequired = aGenericName + "_required";
1626
1627   std::string aLogFileName    = aGenericName + ".log";    // log
1628   std::string aResultFileName;
1629   bool Ok;
1630
1631   std::string aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1632   aGMFFileName              = aGenericName + ".mesh"; // GMF mesh file
1633   aResultFileName           = aGenericName + "Vol.mesh"; // GMF mesh file
1634   aResSolFileName           = aGenericName + "Vol.sol"; // GMF mesh file
1635   aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1636   aSolFileName              = aGenericNameRequired + ".sol"; // GMF solution file
1637
1638   std::map <int, int> nodeID2nodeIndexMap;
1639   std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1640   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues coordsSizeMap;
1641   TopoDS_Shape GeomShape;
1642   std::vector<double> coords;
1643   gp_Pnt aPnt;
1644   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertex* enfVertex;
1645
1646   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList enfVertices = HYBRIDPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1647   HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1648
1649   for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1650   {
1651     enfVertex = (*enfVerIt);
1652     if (enfVertex->coords.size()) {
1653       coordsSizeMap.insert(std::make_pair(enfVertex->coords,enfVertex->size));
1654       enfVerticesWithGroup.insert(std::make_pair(enfVertex->coords,enfVertex->groupName));
1655     }
1656     else {
1657       GeomShape = entryToShape(enfVertex->geomEntry);
1658       for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1659         coords.clear();
1660         if (it.Value().ShapeType() == TopAbs_VERTEX){
1661           aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1662           coords.push_back(aPnt.X());
1663           coords.push_back(aPnt.Y());
1664           coords.push_back(aPnt.Z());
1665           if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1666             coordsSizeMap.insert(std::make_pair(coords,enfVertex->size));
1667             enfVerticesWithGroup.insert(std::make_pair(coords,enfVertex->groupName));
1668           }
1669         }
1670       }
1671     }
1672   }
1673
1674   HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = HYBRIDPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1675   HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = HYBRIDPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1676   HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = HYBRIDPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1677   HYBRIDPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = HYBRIDPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1678
1679   std::string tmpStr;
1680
1681   int nbEnforcedVertices = coordsSizeMap.size();
1682   int nbEnforcedNodes = enforcedNodes.size();
1683   (nbEnforcedNodes <= 1) ? tmpStr = "node" : tmpStr = "nodes";
1684   std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1685   (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : tmpStr = "vertices";
1686   std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1687
1688   std::vector <const SMDS_MeshNode*> aNodeByHybridId, anEnforcedNodeByHybridId;
1689   std::vector <const SMDS_MeshElement*> aFaceByHybridId;
1690   std::map<const SMDS_MeshNode*,int> aNodeToHybridIdMap;
1691   std::vector<std::string> aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId;
1692
1693   SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1694
1695   MG_HYBRID_API mgHybrid( _computeCanceled, _progress );
1696
1697   Ok = writeGMFFile(&mgHybrid,
1698                     aGMFFileName.c_str(),
1699                     aRequiredVerticesFileName.c_str(), aSolFileName.c_str(),
1700                     *proxyMesh, *theHelper,
1701                     aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1702                     aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1703                     enforcedNodes, enforcedEdges, enforcedTriangles,
1704                     enfVerticesWithGroup, coordsSizeMap);
1705
1706   // -----------------
1707   // run hybrid mesher
1708   // -----------------
1709
1710   std::string cmd = HYBRIDPlugin_Hypothesis::CommandToRun( _hyp, theMesh );
1711
1712   if ( mgHybrid.IsExecutable() )
1713   {
1714     cmd += " --in "  + aGMFFileName;
1715     cmd += " --out " + aResultFileName;
1716   }
1717   if ( !_logInStandardOutput )
1718   {
1719     cmd += " 1> " + aLogFileName;  // dump into file
1720     mgHybrid.SetLogFile( aLogFileName );
1721   }
1722   std::cout << std::endl;
1723   std::cout << "Hybrid execution w/o geometry..." << std::endl;
1724   std::cout << cmd << std::endl;
1725
1726   _computeCanceled = false;
1727
1728   std::string errStr;
1729   Ok = mgHybrid.Compute( cmd, errStr ); // run
1730
1731   if ( _logInStandardOutput && mgHybrid.IsLibrary() )
1732     std::cout << std::endl << mgHybrid.GetLog() << std::endl;
1733   if ( Ok )
1734     std::cout << "End of Hybrid execution !" << std::endl;
1735
1736   // --------------
1737   // read a result
1738   // --------------
1739   HYBRIDPlugin_Hypothesis::TSetStrings groupsToRemove = HYBRIDPlugin_Hypothesis::GetGroupsToRemove(_hyp);
1740   const bool toMakeGroupsOfDomains = HYBRIDPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
1741
1742   Ok = Ok && readGMFFile(&mgHybrid,
1743                          aResultFileName.c_str(),
1744                          this,
1745                          theHelper, aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1746                          aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1747                          groupsToRemove, toMakeGroupsOfDomains);
1748
1749   updateMeshGroups(theHelper->GetMesh(), groupsToRemove);
1750   removeEmptyGroupsOfDomains( theHelper->GetMesh(), !toMakeGroupsOfDomains );
1751
1752   if ( Ok ) {
1753     HYBRIDPlugin_Hypothesis* that = (HYBRIDPlugin_Hypothesis*)this->_hyp;
1754     if (that)
1755       that->ClearGroupsToRemove();
1756   }
1757   // ---------------------
1758   // remove working files
1759   // ---------------------
1760
1761   if ( Ok )
1762   {
1763     if ( _removeLogOnSuccess )
1764       removeFile( aLogFileName );
1765   }
1766   else if ( mgHybrid.HasLog() )
1767   {
1768     // get problem description from the log file
1769     _Ghs2smdsConvertor conv( aNodeByHybridId );
1770     storeErrorDescription( _logInStandardOutput ? 0 : aLogFileName.c_str(),
1771                            mgHybrid.GetLog(), conv );
1772   }
1773   else {
1774     // the log file is empty
1775     removeFile( aLogFileName );
1776     INFOS( "HYBRID Error, command '" << cmd << "' failed" );
1777     error(COMPERR_ALGO_FAILED, "hybrid: command not found" );
1778   }
1779
1780   if ( !_keepFiles )
1781   {
1782     if (! Ok && _computeCanceled)
1783       removeFile( aLogFileName );
1784     removeFile( aGMFFileName );
1785     removeFile( aResultFileName );
1786     removeFile( aRequiredVerticesFileName );
1787     removeFile( aSolFileName );
1788     removeFile( aResSolFileName );
1789   }
1790   return Ok;
1791 }
1792
1793 void HYBRIDPlugin_HYBRID::CancelCompute()
1794 {
1795   _computeCanceled = true;
1796 #if !defined(WIN32) && !defined(__APPLE__)
1797   std::string cmd = "ps xo pid,args | grep " + _genericName;
1798   //cmd += " | grep -e \"^ *[0-9]\\+ \\+" + HYBRIDPlugin_Hypothesis::GetExeName() + "\"";
1799   cmd += " | awk '{print $1}' | xargs kill -9 > /dev/null 2>&1";
1800   system( cmd.c_str() );
1801 #endif
1802 }
1803
1804 //================================================================================
1805 /*!
1806  * \brief Provide human readable text by error code reported by hybrid
1807  */
1808 //================================================================================
1809
1810 static const char* translateError(const int errNum)
1811 {
1812   switch ( errNum ) {
1813   case 0:
1814     return "error distene 0";
1815   case 1:
1816     return "error distene 1";
1817   }
1818   return "unknown distene error";
1819 }
1820
1821 //================================================================================
1822 /*!
1823  * \brief Retrieve from a string given number of integers
1824  */
1825 //================================================================================
1826
1827 static char* getIds( char* ptr, int nbIds, std::vector<int>& ids )
1828 {
1829   ids.clear();
1830   ids.reserve( nbIds );
1831   while ( nbIds )
1832   {
1833     while ( !isdigit( *ptr )) ++ptr;
1834     if ( ptr[-1] == '-' ) --ptr;
1835     ids.push_back( strtol( ptr, &ptr, 10 ));
1836     --nbIds;
1837   }
1838   return ptr;
1839 }
1840
1841 //================================================================================
1842 /*!
1843  * \brief Retrieve problem description form a log file
1844  *  \retval bool - always false
1845  */
1846 //================================================================================
1847
1848 bool HYBRIDPlugin_HYBRID::storeErrorDescription(const char*                logFile,
1849                                                 const std::string&         log,
1850                                                 const _Ghs2smdsConvertor & toSmdsConvertor )
1851 {
1852   if(_computeCanceled)
1853     return error(SMESH_Comment("interruption initiated by user"));
1854
1855   char* ptr = const_cast<char*>( log.c_str() );
1856   char* buf = ptr, * bufEnd = ptr + log.size();
1857
1858   SMESH_Comment errDescription;
1859
1860   enum { NODE = 1, EDGE, TRIA, VOL, SKIP_ID = 1 };
1861
1862   // look for MeshGems version
1863   // Since "MG-TETRA -- MeshGems 1.1-3 (January, 2013)" error codes change.
1864   // To discriminate old codes from new ones we add 1000000 to the new codes.
1865   // This way value of the new codes is same as absolute value of codes printed
1866   // in the log after "MGMESSAGE" string.
1867   int versionAddition = 0;
1868   {
1869     char* verPtr = ptr;
1870     while ( ++verPtr < bufEnd )
1871     {
1872       if ( strncmp( verPtr, "MG-TETRA -- MeshGems ", 21 ) != 0 )
1873         continue;
1874       if ( strcmp( verPtr, "MG-TETRA -- MeshGems 1.1-3 " ) >= 0 )
1875         versionAddition = 1000000;
1876       ptr = verPtr;
1877       break;
1878     }
1879   }
1880
1881   // look for errors "ERR #"
1882
1883   std::set<std::string> foundErrorStr; // to avoid reporting same error several times
1884   std::set<int>         elemErrorNums; // not to report different types of errors with bad elements
1885   while ( ++ptr < bufEnd )
1886   {
1887     if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1888       continue;
1889
1890     std::list<const SMDS_MeshElement*> badElems;
1891     std::vector<int> nodeIds;
1892
1893     ptr += 4;
1894     char* errBeg = ptr;
1895     int   errNum = strtol(ptr, &ptr, 10) + versionAddition;
1896     // we treat errors enumerated in [SALOME platform 0019316] issue
1897     // and all errors from a new (Release 1.1) MeshGems User Manual
1898     switch ( errNum ) {
1899     case 0015: // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1900     case 1005620 : // a too bad quality face is detected. This face is considered degenerated.
1901       ptr = getIds(ptr, SKIP_ID, nodeIds);
1902       ptr = getIds(ptr, TRIA, nodeIds);
1903       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1904       break;
1905     case 1005621 : // a too bad quality face is detected. This face is degenerated.
1906       // hence the is degenerated it is invisible, add its edges in addition
1907       ptr = getIds(ptr, SKIP_ID, nodeIds);
1908       ptr = getIds(ptr, TRIA, nodeIds);
1909       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1910       {
1911         std::vector<int> edgeNodes( nodeIds.begin(), --nodeIds.end() ); // 01
1912         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1913         edgeNodes[1] = nodeIds[2]; // 02
1914         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1915         edgeNodes[0] = nodeIds[1]; // 12
1916       }      
1917       break;
1918     case 1000: // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1919       // ERR  1000 :  1 3 2
1920     case 1002: // Face (f 1, f 2, f 3) has a vertex negative or null
1921     case 3019: // Constrained face (f 1, f 2, f 3) cannot be enforced
1922     case 1002211: // a face has a vertex negative or null.
1923     case 1005200 : // a surface mesh appears more than once in the input surface mesh.
1924     case 1008423 : // a constrained face cannot be enforced (regeneration phase failed).
1925       ptr = getIds(ptr, TRIA, nodeIds);
1926       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1927       break;
1928     case 1001: // Edge (e1, e2) appears more than once in the input surface mesh
1929     case 3009: // Constrained edge (e1, e2) cannot be enforced (warning).
1930       // ERR  3109 :  EDGE  5 6 UNIQUE
1931     case 3109: // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1932     case 1005210 : // an edge appears more than once in the input surface mesh.
1933     case 1005820 : // an edge is unique (i.e., bounds a hole in the surface).
1934     case 1008441 : // a constrained edge cannot be enforced.
1935       ptr = getIds(ptr, EDGE, nodeIds);
1936       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1937       break;
1938     case 2004: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1939     case 2014: // at least two points whose distance is dist, i.e., considered as coincident
1940     case 2103: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1941       // ERR  2103 :  16 WITH  3
1942     case 1005105 : // two vertices are too close to one another or coincident.
1943     case 1005107: // Two vertices are too close to one another or coincident.
1944       ptr = getIds(ptr, NODE, nodeIds);
1945       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1946       ptr = getIds(ptr, NODE, nodeIds);
1947       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1948       break;
1949     case 2012: // Vertex v1 cannot be inserted (warning).
1950     case 1005106 : // a vertex cannot be inserted.
1951       ptr = getIds(ptr, NODE, nodeIds);
1952       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1953       break;
1954     case 3103: // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1955     case 1005110 : // two surface edges are intersecting.
1956       // ERR  3103 :  1 2 WITH  7 3
1957       ptr = getIds(ptr, EDGE, nodeIds);
1958       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1959       ptr = getIds(ptr, EDGE, nodeIds);
1960       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1961       break;
1962     case 3104: // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1963       // ERR  3104 :  9 10 WITH  1 2 3
1964     case 3106: // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1965     case 1005120 : // a surface edge intersects a surface face.
1966       ptr = getIds(ptr, EDGE, nodeIds);
1967       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1968       ptr = getIds(ptr, TRIA, nodeIds);
1969       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1970       break;
1971     case 3105: // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1972       // ERR  3105 :  8 IN  2 3 5
1973     case 1005150 : // a boundary point lies within a surface face.
1974       ptr = getIds(ptr, NODE, nodeIds);
1975       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1976       ptr = getIds(ptr, TRIA, nodeIds);
1977       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1978       break;
1979     case 3107: // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1980       // ERR  3107 :  2 IN  4 1
1981     case 1005160 : // a boundary point lies within a surface edge.
1982       ptr = getIds(ptr, NODE, nodeIds);
1983       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1984       ptr = getIds(ptr, EDGE, nodeIds);
1985       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1986       break;
1987     case 9000: // ERR  9000
1988       //  ELEMENT  261 WITH VERTICES :  7 396 -8 242
1989       //  VOLUME   : -1.11325045E+11 W.R.T. EPSILON   0.
1990       // A too small volume element is detected. Are reported the index of the element,
1991       // its four vertex indices, its volume and the tolerance threshold value
1992       ptr = getIds(ptr, SKIP_ID, nodeIds);
1993       ptr = getIds(ptr, VOL, nodeIds);
1994       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1995       // even if all nodes found, volume it most probably invisible,
1996       // add its faces to demonstrate it anyhow
1997       {
1998         std::vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1999         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2000         faceNodes[2] = nodeIds[3]; // 013
2001         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2002         faceNodes[1] = nodeIds[2]; // 023
2003         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2004         faceNodes[0] = nodeIds[1]; // 123
2005         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2006       }
2007       break;
2008     case 9001: // ERR  9001
2009       //  %% NUMBER OF NEGATIVE VOLUME TETS  :  1
2010       //  %% THE LARGEST NEGATIVE TET        :   1.75376581E+11
2011       //  %%  NUMBER OF NULL VOLUME TETS     :  0
2012       // There exists at least a null or negative volume element
2013       break;
2014     case 9002:
2015       // There exist n null or negative volume elements
2016       break;
2017     case 9003:
2018       // A too small volume element is detected
2019       break;
2020     case 9102:
2021       // A too bad quality face is detected. This face is considered degenerated,
2022       // its index, its three vertex indices together with its quality value are reported
2023       break; // same as next
2024     case 9112: // ERR  9112
2025       //  FACE   2 WITH VERTICES :  4 2 5
2026       //  SMALL INRADIUS :   0.
2027       // A too bad quality face is detected. This face is degenerated,
2028       // its index, its three vertex indices together with its inradius are reported
2029       ptr = getIds(ptr, SKIP_ID, nodeIds);
2030       ptr = getIds(ptr, TRIA, nodeIds);
2031       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2032       // add triangle edges as it most probably has zero area and hence invisible
2033       {
2034         std::vector<int> edgeNodes(2);
2035         edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
2036         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2037         edgeNodes[1] = nodeIds[2]; // 0-2
2038         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2039         edgeNodes[0] = nodeIds[1]; // 1-2
2040         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2041       }
2042       break;
2043     case 1005103 : // the vertices of an element are too close to one another or coincident.
2044       ptr = getIds(ptr, TRIA, nodeIds);
2045       if ( nodeIds.back() == 0 ) // index of the third vertex of the element (0 for an edge)
2046         nodeIds.resize( EDGE );
2047       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2048       break;
2049     }
2050
2051     bool isNewError = foundErrorStr.insert( std::string( errBeg, ptr )).second;
2052     if ( !isNewError )
2053       continue; // not to report same error several times
2054
2055 //     const SMDS_MeshElement* nullElem = 0;
2056 //     bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
2057
2058 //     if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
2059 //       bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
2060 //       if ( oneMoreErrorType )
2061 //         continue; // not to report different types of errors with bad elements
2062 //     }
2063
2064     // store bad elements
2065     //if ( allElemsOk ) {
2066       std::list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
2067       for ( ; elem != badElems.end(); ++elem )
2068         addBadInputElement( *elem );
2069       //}
2070
2071     // make error text
2072     std::string text = translateError( errNum );
2073     if ( errDescription.find( text ) == text.npos ) {
2074       if ( !errDescription.empty() )
2075         errDescription << "\n";
2076       errDescription << text;
2077     }
2078
2079   } // end while
2080
2081   if ( errDescription.empty() ) { // no errors found
2082     char msgLic1[] = "connection to server failed";
2083     char msgLic2[] = " Dlim ";
2084     if ( std::search( &buf[0], bufEnd, msgLic1, msgLic1 + strlen(msgLic1)) != bufEnd ||
2085          std::search( &buf[0], bufEnd, msgLic2, msgLic2 + strlen(msgLic2)) != bufEnd )
2086       errDescription << "Licence problems.";
2087     else
2088     {
2089       char msg2[] = "SEGMENTATION FAULT";
2090       if ( std::search( &buf[0], bufEnd, msg2, msg2 + strlen(msg2)) != bufEnd )
2091         errDescription << "hybrid: SEGMENTATION FAULT. ";
2092     }
2093   }
2094
2095   if ( logFile && logFile[0] )
2096   {
2097     if ( errDescription.empty() )
2098       errDescription << "See " << logFile << " for problem description";
2099     else
2100       errDescription << "\nSee " << logFile << " for more information";
2101   }
2102   return error( errDescription );
2103 }
2104
2105 //================================================================================
2106 /*!
2107  * \brief Creates _Ghs2smdsConvertor
2108  */
2109 //================================================================================
2110
2111 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const std::map <int,const SMDS_MeshNode*> & ghs2NodeMap)
2112   :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
2113 {
2114 }
2115
2116 //================================================================================
2117 /*!
2118  * \brief Creates _Ghs2smdsConvertor
2119  */
2120 //================================================================================
2121
2122 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const std::vector <const SMDS_MeshNode*> &  nodeByGhsId)
2123   : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
2124 {
2125 }
2126
2127 //================================================================================
2128 /*!
2129  * \brief Return SMDS element by ids of HYBRID nodes
2130  */
2131 //================================================================================
2132
2133 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const std::vector<int>& ghsNodes) const
2134 {
2135   size_t nbNodes = ghsNodes.size();
2136   std::vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
2137   for ( size_t i = 0; i < nbNodes; ++i ) {
2138     int ghsNode = ghsNodes[ i ];
2139     if ( _ghs2NodeMap ) {
2140       std::map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
2141       if ( in == _ghs2NodeMap->end() )
2142         return 0;
2143       nodes[ i ] = in->second;
2144     }
2145     else {
2146       if ( ghsNode < 1 || ghsNode > (int)_nodeByGhsId->size() )
2147         return 0;
2148       nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
2149     }
2150   }
2151   if ( nbNodes == 1 )
2152     return nodes[0];
2153
2154   if ( nbNodes == 2 ) {
2155     const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
2156     if ( !edge )
2157       edge = new SMDS_LinearEdge( nodes[0], nodes[1] );
2158     return edge;
2159   }
2160   if ( nbNodes == 3 ) {
2161     const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
2162     if ( !face )
2163       face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
2164     return face;
2165   }
2166   if ( nbNodes == 4 )
2167     return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
2168
2169   return 0;
2170 }
2171
2172
2173 //=============================================================================
2174 /*!
2175  *
2176  */
2177 //=============================================================================
2178 bool HYBRIDPlugin_HYBRID::Evaluate(SMESH_Mesh& aMesh,
2179                                  const TopoDS_Shape& aShape,
2180                                  MapShapeNbElems& aResMap)
2181 {
2182   int nbtri = 0, nbqua = 0;
2183   double fullArea = 0.0;
2184   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2185     TopoDS_Face F = TopoDS::Face( exp.Current() );
2186     SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2187     MapShapeNbElemsItr anIt = aResMap.find(sm);
2188     if( anIt==aResMap.end() ) {
2189       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2190       smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
2191                                             "Submesh can not be evaluated",this));
2192       return false;
2193     }
2194     std::vector<int> aVec = (*anIt).second;
2195     nbtri += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
2196     nbqua += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
2197     GProp_GProps G;
2198     BRepGProp::SurfaceProperties(F,G);
2199     double anArea = G.Mass();
2200     fullArea += anArea;
2201   }
2202
2203   // collect info from edges
2204   int nb0d_e = 0, nb1d_e = 0;
2205   bool IsQuadratic = false;
2206   bool IsFirst = true;
2207   TopTools_MapOfShape tmpMap;
2208   for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2209     TopoDS_Edge E = TopoDS::Edge(exp.Current());
2210     if( tmpMap.Contains(E) )
2211       continue;
2212     tmpMap.Add(E);
2213     SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
2214     MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
2215     std::vector<int> aVec = (*anIt).second;
2216     nb0d_e += aVec[SMDSEntity_Node];
2217     nb1d_e += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
2218     if(IsFirst) {
2219       IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
2220       IsFirst = false;
2221     }
2222   }
2223   tmpMap.Clear();
2224
2225   double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
2226
2227   GProp_GProps G;
2228   BRepGProp::VolumeProperties(aShape,G);
2229   double aVolume = G.Mass();
2230   double tetrVol = 0.1179*ELen*ELen*ELen;
2231   double CoeffQuality = 0.9;
2232   int nbVols = int(aVolume/tetrVol/CoeffQuality);
2233   int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2234   int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2235   std::vector<int> aVec(SMDSEntity_Last);
2236   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2237   if( IsQuadratic ) {
2238     aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2239     aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2240     aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2241   }
2242   else {
2243     aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2244     aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2245     aVec[SMDSEntity_Pyramid] = nbqua;
2246   }
2247   SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2248   aResMap.insert(std::make_pair(sm,aVec));
2249
2250   return true;
2251 }
2252
2253 bool HYBRIDPlugin_HYBRID::importGMFMesh(const char* theGMFFileName, SMESH_Mesh& theMesh)
2254 {
2255   SMESH_ComputeErrorPtr err = theMesh.GMFToMesh( theGMFFileName, /*makeRequiredGroups =*/ true );
2256
2257   theMesh.GetMeshDS()->Modified();
2258
2259   return ( !err || err->IsOK());
2260 }
2261
2262 namespace
2263 {
2264   //================================================================================
2265   /*!
2266    * \brief Sub-mesh event listener setting enforced elements as soon as an enforced
2267    *        mesh is loaded
2268    */
2269   struct _EnforcedMeshRestorer : public SMESH_subMeshEventListener
2270   {
2271     _EnforcedMeshRestorer():
2272       SMESH_subMeshEventListener( /*isDeletable = */true, Name() )
2273     {}
2274
2275     //================================================================================
2276     /*!
2277      * \brief Returns an ID of listener
2278      */
2279     static const char* Name() { return "HYBRIDPlugin_HYBRID::_EnforcedMeshRestorer"; }
2280
2281     //================================================================================
2282     /*!
2283      * \brief Treat events of the subMesh
2284      */
2285     void ProcessEvent(const int                       event,
2286                       const int                       eventType,
2287                       SMESH_subMesh*                  subMesh,
2288                       SMESH_subMeshEventListenerData* data,
2289                       const SMESH_Hypothesis*         hyp)
2290     {
2291       if ( SMESH_subMesh::SUBMESH_LOADED == event &&
2292            SMESH_subMesh::COMPUTE_EVENT  == eventType &&
2293            data &&
2294            !data->mySubMeshes.empty() )
2295       {
2296         // An enforced mesh (subMesh->_father) has been loaded from hdf file
2297         if ( HYBRIDPlugin_Hypothesis* hyp = GetGHSHypothesis( data->mySubMeshes.front() ))
2298           hyp->RestoreEnfElemsByMeshes();
2299       }
2300     }
2301     //================================================================================
2302     /*!
2303      * \brief Returns HYBRIDPlugin_Hypothesis used to compute a subMesh
2304      */
2305     static HYBRIDPlugin_Hypothesis* GetGHSHypothesis( SMESH_subMesh* subMesh )
2306     {
2307       SMESH_HypoFilter ghsHypFilter( SMESH_HypoFilter::HasName( "HYBRID_Parameters" ));
2308       return (HYBRIDPlugin_Hypothesis* )
2309         subMesh->GetFather()->GetHypothesis( subMesh->GetSubShape(),
2310                                              ghsHypFilter,
2311                                              /*visitAncestors=*/true);
2312     }
2313   };
2314
2315   //================================================================================
2316   /*!
2317    * \brief Sub-mesh event listener removing empty groups created due to "To make
2318    *        groups of domains".
2319    */
2320   struct _GroupsOfDomainsRemover : public SMESH_subMeshEventListener
2321   {
2322     _GroupsOfDomainsRemover():
2323       SMESH_subMeshEventListener( /*isDeletable = */true,
2324                                   "HYBRIDPlugin_HYBRID::_GroupsOfDomainsRemover" ) {}
2325     /*!
2326      * \brief Treat events of the subMesh
2327      */
2328     void ProcessEvent(const int                       event,
2329                       const int                       eventType,
2330                       SMESH_subMesh*                  subMesh,
2331                       SMESH_subMeshEventListenerData* data,
2332                       const SMESH_Hypothesis*         hyp)
2333     {
2334       if (SMESH_subMesh::ALGO_EVENT == eventType &&
2335           !subMesh->GetAlgo() )
2336       {
2337         removeEmptyGroupsOfDomains( subMesh->GetFather(), /*notEmptyAsWell=*/true );
2338       }
2339     }
2340   };
2341 }
2342
2343 //================================================================================
2344 /*!
2345  * \brief Set an event listener to set enforced elements as soon as an enforced
2346  *        mesh is loaded
2347  */
2348 //================================================================================
2349
2350 void HYBRIDPlugin_HYBRID::SubmeshRestored(SMESH_subMesh* subMesh)
2351 {
2352   if ( HYBRIDPlugin_Hypothesis* hyp = _EnforcedMeshRestorer::GetGHSHypothesis( subMesh ))
2353   {
2354     HYBRIDPlugin_Hypothesis::THYBRIDEnforcedMeshList enfMeshes = hyp->_GetEnforcedMeshes();
2355     HYBRIDPlugin_Hypothesis::THYBRIDEnforcedMeshList::iterator it = enfMeshes.begin();
2356     for(;it != enfMeshes.end();++it) {
2357       HYBRIDPlugin_Hypothesis::THYBRIDEnforcedMesh* enfMesh = *it;
2358       if ( SMESH_Mesh* mesh = GetMeshByPersistentID( enfMesh->persistID ))
2359       {
2360         SMESH_subMesh* smToListen = mesh->GetSubMesh( mesh->GetShapeToMesh() );
2361         // a listener set to smToListen will care of hypothesis stored in SMESH_EventListenerData
2362         subMesh->SetEventListener( new _EnforcedMeshRestorer(),
2363                                    SMESH_subMeshEventListenerData::MakeData( subMesh ),
2364                                    smToListen);
2365       }
2366     }
2367   }
2368 }
2369
2370 //================================================================================
2371 /*!
2372  * \brief Sets an event listener removing empty groups created due to "To make
2373  *        groups of domains".
2374  * \param subMesh - submesh where algo is set
2375  *
2376  * This method is called when a submesh gets HYP_OK algo_state.
2377  * After being set, event listener is notified on each event of a submesh.
2378  */
2379 //================================================================================
2380
2381 void HYBRIDPlugin_HYBRID::SetEventListener(SMESH_subMesh* subMesh)
2382 {
2383   subMesh->SetEventListener( new _GroupsOfDomainsRemover(), 0, subMesh );
2384 }