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