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