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