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