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