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