Salome HOME
Merge from V6_main_20120808 08Aug12
[plugins/blsurfplugin.git] / src / BLSURFPlugin / BLSURFPlugin_BLSURF.cxx
1 // Copyright (C) 2007-2012  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.
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 // File    : BLSURFPlugin_BLSURF.cxx
21 // Authors : Francis KLOSS (OCC) & Patrick LAUG (INRIA) & Lioka RAZAFINDRAZAKA (CEA)
22 //           & Aurelien ALLEAUME (DISTENE)
23 //           Size maps developement: Nicolas GEIMER (OCC) & Gilles DAVID (EURIWARE)
24 // ---
25
26 #include "BLSURFPlugin_BLSURF.hxx"
27 #include "BLSURFPlugin_Hypothesis.hxx"
28 #include "BLSURFPlugin_Attractor.hxx"
29
30 extern "C"{
31 #include <distene/api.h>
32 #include <distene/blsurf.h>
33 #include <distene/precad.h>
34 }
35
36 #include <structmember.h>
37
38
39 #include <Basics_Utils.hxx>
40 #include <Basics_OCCTVersion.hxx>
41
42 #include <SMESHDS_Group.hxx>
43 #include <SMESH_Gen.hxx>
44 #include <SMESH_Group.hxx>
45 #include <SMESH_Mesh.hxx>
46 #include <SMESH_MeshEditor.hxx>
47 #include <SMESH_MesherHelper.hxx>
48 #include <StdMeshers_FaceSide.hxx>
49
50 #include <utilities.h>
51
52 #include <limits>
53 #include <list>
54 #include <vector>
55 #include <set>
56 #include <cstdlib>
57
58 // OPENCASCADE includes
59 #include <BRep_Tool.hxx>
60 #include <TopExp.hxx>
61 #include <TopExp_Explorer.hxx>
62 #include <TopoDS.hxx>
63 #include <NCollection_Map.hxx>
64
65 #include <Geom_Surface.hxx>
66 #include <Handle_Geom_Surface.hxx>
67 #include <Geom2d_Curve.hxx>
68 #include <Handle_Geom2d_Curve.hxx>
69 #include <Geom_Curve.hxx>
70 #include <Handle_Geom_Curve.hxx>
71 #include <Handle_AIS_InteractiveObject.hxx>
72 #include <TopoDS_Vertex.hxx>
73 #include <TopoDS_Edge.hxx>
74 #include <TopoDS_Wire.hxx>
75 #include <TopoDS_Face.hxx>
76
77 #include <gp_Pnt2d.hxx>
78 #include <gp_Pnt.hxx>
79 #include <TopTools_IndexedMapOfShape.hxx>
80 #include <TopoDS_Shape.hxx>
81 #include <BRep_Builder.hxx>
82 #include <BRepBuilderAPI_MakeVertex.hxx>
83 #include <BRepTools.hxx>
84
85 #include <TopTools_DataMapOfShapeInteger.hxx>
86 #include <GProp_GProps.hxx>
87 #include <BRepGProp.hxx>
88
89 #ifndef WNT
90 #include <fenv.h>
91 #endif
92
93 #include <Standard_ErrorHandler.hxx>
94 #include <GeomAPI_ProjectPointOnCurve.hxx>
95 #include <GeomAPI_ProjectPointOnSurf.hxx>
96 #include <gp_XY.hxx>
97 #include <gp_XYZ.hxx>
98 // #include <BRepClass_FaceClassifier.hxx>
99 #include <TopTools_MapOfShape.hxx>
100
101 /* ==================================
102  * ===========  PYTHON ==============
103  * ==================================*/
104
105 typedef struct {
106   PyObject_HEAD
107   int softspace;
108   std::string *out;
109   } PyStdOut;
110
111 static void
112 PyStdOut_dealloc(PyStdOut *self)
113 {
114   PyObject_Del(self);
115 }
116
117 static PyObject *
118 PyStdOut_write(PyStdOut *self, PyObject *args)
119 {
120   char *c;
121   int l;
122   if (!PyArg_ParseTuple(args, "t#:write",&c, &l))
123     return NULL;
124
125   //std::cerr << c ;
126   *(self->out)=*(self->out)+c;
127
128   Py_INCREF(Py_None);
129   return Py_None;
130 }
131
132 static PyMethodDef PyStdOut_methods[] = {
133   {"write",  (PyCFunction)PyStdOut_write,  METH_VARARGS,
134     PyDoc_STR("write(string) -> None")},
135   {NULL,    NULL}   /* sentinel */
136 };
137
138 static PyMemberDef PyStdOut_memberlist[] = {
139   {(char*)"softspace", T_INT,  offsetof(PyStdOut, softspace), 0,
140    (char*)"flag indicating that a space needs to be printed; used by print"},
141   {NULL} /* Sentinel */
142 };
143
144 static PyTypeObject PyStdOut_Type = {
145   /* The ob_type field must be initialized in the module init function
146    * to be portable to Windows without using C++. */
147   PyObject_HEAD_INIT(NULL)
148   0,                            /*ob_size*/
149   "PyOut",                      /*tp_name*/
150   sizeof(PyStdOut),             /*tp_basicsize*/
151   0,                            /*tp_itemsize*/
152   /* methods */
153   (destructor)PyStdOut_dealloc, /*tp_dealloc*/
154   0,                            /*tp_print*/
155   0,                            /*tp_getattr*/
156   0,                            /*tp_setattr*/
157   0,                            /*tp_compare*/
158   0,                            /*tp_repr*/
159   0,                            /*tp_as_number*/
160   0,                            /*tp_as_sequence*/
161   0,                            /*tp_as_mapping*/
162   0,                            /*tp_hash*/
163   0,                            /*tp_call*/
164   0,                            /*tp_str*/
165   PyObject_GenericGetAttr,      /*tp_getattro*/
166   /* softspace is writable:  we must supply tp_setattro */
167   PyObject_GenericSetAttr,      /* tp_setattro */
168   0,                            /*tp_as_buffer*/
169   Py_TPFLAGS_DEFAULT,           /*tp_flags*/
170   0,                            /*tp_doc*/
171   0,                            /*tp_traverse*/
172   0,                            /*tp_clear*/
173   0,                            /*tp_richcompare*/
174   0,                            /*tp_weaklistoffset*/
175   0,                            /*tp_iter*/
176   0,                            /*tp_iternext*/
177   PyStdOut_methods,             /*tp_methods*/
178   PyStdOut_memberlist,          /*tp_members*/
179   0,                            /*tp_getset*/
180   0,                            /*tp_base*/
181   0,                            /*tp_dict*/
182   0,                            /*tp_descr_get*/
183   0,                            /*tp_descr_set*/
184   0,                            /*tp_dictoffset*/
185   0,                            /*tp_init*/
186   0,                            /*tp_alloc*/
187   0,                            /*tp_new*/
188   0,                            /*tp_free*/
189   0,                            /*tp_is_gc*/
190 };
191
192 PyObject * newPyStdOut( std::string& out )
193 {
194   PyStdOut *self;
195   self = PyObject_New(PyStdOut, &PyStdOut_Type);
196   if (self == NULL)
197     return NULL;
198   self->softspace = 0;
199   self->out=&out;
200   return (PyObject*)self;
201 }
202
203
204 ////////////////////////END PYTHON///////////////////////////
205
206 //////////////////MY MAPS////////////////////////////////////////
207 TopTools_IndexedMapOfShape FacesWithSizeMap;
208 std::map<int,string> FaceId2SizeMap;
209 TopTools_IndexedMapOfShape EdgesWithSizeMap;
210 std::map<int,string> EdgeId2SizeMap;
211 TopTools_IndexedMapOfShape VerticesWithSizeMap;
212 std::map<int,string> VertexId2SizeMap;
213
214 std::map<int,PyObject*> FaceId2PythonSmp;
215 std::map<int,PyObject*> EdgeId2PythonSmp;
216 std::map<int,PyObject*> VertexId2PythonSmp;
217
218 std::map<int,std::vector<double> > FaceId2AttractorCoords;
219 std::map<int,BLSURFPlugin_Attractor*> FaceId2ClassAttractor;
220 std::map<int,BLSURFPlugin_Attractor*> FaceIndex2ClassAttractor;
221
222 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
223 std::map< int, BLSURFPlugin_Hypothesis::TEnfVertexCoordsList > FaceId2EnforcedVertexCoords;
224 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexCoords > EnfVertexCoords2ProjVertex;
225 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList > EnfVertexCoords2EnfVertexList;
226
227 bool HasSizeMapOnFace=false;
228 bool HasSizeMapOnEdge=false;
229 bool HasSizeMapOnVertex=false;
230 //bool HasAttractorOnFace=false;
231
232 //=============================================================================
233 /*!
234  *
235  */
236 //=============================================================================
237
238 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId, int studyId,
239                                                SMESH_Gen* gen)
240   : SMESH_2D_Algo(hypId, studyId, gen)
241 {
242   MESSAGE("BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF");
243
244   _name = "BLSURF";
245   _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
246   _compatibleHypothesis.push_back(BLSURFPlugin_Hypothesis::GetHypType());
247   _requireDiscreteBoundary = false;
248   _onlyUnaryInput = false;
249   _hypothesis = NULL;
250   _supportSubmeshes = true;
251
252   smeshGen_i = SMESH_Gen_i::GetSMESHGen();
253   CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
254   SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
255
256   MESSAGE("studyid = " << _studyId);
257
258   myStudy = NULL;
259   myStudy = aStudyMgr->GetStudyByID(_studyId);
260   if (myStudy)
261     MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
262
263   /* Initialize the Python interpreter */
264   assert(Py_IsInitialized());
265   PyGILState_STATE gstate;
266   gstate = PyGILState_Ensure();
267
268   main_mod = NULL;
269   main_mod = PyImport_AddModule("__main__");
270
271   main_dict = NULL;
272   main_dict = PyModule_GetDict(main_mod);
273
274   PyRun_SimpleString("from math import *");
275   PyGILState_Release(gstate);
276
277   FacesWithSizeMap.Clear();
278   FaceId2SizeMap.clear();
279   EdgesWithSizeMap.Clear();
280   EdgeId2SizeMap.clear();
281   VerticesWithSizeMap.Clear();
282   VertexId2SizeMap.clear();
283   FaceId2PythonSmp.clear();
284   EdgeId2PythonSmp.clear();
285   VertexId2PythonSmp.clear();
286   FaceId2AttractorCoords.clear();
287   FaceId2ClassAttractor.clear();
288   FaceIndex2ClassAttractor.clear();
289   FacesWithEnforcedVertices.Clear();
290   FaceId2EnforcedVertexCoords.clear();
291   EnfVertexCoords2ProjVertex.clear();
292   EnfVertexCoords2EnfVertexList.clear();
293
294 #ifdef WITH_SMESH_CANCEL_COMPUTE
295   _compute_canceled = false;
296 #endif
297 }
298
299 //=============================================================================
300 /*!
301  *
302  */
303 //=============================================================================
304
305 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
306 {
307   MESSAGE("BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF");
308 }
309
310
311 //=============================================================================
312 /*!
313  *
314  */
315 //=============================================================================
316
317 bool BLSURFPlugin_BLSURF::CheckHypothesis
318                          (SMESH_Mesh&                          aMesh,
319                           const TopoDS_Shape&                  aShape,
320                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
321 {
322   _hypothesis = NULL;
323
324   list<const SMESHDS_Hypothesis*>::const_iterator itl;
325   const SMESHDS_Hypothesis* theHyp;
326
327   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
328   int nbHyp = hyps.size();
329   if (!nbHyp)
330   {
331     aStatus = SMESH_Hypothesis::HYP_OK;
332     return true;  // can work with no hypothesis
333   }
334
335   itl = hyps.begin();
336   theHyp = (*itl); // use only the first hypothesis
337
338   string hypName = theHyp->GetName();
339
340   if (hypName == "BLSURF_Parameters")
341   {
342     _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
343     ASSERT(_hypothesis);
344     if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
345          _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
346       //  hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
347       aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
348     else
349       aStatus = SMESH_Hypothesis::HYP_OK;
350   }
351   else
352     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
353
354   return aStatus == SMESH_Hypothesis::HYP_OK;
355 }
356
357 //=============================================================================
358 /*!
359  * Pass parameters to BLSURF
360  */
361 //=============================================================================
362
363 inline std::string to_string(double d)
364 {
365    std::ostringstream o;
366    o << d;
367    return o.str();
368 }
369
370 inline std::string to_string(int i)
371 {
372    std::ostringstream o;
373    o << i;
374    return o.str();
375 }
376
377 double _smp_phy_size;
378 // #if BLSURF_VERSION_LONG >= "3.1.1"
379 // //   sizemap_t *geo_sizemap_e, *geo_sizemap_f;
380 //   sizemap_t *iso_sizemap_p, *iso_sizemap_e, *iso_sizemap_f;
381 // //   sizemap_t *clean_geo_sizemap_e, *clean_geo_sizemap_f;
382 //   sizemap_t *clean_iso_sizemap_p, *clean_iso_sizemap_e, *clean_iso_sizemap_f;
383 // #endif
384 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
385 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
386 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
387
388 double my_u_min=1e6,my_v_min=1e6,my_u_max=-1e6,my_v_max=-1e6;
389
390 typedef struct {
391         gp_XY uv;
392         gp_XYZ xyz;
393 } projectionPoint;
394 /////////////////////////////////////////////////////////
395 projectionPoint getProjectionPoint(const TopoDS_Face& face, const gp_Pnt& point)
396 {
397   projectionPoint myPoint;
398   Handle(Geom_Surface) surface = BRep_Tool::Surface(face);
399   GeomAPI_ProjectPointOnSurf projector( point, surface );
400   if ( !projector.IsDone() || projector.NbPoints()==0 )
401     throw "getProjectionPoint: Can't project";
402
403   Quantity_Parameter u,v;
404   projector.LowerDistanceParameters(u,v);
405   myPoint.uv = gp_XY(u,v);
406   gp_Pnt aPnt = projector.NearestPoint();
407   myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
408   //return gp_XY(u,v);
409   return myPoint;
410 }
411 /////////////////////////////////////////////////////////
412
413 /////////////////////////////////////////////////////////
414 double getT(const TopoDS_Edge& edge, const gp_Pnt& point)
415 {
416   Standard_Real f,l;
417   Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, f,l);
418   GeomAPI_ProjectPointOnCurve projector( point, curve);
419   if ( projector.NbPoints() == 0 )
420     throw;
421   return projector.LowerDistanceParameter();
422 }
423
424 /////////////////////////////////////////////////////////
425 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
426 {
427   MESSAGE("BLSURFPlugin_BLSURF::entryToShape "<<entry );
428   GEOM::GEOM_Object_var aGeomObj;
429   TopoDS_Shape S = TopoDS_Shape();
430   SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
431   SALOMEDS::GenericAttribute_var anAttr;
432
433   if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR")) {
434     SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
435     CORBA::String_var aVal = anIOR->Value();
436     CORBA::Object_var obj = myStudy->ConvertIORToObject(aVal);
437     aGeomObj = GEOM::GEOM_Object::_narrow(obj);
438   }
439   if ( !aGeomObj->_is_nil() )
440     S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
441   return S;
442 }
443
444 void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
445 {
446   BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
447   enf_coords.clear();
448   coords.clear();
449   s_coords.clear();
450
451   // Get the (u,v) values of the enforced vertex on the face
452   projectionPoint myPoint = getProjectionPoint(faceShape,aPnt);
453
454   MESSAGE("Enforced Vertex: " << aPnt.X() << ", " << aPnt.Y() << ", " << aPnt.Z());
455   MESSAGE("Projected Vertex: " << myPoint.xyz.X() << ", " << myPoint.xyz.Y() << ", " << myPoint.xyz.Z());
456   MESSAGE("Parametric coordinates: " << myPoint.uv.X() << ", " << myPoint.uv.Y() );
457
458   enf_coords.push_back(aPnt.X());
459   enf_coords.push_back(aPnt.Y());
460   enf_coords.push_back(aPnt.Z());
461
462   coords.push_back(myPoint.uv.X());
463   coords.push_back(myPoint.uv.Y());
464   coords.push_back(myPoint.xyz.X());
465   coords.push_back(myPoint.xyz.Y());
466   coords.push_back(myPoint.xyz.Z());
467
468   s_coords.push_back(myPoint.xyz.X());
469   s_coords.push_back(myPoint.xyz.Y());
470   s_coords.push_back(myPoint.xyz.Z());
471
472   // Save pair projected vertex / enf vertex
473   MESSAGE("Storing pair projected vertex / enf vertex:");
474   MESSAGE("("<< myPoint.xyz.X() << ", " << myPoint.xyz.Y() << ", " << myPoint.xyz.Z() <<") / (" << aPnt.X() << ", " << aPnt.Y() << ", " << aPnt.Z()<<")");
475   EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
476   MESSAGE("Group name is: \"" << enfVertex->grpName << "\"");
477   pair<BLSURFPlugin_Hypothesis::TEnfVertexList::iterator,bool> ret;
478   BLSURFPlugin_Hypothesis::TEnfVertexList::iterator it;
479   ret = EnfVertexCoords2EnfVertexList[s_coords].insert(enfVertex);
480   if (ret.second == false) {
481     it = ret.first;
482     (*it)->grpName = enfVertex->grpName;
483   }
484
485   int key = 0;
486   if (! FacesWithEnforcedVertices.Contains(faceShape)) {
487     key = FacesWithEnforcedVertices.Add(faceShape);
488   }
489   else {
490     key = FacesWithEnforcedVertices.FindIndex(faceShape);
491   }
492
493   // If a node is already created by an attractor, do not create enforced vertex
494   int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
495   bool sameAttractor = false;
496   if (attractorKey >= 0)
497     if (FaceId2AttractorCoords.count(attractorKey) > 0)
498       if (FaceId2AttractorCoords[attractorKey] == coords)
499         sameAttractor = true;
500
501   if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
502     MESSAGE("Map of enf. vertex has key " << key)
503     MESSAGE("Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
504     if (! sameAttractor)
505       FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
506     else
507       MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
508     MESSAGE("New Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
509   }
510   else {
511     MESSAGE("Map of enf. vertex has not key " << key << ": creating it")
512     if (! sameAttractor) {
513       BLSURFPlugin_Hypothesis::TEnfVertexCoordsList ens;
514       ens.insert(coords);
515       FaceId2EnforcedVertexCoords[key] = ens;
516     }
517     else
518       MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
519   }
520 }
521
522 /////////////////////////////////////////////////////////
523 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
524 {
525   BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
526   gp_Pnt aPnt;
527
528   BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
529
530   for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
531     enfVertex = *enfVertexListIt;
532     // Case of manual coords
533     if (enfVertex->coords.size() != 0) {
534       aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
535       _createEnforcedVertexOnFace( TopoDS::Face(faceShape),  aPnt, enfVertex);
536     }
537
538     // Case of geom vertex coords
539     if (enfVertex->geomEntry != "") {
540       TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
541       TopAbs_ShapeEnum GeomType  = GeomShape.ShapeType();
542        if (GeomType == TopAbs_VERTEX){
543          aPnt = BRep_Tool::Pnt(TopoDS::Vertex(GeomShape));
544          _createEnforcedVertexOnFace( TopoDS::Face(faceShape),  aPnt, enfVertex);
545        }
546        // Group Management
547        if (GeomType == TopAbs_COMPOUND){
548          for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
549            if (it.Value().ShapeType() == TopAbs_VERTEX){
550              aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
551              _createEnforcedVertexOnFace( TopoDS::Face(faceShape),  aPnt, enfVertex);
552            }
553          }
554        }
555     }
556   }
557 }
558
559 /////////////////////////////////////////////////////////
560 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction)
561 {
562   MESSAGE("Attractor function: "<< AttractorFunction);
563   double xa, ya, za; // Coordinates of attractor point
564   double a, b;       // Attractor parameter
565   double d = 0.;
566   bool createNode=false; // To create a node on attractor projection
567   int pos1, pos2;
568   const char *sep = ";";
569   // atIt->second has the following pattern:
570   // ATTRACTOR(xa;ya;za;a;b;True|False;d)
571   // where:
572   // xa;ya;za : coordinates of  attractor
573   // a        : desired size on attractor
574   // b        : distance of influence of attractor
575   // d        : distance until which the size remains constant
576   //
577   // We search the parameters in the string
578   // xa
579   pos1 = AttractorFunction.find(sep);
580   if (pos1!=string::npos)
581   xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
582   // ya
583   pos2 = AttractorFunction.find(sep, pos1+1);
584   if (pos2!=string::npos) {
585   ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
586   pos1 = pos2;
587   }
588   // za
589   pos2 = AttractorFunction.find(sep, pos1+1);
590   if (pos2!=string::npos) {
591   za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
592   pos1 = pos2;
593   }
594   // a
595   pos2 = AttractorFunction.find(sep, pos1+1);
596   if (pos2!=string::npos) {
597   a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
598   pos1 = pos2;
599   }
600   // b
601   pos2 = AttractorFunction.find(sep, pos1+1);
602   if (pos2!=string::npos) {
603   b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
604   pos1 = pos2;
605   }
606   // createNode
607   pos2 = AttractorFunction.find(sep, pos1+1);
608   if (pos2!=string::npos) {
609     string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
610     MESSAGE("createNode: " << createNodeStr);
611     createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
612     pos1=pos2;
613   }
614   // d
615   pos2 = AttractorFunction.find(")");
616   if (pos2!=string::npos) {
617   d = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
618   }
619
620   // Get the (u,v) values of the attractor on the face
621   projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
622   gp_XY uvPoint = myPoint.uv;
623   gp_XYZ xyzPoint = myPoint.xyz;
624   Standard_Real u0 = uvPoint.X();
625   Standard_Real v0 = uvPoint.Y();
626   Standard_Real x0 = xyzPoint.X();
627   Standard_Real y0 = xyzPoint.Y();
628   Standard_Real z0 = xyzPoint.Z();
629   std::vector<double> coords;
630   coords.push_back(u0);
631   coords.push_back(v0);
632   coords.push_back(x0);
633   coords.push_back(y0);
634   coords.push_back(z0);
635   // We construct the python function
636   ostringstream attractorFunctionStream;
637   attractorFunctionStream << "def f(u,v): return ";
638   attractorFunctionStream << _smp_phy_size << "-(" << _smp_phy_size <<"-" << a << ")";
639   //attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
640   // rnc: make possible to keep the size constant until
641   // a defined distance. Distance is expressed as the positiv part
642   // of r-d where r is the distance to (u0,v0)
643   attractorFunctionStream << "*exp(-(0.5*(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"+abs(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"))/(" << b << "))**2)";
644
645   MESSAGE("Python function for attractor:" << std::endl << attractorFunctionStream.str());
646
647   int key;
648   if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
649     key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
650   }
651   else {
652     key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
653   }
654   FaceId2SizeMap[key] =attractorFunctionStream.str();
655   if (createNode) {
656     MESSAGE("Creating node on ("<<x0<<","<<y0<<","<<z0<<")");
657     FaceId2AttractorCoords[key] = coords;
658   }
659 //   // Test for new attractors
660 //   gp_Pnt myP(xyzPoint);
661 //   TopoDS_Vertex myV = BRepBuilderAPI_MakeVertex(myP);
662 //   BLSURFPlugin_Attractor myAttractor(TopoDS::Face(GeomShape),myV,200);
663 //   myAttractor.SetParameters(a, _smp_phy_size, b, d);
664 //   myAttractor.SetType(1);
665 //   FaceId2ClassAttractor[key] = myAttractor;
666 //   if(FaceId2ClassAttractor[key].GetFace().IsNull()){
667 //     MESSAGE("face nulle ");
668 //   }
669 //   else
670 //     MESSAGE("face OK");
671 //
672 //   if (FaceId2ClassAttractor[key].GetAttractorShape().IsNull()){
673 //     MESSAGE("pas de point");
674 //   }
675 //   else
676 //     MESSAGE("point OK");
677 }
678
679 /////////////////////////////////////////////////////////
680
681 void BLSURFPlugin_BLSURF::SetParameters(
682 // #if BLSURF_VERSION_LONG >= "3.1.1"
683 //                                         cad_t *                          c,
684 // #endif
685                                         const BLSURFPlugin_Hypothesis* hyp,
686                                         blsurf_session_t *             bls,
687                                         precad_session_t *             pcs,
688                                         SMESH_Mesh&                   mesh,
689                                         bool *                  use_precad
690                                        )
691 {
692   int    _topology      = BLSURFPlugin_Hypothesis::GetDefaultTopology();
693   int    _physicalMesh  = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
694   double _phySize       = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
695   int    _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
696   double _angleMeshS    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
697   double _angleMeshC    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
698   double _gradation     = BLSURFPlugin_Hypothesis::GetDefaultGradation();
699   bool   _quadAllowed   = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
700   bool   _decimesh      = BLSURFPlugin_Hypothesis::GetDefaultDecimesh();
701   int    _verb          = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
702
703   // PreCAD
704   int _precadMergeEdges      = BLSURFPlugin_Hypothesis::GetDefaultPreCADMergeEdges();
705   int _precadRemoveNanoEdges = BLSURFPlugin_Hypothesis::GetDefaultPreCADRemoveNanoEdges();
706   int _precadDiscardInput    = BLSURFPlugin_Hypothesis::GetDefaultPreCADDiscardInput();
707   double _precadEpsNano      = BLSURFPlugin_Hypothesis::GetDefaultPreCADEpsNano();
708
709
710   if (hyp) {
711     MESSAGE("BLSURFPlugin_BLSURF::SetParameters");
712     _topology      = (int) hyp->GetTopology();
713     _physicalMesh  = (int) hyp->GetPhysicalMesh();
714     _phySize       = hyp->GetPhySize();
715     _geometricMesh = (int) hyp->GetGeometricMesh();
716     _angleMeshS    = hyp->GetAngleMeshS();
717     _angleMeshC    = hyp->GetAngleMeshC();
718     _gradation     = hyp->GetGradation();
719     _quadAllowed   = hyp->GetQuadAllowed();
720     _decimesh      = hyp->GetDecimesh();
721     _verb          = hyp->GetVerbosity();
722     if ( hyp->GetPhyMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
723       blsurf_set_param(bls, "hphymin", to_string(hyp->GetPhyMin()).c_str());
724     if ( hyp->GetPhyMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
725       blsurf_set_param(bls, "hphymax", to_string(hyp->GetPhyMax()).c_str());
726     if ( hyp->GetGeoMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
727       blsurf_set_param(bls, "hgeomin", to_string(hyp->GetGeoMin()).c_str());
728     if ( hyp->GetGeoMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
729       blsurf_set_param(bls, "hgeomax", to_string(hyp->GetGeoMax()).c_str());
730
731     const BLSURFPlugin_Hypothesis::TOptionValues & opts = hyp->GetOptionValues();
732     BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
733     for ( opIt = opts.begin(); opIt != opts.end(); ++opIt )
734       if ( !opIt->second.empty() ) {
735         MESSAGE("blsurf_set_param(): " << opIt->first << " = " << opIt->second);
736         blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
737       }
738       
739     const BLSURFPlugin_Hypothesis::TOptionValues & preCADopts = hyp->GetPreCADOptionValues();
740     for ( opIt = preCADopts.begin(); opIt != preCADopts.end(); ++opIt )
741       if ( !opIt->second.empty() ) {
742         if (_topology == BLSURFPlugin_Hypothesis::PreCAD) {
743           MESSAGE("precad_set_param(): " << opIt->first << " = " << opIt->second);
744           blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
745         }
746       }
747
748     // PreCAD
749     _precadMergeEdges = hyp->GetPreCADMergeEdges();
750     _precadRemoveNanoEdges = hyp->GetPreCADRemoveNanoEdges();
751     _precadDiscardInput = hyp->GetPreCADDiscardInput();
752     _precadEpsNano = hyp->GetPreCADEpsNano();
753
754   } else {
755     //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
756     // GetDefaultPhySize() sometimes leads to computation failure
757     _phySize = mesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
758     MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
759   }
760
761   // PreCAD
762   if (_topology == BLSURFPlugin_Hypothesis::PreCAD) {
763     *use_precad = true;
764     precad_set_param(pcs, "verbose",                to_string(_verb).c_str());
765     precad_set_param(pcs, "merge_edges",            _precadMergeEdges ? "1" : "0");
766     precad_set_param(pcs, "remove_nano_edges",      _precadRemoveNanoEdges ? "1" : "0");
767     precad_set_param(pcs, "discard_input_topology", _precadDiscardInput ? "1" : "0");
768     if ( _precadEpsNano != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
769       precad_set_param(pcs, "eps_nano",               to_string(_precadEpsNano).c_str());
770   }
771
772   _smp_phy_size = _phySize;
773   blsurf_set_param(bls, "topo_points",       _topology == BLSURFPlugin_Hypothesis::Process || _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
774   blsurf_set_param(bls, "topo_curves",       _topology == BLSURFPlugin_Hypothesis::Process || _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
775   blsurf_set_param(bls, "topo_project",      _topology == BLSURFPlugin_Hypothesis::Process || _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
776   blsurf_set_param(bls, "clean_boundary",    _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
777   blsurf_set_param(bls, "close_boundary",    _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
778   blsurf_set_param(bls, "hphy_flag",         to_string(_physicalMesh).c_str());
779   blsurf_set_param(bls, "hphydef",           to_string(_phySize).c_str());
780   blsurf_set_param(bls, "hgeo_flag",         to_string(_geometricMesh).c_str());
781   blsurf_set_param(bls, "relax_size",        _decimesh ? "0": to_string(_geometricMesh).c_str());
782   blsurf_set_param(bls, "angle_meshs",       to_string(_angleMeshS).c_str());
783   blsurf_set_param(bls, "angle_meshc",       to_string(_angleMeshC).c_str());
784   blsurf_set_param(bls, "gradation",         to_string(_gradation).c_str());
785   blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
786   blsurf_set_param(bls, "element",           _quadAllowed ? "q1.0" : "p1");
787   blsurf_set_param(bls, "verb",              to_string(_verb).c_str());
788   
789   if (_physicalMesh == BLSURFPlugin_Hypothesis::SizeMap){
790     TopoDS_Shape GeomShape;
791     TopoDS_Shape AttShape;
792     TopAbs_ShapeEnum GeomType;
793     //
794     // Standard Size Maps
795     //
796     MESSAGE("Setting a Size Map");
797     const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
798     BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
799     for ( ; smIt != sizeMaps.end(); ++smIt ) {
800       if ( !smIt->second.empty() ) {
801         MESSAGE("blsurf_set_sizeMap(): " << smIt->first << " = " << smIt->second);
802         GeomShape = entryToShape(smIt->first);
803         GeomType  = GeomShape.ShapeType();
804         MESSAGE("Geomtype is " << GeomType);
805         int key = -1;
806         // Group Management
807         if (GeomType == TopAbs_COMPOUND){
808           for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
809             // Group of faces
810             if (it.Value().ShapeType() == TopAbs_FACE){
811               HasSizeMapOnFace = true;
812               if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
813                 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
814               }
815               else {
816                 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
817 //                 MESSAGE("Face with key " << key << " already in map");
818               }
819               FaceId2SizeMap[key] = smIt->second;
820             }
821             // Group of edges
822             if (it.Value().ShapeType() == TopAbs_EDGE){
823               HasSizeMapOnEdge = true;
824               HasSizeMapOnFace = true;
825               if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
826                 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
827               }
828               else {
829                 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
830 //                 MESSAGE("Edge with key " << key << " already in map");
831               }
832               EdgeId2SizeMap[key] = smIt->second;
833             }
834             // Group of vertices
835             if (it.Value().ShapeType() == TopAbs_VERTEX){
836               HasSizeMapOnVertex = true;
837               HasSizeMapOnEdge = true;
838               HasSizeMapOnFace = true;
839               if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
840                 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
841               }
842               else {
843                 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
844                 MESSAGE("Group of vertices with key " << key << " already in map");
845               }
846               MESSAGE("Group of vertices with key " << key << " has a size map: " << smIt->second);
847               VertexId2SizeMap[key] = smIt->second;
848             }
849           }
850         }
851         // Single face
852         if (GeomType == TopAbs_FACE){
853           HasSizeMapOnFace = true;
854           if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
855             key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
856           }
857           else {
858             key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
859 //             MESSAGE("Face with key " << key << " already in map");
860           }
861           FaceId2SizeMap[key] = smIt->second;
862         }
863         // Single edge
864         if (GeomType == TopAbs_EDGE){
865           HasSizeMapOnEdge = true;
866           HasSizeMapOnFace = true;
867           if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
868             key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
869           }
870           else {
871             key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
872 //             MESSAGE("Edge with key " << key << " already in map");
873           }
874           EdgeId2SizeMap[key] = smIt->second;
875         }
876         // Single vertex
877         if (GeomType == TopAbs_VERTEX){
878           HasSizeMapOnVertex = true;
879           HasSizeMapOnEdge   = true;
880           HasSizeMapOnFace   = true;
881           if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
882             key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
883           }
884           else {
885             key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
886              MESSAGE("Vertex with key " << key << " already in map");
887           }
888           MESSAGE("Vertex with key " << key << " has a size map: " << smIt->second);
889           VertexId2SizeMap[key] = smIt->second;
890         }
891       }
892     }
893
894     //
895     // Attractors
896     //
897     // TODO appeler le constructeur des attracteurs directement ici
898     MESSAGE("Setting Attractors");
899     const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
900     BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
901     for ( ; atIt != attractors.end(); ++atIt ) {
902       if ( !atIt->second.empty() ) {
903         MESSAGE("blsurf_set_attractor(): " << atIt->first << " = " << atIt->second);
904         GeomShape = entryToShape(atIt->first);
905         GeomType  = GeomShape.ShapeType();
906         // Group Management
907         if (GeomType == TopAbs_COMPOUND){
908           for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
909             if (it.Value().ShapeType() == TopAbs_FACE){
910               HasSizeMapOnFace = true;
911               createAttractorOnFace(it.Value(), atIt->second);
912             }
913           }
914         }
915
916         if (GeomType == TopAbs_FACE){
917           HasSizeMapOnFace = true;
918           createAttractorOnFace(GeomShape, atIt->second);
919         }
920 /*
921         if (GeomType == TopAbs_EDGE){
922           HasSizeMapOnEdge = true;
923           HasSizeMapOnFace = true;
924         EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
925         }
926         if (GeomType == TopAbs_VERTEX){
927           HasSizeMapOnVertex = true;
928           HasSizeMapOnEdge   = true;
929           HasSizeMapOnFace   = true;
930         VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
931         }
932 */
933       }
934     }
935
936     // Class Attractors
937     // temporary commented out for testing
938     // TODO
939     //  - Fill in the BLSURFPlugin_Hypothesis::TAttractorMap map in the hypothesis
940     //  - Uncomment and complete this part to construct the attractors from the attractor shape and the passed parameters on each face of the map
941     //  - To do this use the public methodss: SetParameters(several double parameters) and SetType(int type)
942     //  OR, even better:
943     //  - Construct the attractors with an empty dist. map in the hypothesis
944     //  - build the map here for each face with an attractor set and only if the attractor shape as changed since the last call to _buildmap()
945     //  -> define a bool _mapbuilt in the class that is set to false by default and set to true when calling _buildmap()  OK
946
947     const BLSURFPlugin_Hypothesis::TAttractorMap class_attractors = BLSURFPlugin_Hypothesis::GetClassAttractorEntries(hyp);
948     int key=-1;
949     BLSURFPlugin_Hypothesis::TAttractorMap::const_iterator AtIt = class_attractors.begin();
950     for ( ; AtIt != class_attractors.end(); ++AtIt ) {
951       if ( !AtIt->second->Empty() ) {
952        // MESSAGE("blsurf_set_attractor(): " << AtIt->first << " = " << AtIt->second);
953         GeomShape = entryToShape(AtIt->first);
954         AttShape = AtIt->second->GetAttractorShape();
955         GeomType  = GeomShape.ShapeType();
956         // Group Management
957 //         if (GeomType == TopAbs_COMPOUND){
958 //           for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
959 //             if (it.Value().ShapeType() == TopAbs_FACE){
960 //               HasAttractorOnFace = true;
961 //               myAttractor = BLSURFPluginAttractor(GeomShape, AttShape);
962 //             }
963 //           }
964 //         }
965
966         if (GeomType == TopAbs_FACE
967           && (AttShape.ShapeType() == TopAbs_VERTEX
968            || AttShape.ShapeType() == TopAbs_EDGE
969            || AttShape.ShapeType() == TopAbs_WIRE
970            || AttShape.ShapeType() == TopAbs_COMPOUND) ){
971             HasSizeMapOnFace = true;
972
973             if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape)) ) {
974                 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
975             }
976             else {
977               key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
978 //                 MESSAGE("Face with key " << key << " already in map");
979             }
980
981             FaceId2ClassAttractor[key] = AtIt->second;
982         }
983         else{
984           MESSAGE("Wrong shape type !!")
985         }
986
987       }
988     }
989
990
991     //
992     // Enforced Vertices
993     //
994     MESSAGE("Setting Enforced Vertices");
995     const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
996     BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
997     for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
998       if ( !enfIt->second.empty() ) {
999         GeomShape = entryToShape(enfIt->first);
1000         GeomType  = GeomShape.ShapeType();
1001         // Group Management
1002         if (GeomType == TopAbs_COMPOUND){
1003           for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1004             if (it.Value().ShapeType() == TopAbs_FACE){
1005               HasSizeMapOnFace = true;
1006               createEnforcedVertexOnFace(it.Value(), enfIt->second);
1007             }
1008           }
1009         }
1010
1011         if (GeomType == TopAbs_FACE){
1012           HasSizeMapOnFace = true;
1013           createEnforcedVertexOnFace(GeomShape, enfIt->second);
1014         }
1015       }
1016     }
1017
1018     // Internal vertices
1019     bool useInternalVertexAllFaces = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFaces(hyp);
1020     if (useInternalVertexAllFaces) {
1021       std::string grpName = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFacesGroup(hyp);
1022       MESSAGE("Setting Internal Enforced Vertices");
1023       GeomShape = mesh.GetShapeToMesh();
1024       gp_Pnt aPnt;
1025       TopExp_Explorer exp (GeomShape, TopAbs_FACE);
1026       for (; exp.More(); exp.Next()){
1027         MESSAGE("Iterating shapes. Shape type is " << exp.Current().ShapeType());
1028         TopExp_Explorer exp_face (exp.Current(), TopAbs_VERTEX);
1029         for (; exp_face.More(); exp_face.Next())
1030         {
1031           // Get coords of vertex
1032           // Check if current coords is already in enfVertexList
1033           // If coords not in enfVertexList, add new enfVertex
1034           aPnt = BRep_Tool::Pnt(TopoDS::Vertex(exp_face.Current()));
1035           MESSAGE("Found vertex on face at " << aPnt.X() <<", "<<aPnt.Y()<<", "<<aPnt.Z());
1036           BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex = new BLSURFPlugin_Hypothesis::TEnfVertex();
1037           enfVertex->coords.push_back(aPnt.X());
1038           enfVertex->coords.push_back(aPnt.Y());
1039           enfVertex->coords.push_back(aPnt.Z());
1040           enfVertex->name = "";
1041           enfVertex->faceEntries.clear();
1042           enfVertex->geomEntry = "";
1043           enfVertex->grpName = grpName;
1044           _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()),  aPnt, enfVertex);
1045           HasSizeMapOnFace = true;
1046         }
1047       }
1048     }
1049
1050     MESSAGE("Setting Size Map on FACES ");
1051 // #if BLSURF_VERSION_LONG < "3.1.1"
1052     blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
1053 // #else
1054 //     if (*use_precad)
1055 //       iso_sizemap_f = sizemap_new(c, distene_sizemap_type_iso_cad_face, (void *)size_on_surface, NULL);
1056 //     else
1057 //       clean_iso_sizemap_f = sizemap_new(c, distene_sizemap_type_iso_cad_face, (void *)size_on_surface, NULL);
1058 // #endif
1059
1060     if (HasSizeMapOnEdge){
1061       MESSAGE("Setting Size Map on EDGES ");
1062 // #if BLSURF_VERSION_LONG < "3.1.1"
1063       blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
1064 // #else
1065 //       if (*use_precad)
1066 //         iso_sizemap_e = sizemap_new(c, distene_sizemap_type_iso_cad_edge, (void *)size_on_edge, NULL);
1067 //       else
1068 //         clean_iso_sizemap_e = sizemap_new(c, distene_sizemap_type_iso_cad_edge, (void *)size_on_edge, NULL);
1069 // #endif
1070     }
1071     if (HasSizeMapOnVertex){
1072       MESSAGE("Setting Size Map on VERTICES ");
1073 // #if BLSURF_VERSION_LONG < "3.1.1"
1074       blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
1075 // #else
1076 //       if (*use_precad)
1077 //         iso_sizemap_p = sizemap_new(c, distene_sizemap_type_iso_cad_point, (void *)size_on_vertex, NULL);
1078 //       else
1079 //         clean_iso_sizemap_p = sizemap_new(c, distene_sizemap_type_iso_cad_point, (void *)size_on_vertex, NULL);
1080 // #endif
1081     }
1082   }
1083 }
1084
1085 namespace
1086 {
1087   /*!
1088    * \brief Class correctly terminating usage of BLSURF library at destruction
1089    */
1090   class BLSURF_Cleaner
1091   {
1092     context_t *       _ctx;
1093     blsurf_session_t* _bls;
1094     cad_t *           _cad;
1095     dcad_t *          _dcad;
1096   public:
1097     BLSURF_Cleaner(context_t *       ctx,
1098                    blsurf_session_t* bls,
1099                    cad_t *           cad,
1100                    dcad_t *          dcad)
1101       : _ctx ( ctx  ),
1102         _bls ( bls  ),
1103         _cad ( cad  ),
1104         _dcad( dcad )
1105     {
1106     }
1107     ~BLSURF_Cleaner()
1108     {
1109       blsurf_session_delete(_bls);
1110
1111       // #if BLSURF_VERSION_LONG >= "3.1.1"
1112       // //     if(geo_sizemap_e)
1113       // //       distene_sizemap_delete(geo_sizemap_e);
1114       // //     if(geo_sizemap_f)
1115       // //       distene_sizemap_delete(geo_sizemap_f);
1116       //     if(iso_sizemap_p)
1117       //       distene_sizemap_delete(iso_sizemap_p);
1118       //     if(iso_sizemap_e)
1119       //       distene_sizemap_delete(iso_sizemap_e);
1120       //     if(iso_sizemap_f)
1121       //       distene_sizemap_delete(iso_sizemap_f);
1122       // 
1123       // //     if(clean_geo_sizemap_e)
1124       // //       distene_sizemap_delete(clean_geo_sizemap_e);
1125       // //     if(clean_geo_sizemap_f)
1126       // //       distene_sizemap_delete(clean_geo_sizemap_f);
1127       //     if(clean_iso_sizemap_p)
1128       //       distene_sizemap_delete(clean_iso_sizemap_p);
1129       //     if(clean_iso_sizemap_e)
1130       //       distene_sizemap_delete(clean_iso_sizemap_e);
1131       //     if(clean_iso_sizemap_f)
1132       //       distene_sizemap_delete(clean_iso_sizemap_f);
1133       // #endif
1134       
1135       cad_delete(_cad);
1136       dcad_delete(_dcad);
1137       context_delete(_ctx);
1138     }
1139   };
1140 } // namespace
1141
1142 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
1143 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1144                   real *duu, real *duv, real *dvv, void *user_data);
1145 status_t message_cb(message_t *msg, void *user_data);
1146 status_t interrupt_cb(integer *interrupt_status, void *user_data);
1147
1148 //=============================================================================
1149 /*!
1150  *
1151  */
1152 //=============================================================================
1153
1154 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
1155
1156   MESSAGE("BLSURFPlugin_BLSURF::Compute");
1157
1158   // Fix problem with locales
1159   Kernel_Utils::Localizer aLocalizer;
1160
1161   /* create a distene context (generic object) */
1162   status_t status = STATUS_ERROR;
1163
1164   context_t *ctx =  context_new();
1165
1166   /* Set the message callback in the working context */
1167   context_set_message_callback(ctx, message_cb, &_comment);
1168 #ifdef WITH_SMESH_CANCEL_COMPUTE
1169   _compute_canceled = false;
1170   context_set_interrupt_callback(ctx, interrupt_cb, this);
1171 #else
1172   context_set_interrupt_callback(ctx, interrupt_cb, NULL);
1173 #endif
1174
1175   /* create the CAD object we will work on. It is associated to the context ctx. */
1176   cad_t *c = cad_new(ctx);
1177   dcad_t *dcad = dcad_new(c);
1178
1179   FacesWithSizeMap.Clear();
1180   FaceId2SizeMap.clear();
1181   FaceId2ClassAttractor.clear();
1182   FaceIndex2ClassAttractor.clear();
1183   EdgesWithSizeMap.Clear();
1184   EdgeId2SizeMap.clear();
1185   VerticesWithSizeMap.Clear();
1186   VertexId2SizeMap.clear();
1187
1188   SMESH_MesherHelper helper( aMesh );
1189   // do not call helper.IsQuadraticSubMesh() because submeshes
1190   // may be cleaned and helper.myTLinkNodeMap gets invalid in such a case
1191   const bool haveQudraticSubMesh = SMESH_MesherHelper( aMesh ).IsQuadraticSubMesh( aShape );
1192   helper.SetIsQuadratic( haveQudraticSubMesh );
1193   bool needMerge = false;
1194   set< SMESH_subMesh* > edgeSubmeshes;
1195
1196   /* Now fill the CAD object with data from your CAD
1197    * environement. This is the most complex part of a successfull
1198    * integration.
1199    */
1200
1201   // PreCAD
1202   // If user requests it, send the CAD through Distene preprocessor : PreCAD
1203   cad_t *cleanc = NULL;
1204   precad_session_t *pcs = precad_session_new(ctx);
1205   precad_data_set_cad(pcs, c);
1206
1207   blsurf_session_t *bls = blsurf_session_new(ctx);
1208
1209   // an object that correctly deletes all blsurf objects at destruction
1210   BLSURF_Cleaner cleaner( ctx,bls,c,dcad );
1211
1212   MESSAGE("BEGIN SetParameters");
1213   bool use_precad = false;
1214   SetParameters(
1215 // #if BLSURF_VERSION_LONG >= "3.1.1"
1216 //     c,
1217 // #endif
1218     _hypothesis, bls, pcs, aMesh, &use_precad);
1219   MESSAGE("END SetParameters");
1220
1221   // needed to prevent the opencascade memory managmement from freeing things
1222   vector<Handle(Geom2d_Curve)> curves;
1223   vector<Handle(Geom_Surface)> surfaces;
1224
1225   surfaces.resize(0);
1226   curves.resize(0);
1227
1228   TopTools_IndexedMapOfShape fmap;
1229   TopTools_IndexedMapOfShape emap;
1230   TopTools_IndexedMapOfShape pmap;
1231
1232   fmap.Clear();
1233   FaceId2PythonSmp.clear();
1234   emap.Clear();
1235   EdgeId2PythonSmp.clear();
1236   pmap.Clear();
1237   VertexId2PythonSmp.clear();
1238
1239   assert(Py_IsInitialized());
1240   PyGILState_STATE gstate;
1241   gstate = PyGILState_Ensure();
1242
1243   string theSizeMapStr;
1244
1245   /****************************************************************************************
1246                                   FACES
1247   *****************************************************************************************/
1248   int iface = 0;
1249   string bad_end = "return";
1250   int faceKey = -1;
1251   TopTools_IndexedMapOfShape _map;
1252   TopExp::MapShapes(aShape,TopAbs_VERTEX,_map);
1253   int ienf = _map.Extent();
1254
1255   for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
1256     TopoDS_Face f=TopoDS::Face(face_iter.Current());
1257
1258     // make INTERNAL face oriented FORWARD (issue 0020993)
1259     if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
1260       f.Orientation(TopAbs_FORWARD);
1261
1262     if (fmap.FindIndex(f) > 0)
1263       continue;
1264
1265     fmap.Add(f);
1266     iface++;
1267     surfaces.push_back(BRep_Tool::Surface(f));
1268
1269     /* create an object representing the face for blsurf */
1270     /* where face_id is an integer identifying the face.
1271      * surf_function is the function that defines the surface
1272      * (For this face, it will be called by blsurf with your_face_object_ptr
1273      * as last parameter.
1274      */
1275     cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
1276
1277     /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
1278     cad_face_set_tag(fce, iface);
1279
1280     /* Set face orientation (optional if you want a well oriented output mesh)*/
1281     if(f.Orientation() != TopAbs_FORWARD){
1282       cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
1283     } else {
1284       cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
1285     }
1286
1287     if (HasSizeMapOnFace && !use_precad){
1288 //       MESSAGE("A size map is defined on a face")
1289 //       std::cout << "A size map is defined on a face" << std::endl;
1290       // Classic size map
1291       faceKey = FacesWithSizeMap.FindIndex(f);
1292
1293
1294       if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
1295         MESSAGE("A size map is defined on face :"<<faceKey)
1296         theSizeMapStr = FaceId2SizeMap[faceKey];
1297         // check if function ends with "return"
1298         if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1299           continue;
1300         // Expr To Python function, verification is performed at validation in GUI
1301         PyObject * obj = NULL;
1302         obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1303         Py_DECREF(obj);
1304         PyObject * func = NULL;
1305         func = PyObject_GetAttrString(main_mod, "f");
1306         FaceId2PythonSmp[iface]=func;
1307         FaceId2SizeMap.erase(faceKey);
1308       }
1309
1310       // Specific size map = Attractor
1311       std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
1312
1313       for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
1314         if (attractor_iter->first == faceKey) {
1315           MESSAGE("Face indice: " << iface);
1316           MESSAGE("Adding attractor");
1317
1318           double xyzCoords[3]  = {attractor_iter->second[2],
1319                                   attractor_iter->second[3],
1320                                   attractor_iter->second[4]};
1321
1322           MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1323           gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1324           BRepClass_FaceClassifier scl(f,P,1e-7);
1325           // scl.Perform() is bugged. The function was rewritten
1326 //          scl.Perform();
1327           BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1328           TopAbs_State result = scl.State();
1329           MESSAGE("Position of point on face: "<<result);
1330           if ( result == TopAbs_OUT )
1331               MESSAGE("Point is out of face: node is not created");
1332           if ( result == TopAbs_UNKNOWN )
1333               MESSAGE("Point position on face is unknown: node is not created");
1334           if ( result == TopAbs_ON )
1335               MESSAGE("Point is on border of face: node is not created");
1336           if ( result == TopAbs_IN )
1337           {
1338             // Point is inside face and not on border
1339             MESSAGE("Point is in face: node is created");
1340             double uvCoords[2]   = {attractor_iter->second[0],attractor_iter->second[1]};
1341             ienf++;
1342             MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1343             cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1344             cad_point_set_tag(point_p, ienf);
1345           }
1346           FaceId2AttractorCoords.erase(faceKey);
1347         }
1348       }
1349
1350       // Class Attractors
1351       std::map<int,BLSURFPlugin_Attractor* >::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
1352       if (clAttractor_iter != FaceId2ClassAttractor.end()){
1353           MESSAGE("Face indice: " << iface);
1354           MESSAGE("Adding attractor");
1355           FaceIndex2ClassAttractor[iface]=clAttractor_iter->second;
1356           FaceId2ClassAttractor.erase(clAttractor_iter);
1357         }
1358       }
1359
1360       // Enforced Vertices
1361       faceKey = FacesWithEnforcedVertices.FindIndex(f);
1362       std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1363       if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1364         MESSAGE("Some enforced vertices are defined");
1365         BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl;
1366         MESSAGE("Face indice: " << iface);
1367         MESSAGE("Adding enforced vertices");
1368         evl = evmIt->second;
1369         MESSAGE("Number of vertices to add: "<< evl.size());
1370         BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
1371         for (; evlIt != evl.end(); ++evlIt) {
1372           BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
1373           xyzCoords.push_back(evlIt->at(2));
1374           xyzCoords.push_back(evlIt->at(3));
1375           xyzCoords.push_back(evlIt->at(4));
1376           MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1377           gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1378           BRepClass_FaceClassifier scl(f,P,1e-7);
1379           // OCC 6.3sp6 : scl.Perform() is bugged. The function was rewritten
1380 //          BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1381           // OCC 6.5.2: scl.Perform() is not bugged anymore
1382           scl.Perform(f, P, 1e-7);
1383           TopAbs_State result = scl.State();
1384           MESSAGE("Position of point on face: "<<result);
1385           if ( result == TopAbs_OUT ) {
1386             MESSAGE("Point is out of face: node is not created");
1387             if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1388               EnfVertexCoords2ProjVertex.erase(xyzCoords);
1389               EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1390             }
1391           }
1392           if ( result == TopAbs_UNKNOWN ) {
1393             MESSAGE("Point position on face is unknown: node is not created");
1394             if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1395               EnfVertexCoords2ProjVertex.erase(xyzCoords);
1396               EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1397             }
1398           }
1399           if ( result == TopAbs_ON ) {
1400             MESSAGE("Point is on border of face: node is not created");
1401             if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1402               EnfVertexCoords2ProjVertex.erase(xyzCoords);
1403               EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1404             }
1405           }
1406           if ( result == TopAbs_IN )
1407           {
1408             // Point is inside face and not on border
1409             MESSAGE("Point is in face: node is created");
1410             double uvCoords[2]   = {evlIt->at(0),evlIt->at(1)};
1411             ienf++;
1412             MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1413             cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1414             cad_point_set_tag(point_p, ienf);
1415           }
1416         }
1417         FaceId2EnforcedVertexCoords.erase(faceKey);
1418       }
1419 //       else
1420 //         std::cout << "No enforced vertex defined" << std::endl;
1421 //     }
1422
1423
1424     /****************************************************************************************
1425                                     EDGES
1426                    now create the edges associated to this face
1427     *****************************************************************************************/
1428     int edgeKey = -1;
1429     for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1430       TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1431       int ic = emap.FindIndex(e);
1432       if (ic <= 0)
1433         ic = emap.Add(e);
1434
1435       double tmin,tmax;
1436       curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1437
1438       if (HasSizeMapOnEdge){
1439         edgeKey = EdgesWithSizeMap.FindIndex(e);
1440         if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1441           MESSAGE("Sizemap defined on edge with index " << edgeKey);
1442           theSizeMapStr = EdgeId2SizeMap[edgeKey];
1443           if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1444             continue;
1445           // Expr To Python function, verification is performed at validation in GUI
1446           PyObject * obj = NULL;
1447           obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1448           Py_DECREF(obj);
1449           PyObject * func = NULL;
1450           func = PyObject_GetAttrString(main_mod, "f");
1451           EdgeId2PythonSmp[ic]=func;
1452           EdgeId2SizeMap.erase(edgeKey);
1453         }
1454       }
1455
1456       /* attach the edge to the current blsurf face */
1457       cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1458
1459       /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1460       cad_edge_set_tag(edg, ic);
1461
1462       /* by default, an edge does not necessalry appear in the resulting mesh,
1463      unless the following property is set :
1464       */
1465       cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1466
1467       /* by default an edge is a boundary edge */
1468       if (e.Orientation() == TopAbs_INTERNAL)
1469         cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1470
1471       // pass existing nodes of sub-meshes to BLSURF
1472       SMESH_subMesh* sm = aMesh.GetSubMesh( e );
1473       if ( !sm->IsEmpty() )
1474       {
1475         edgeSubmeshes.insert( sm );
1476
1477         StdMeshers_FaceSide edgeOfFace( f, e, &aMesh, e.Orientation() == TopAbs_FORWARD,
1478                                         /*ignoreMedium=*/haveQudraticSubMesh);
1479         if ( edgeOfFace.MissVertexNode() )
1480           return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
1481
1482         const int nbNodes = edgeOfFace.NbPoints();
1483
1484         dcad_edge_discretization_t *dedge;
1485         dcad_get_edge_discretization(dcad, edg, &dedge);
1486         dcad_edge_discretization_set_vertex_count( dedge, nbNodes );
1487
1488         const std::vector<UVPtStruct>& nodeData = edgeOfFace.GetUVPtStruct();
1489         for ( int iN = 0; iN < nbNodes; ++iN )
1490         {
1491           const UVPtStruct& nData = nodeData[ iN ];
1492           double t                = nData.param;
1493           real uv[2]              = { nData.u, nData.v };
1494           SMESH_TNodeXYZ nXYZ( nData.node );
1495           dcad_edge_discretization_set_vertex_coordinates( dedge, iN+1, t, uv, nXYZ._xyz );
1496         }
1497         dcad_edge_discretization_set_property(dedge, DISTENE_DCAD_PROPERTY_REQUIRED);
1498       }
1499
1500       /****************************************************************************************
1501                                       VERTICES
1502       *****************************************************************************************/
1503
1504       int npts = 0;
1505       int ip1, ip2, *ip;
1506       gp_Pnt2d e0 = curves.back()->Value(tmin);
1507       gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1508       Standard_Real d1=0,d2=0;
1509
1510       int vertexKey = -1;
1511       for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1512         TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1513         ++npts;
1514         if (npts == 1){
1515           ip = &ip1;
1516           d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1517         } else {
1518           ip = &ip2;
1519           d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1520         }
1521         *ip = pmap.FindIndex(v);
1522         if(*ip <= 0)
1523           *ip = pmap.Add(v);
1524
1525         if (HasSizeMapOnVertex){
1526           vertexKey = VerticesWithSizeMap.FindIndex(v);
1527           if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1528             theSizeMapStr = VertexId2SizeMap[vertexKey];
1529             //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1530             if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1531               continue;
1532             // Expr To Python function, verification is performed at validation in GUI
1533             PyObject * obj = NULL;
1534             obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1535             Py_DECREF(obj);
1536             PyObject * func = NULL;
1537             func = PyObject_GetAttrString(main_mod, "f");
1538             VertexId2PythonSmp[*ip]=func;
1539             VertexId2SizeMap.erase(vertexKey);   // do not erase if using a vector
1540           }
1541         }
1542       }
1543       if (npts != 2) {
1544         // should not happen
1545         MESSAGE("An edge does not have 2 extremities.");
1546       } else {
1547         if (d1 < d2) {
1548           // This defines the curves extremity connectivity
1549           cad_edge_set_extremities(edg, ip1, ip2);
1550           /* set the tag (color) to the same value as the extremity id : */
1551           cad_edge_set_extremities_tag(edg, ip1, ip2);
1552         }
1553         else {
1554           cad_edge_set_extremities(edg, ip2, ip1);
1555           cad_edge_set_extremities_tag(edg, ip2, ip1);
1556         }
1557       }
1558     } // for edge
1559   } //for face
1560
1561   // Clear mesh from already meshed edges if possible else
1562   // remember that merge is needed
1563   set< SMESH_subMesh* >::iterator smIt = edgeSubmeshes.begin();
1564   for ( ; smIt != edgeSubmeshes.end(); ++smIt )
1565   {
1566     SMESH_subMesh* sm = *smIt;
1567     SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
1568                                                                  /*complexFirst=*/false);
1569     while ( subsmIt->more() )
1570     {
1571       sm = subsmIt->next();
1572       if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() )
1573       {
1574         SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1575         if ( nIt->more() )
1576         {
1577           const SMDS_MeshNode* n = nIt->next();
1578           if ( n->NbInverseElements( SMDSAbs_Face ) > 0 )
1579           {
1580             needMerge = true;
1581             // add existing medium nodes to helper
1582             if ( aMesh.NbEdges( ORDER_QUADRATIC ) > 0 )
1583             {
1584               SMDS_ElemIteratorPtr edgeIt = smDS->GetElements();
1585               while ( edgeIt->more() )
1586                 helper.AddTLinks( static_cast<const SMDS_MeshEdge*>(edgeIt->next()));
1587             }
1588           }
1589           else
1590           {
1591             sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
1592           }
1593         }
1594       }
1595     }
1596   }
1597
1598   PyGILState_Release(gstate);
1599
1600   if (use_precad){
1601     /* Now launch the PreCAD process */
1602     status = precad_process(pcs);
1603     if(status != STATUS_OK){
1604       cout << "PreCAD processing failed with error code " << status << "\n";
1605       // Now we can delete the PreCAD session
1606       precad_session_delete(pcs);
1607     }
1608     else {
1609       // retrieve the pre-processed CAD object
1610       cleanc = precad_new_cad(pcs);
1611       if(!cleanc){
1612         cout << "Unable to retrieve PreCAD result \n";
1613       }
1614       cout << "PreCAD processing successfull \n";
1615
1616 // #if BLSURF_VERSION_LONG >= "3.1.1"
1617 //       /* We can now get the updated sizemaps (if any) */
1618 // //       if(geo_sizemap_e)
1619 // //         clean_geo_sizemap_e = precad_new_sizemap(pcs, geo_sizemap_e);
1620 // // 
1621 // //       if(geo_sizemap_f)
1622 // //         clean_geo_sizemap_f = precad_new_sizemap(pcs, geo_sizemap_f);
1623 // 
1624 //       if(iso_sizemap_p)
1625 //         clean_iso_sizemap_p = precad_new_sizemap(pcs, iso_sizemap_p);
1626 // 
1627 //       if(iso_sizemap_e)
1628 //         clean_iso_sizemap_e = precad_new_sizemap(pcs, iso_sizemap_e);
1629 // 
1630 //       if(iso_sizemap_f)
1631 //         clean_iso_sizemap_f = precad_new_sizemap(pcs, iso_sizemap_f);
1632 // #endif
1633       
1634       // Now we can delete the PreCAD session 
1635       precad_session_delete(pcs);
1636     }
1637   }
1638
1639   blsurf_data_set_dcad(bls, dcad);
1640   if (cleanc) {
1641     // Give the pre-processed CAD object to the current BLSurf session
1642     blsurf_data_set_cad(bls, cleanc);
1643   }
1644   else {
1645     // Use the original one
1646     blsurf_data_set_cad(bls, c);
1647   }
1648
1649 // #if BLSURF_VERSION_LONG >= "3.1.1"
1650 //   blsurf_data_set_sizemap(bls, clean_iso_sizemap_p);
1651 //   blsurf_data_set_sizemap(bls, clean_iso_sizemap_e);
1652 //   blsurf_data_set_sizemap(bls, clean_iso_sizemap_f);
1653 // #endif
1654
1655   std::cout << std::endl;
1656   std::cout << "Beginning of Surface Mesh generation" << std::endl;
1657   std::cout << std::endl;
1658
1659   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1660 #ifndef WNT
1661   feclearexcept( FE_ALL_EXCEPT );
1662   int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1663 #endif
1664
1665   try {
1666     OCC_CATCH_SIGNALS;
1667
1668     status = blsurf_compute_mesh(bls);
1669
1670   }
1671   catch ( std::exception& exc ) {
1672     _comment += exc.what();
1673   }
1674   catch (Standard_Failure& ex) {
1675     _comment += ex.DynamicType()->Name();
1676     if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1677       _comment += ": ";
1678       _comment += ex.GetMessageString();
1679     }
1680   }
1681   catch (...) {
1682     if ( _comment.empty() )
1683       _comment = "Exception in blsurf_compute_mesh()";
1684   }
1685   if ( status != STATUS_OK) {
1686     // There was an error while meshing
1687     //return error(_comment);
1688   }
1689
1690   std::cout << std::endl;
1691   std::cout << "End of Surface Mesh generation" << std::endl;
1692   std::cout << std::endl;
1693
1694   mesh_t *msh = NULL;
1695   blsurf_data_get_mesh(bls, &msh);
1696   if(!msh){
1697     /* release the mesh object */
1698     blsurf_data_regain_mesh(bls, msh);
1699     return error(_comment);
1700   }
1701
1702   std::string GMFFileName = BLSURFPlugin_Hypothesis::GetDefaultGMFFile();
1703   if (_hypothesis)
1704     GMFFileName = _hypothesis->GetGMFFile();
1705   if (GMFFileName != "") {
1706 //     bool GMFFileMode = _hypothesis->GetGMFFileMode();
1707     bool asciiFound = (GMFFileName.find(".mesh",GMFFileName.length()-5) != std::string::npos);
1708     bool binaryFound = (GMFFileName.find(".meshb",GMFFileName.length()-6) != std::string::npos);
1709     if (!asciiFound && !binaryFound)
1710       GMFFileName.append(".mesh");
1711     /*
1712     if (GMFFileMode) {
1713       if (!binaryFound) {
1714         if (asciiFound)
1715           GMFFileName.append("b");
1716         else
1717           GMFFileName.append(".meshb");
1718       }
1719     }
1720     else {
1721       if (!asciiFound)
1722         GMFFileName.append(".mesh");
1723     }
1724     */
1725     mesh_write_mesh(msh, GMFFileName.c_str());
1726   }
1727
1728   /* retrieve mesh data (see distene/mesh.h) */
1729   integer nv, ne, nt, nq, vtx[4], tag;
1730   real xyz[3];
1731
1732   mesh_get_vertex_count(msh, &nv);
1733   mesh_get_edge_count(msh, &ne);
1734   mesh_get_triangle_count(msh, &nt);
1735   mesh_get_quadrangle_count(msh, &nq);
1736
1737
1738   SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1739   SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1740   bool* tags = new bool[nv+1];
1741
1742   /* enumerated vertices */
1743   for(int iv=1;iv<=nv;iv++) {
1744     mesh_get_vertex_coordinates(msh, iv, xyz);
1745     mesh_get_vertex_tag(msh, iv, &tag);
1746     // Issue 0020656. Use vertex coordinates
1747     if ( tag > 0 && tag <= pmap.Extent() ) {
1748       TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
1749       double tol = BRep_Tool::Tolerance( v );
1750       gp_Pnt p = BRep_Tool::Pnt( v );
1751       if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
1752         xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
1753       else
1754         tag = 0; // enforced or attracted vertex
1755     }
1756     nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1757
1758     // Create group of enforced vertices if requested
1759     BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
1760     projVertex.clear();
1761     projVertex.push_back((double)xyz[0]);
1762     projVertex.push_back((double)xyz[1]);
1763     projVertex.push_back((double)xyz[2]);
1764     std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
1765     if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end()) {
1766       MESSAGE("Found enforced vertex @ " << xyz[0] << ", " << xyz[1] << ", " << xyz[2]);
1767       BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
1768       BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
1769       for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
1770         currentEnfVertex = (*enfListIt);
1771         if (currentEnfVertex->grpName != "") {
1772           bool groupDone = false;
1773           SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
1774           MESSAGE("currentEnfVertex->grpName: " << currentEnfVertex->grpName);
1775           MESSAGE("Parsing the groups of the mesh");
1776           while (grIt->more()) {
1777             SMESH_Group * group = grIt->next();
1778             if ( !group ) continue;
1779             MESSAGE("Group: " << group->GetName());
1780             SMESHDS_GroupBase* groupDS = group->GetGroupDS();
1781             if ( !groupDS ) continue;
1782             MESSAGE("group->SMDSGroup().GetType(): " << (groupDS->GetType()));
1783             MESSAGE("group->SMDSGroup().GetType()==SMDSAbs_Node: " << (groupDS->GetType()==SMDSAbs_Node));
1784             MESSAGE("currentEnfVertex->grpName.compare(group->GetStoreName())==0: " << (currentEnfVertex->grpName.compare(group->GetName())==0));
1785             if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
1786               SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
1787               aGroupDS->SMDSGroup().Add(nodes[iv]);
1788               MESSAGE("Node ID: " << nodes[iv]->GetID());
1789               // How can I inform the hypothesis ?
1790 //                 _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
1791               groupDone = true;
1792               MESSAGE("Successfully added enforced vertex to existing group " << currentEnfVertex->grpName);
1793               break;
1794             }
1795           }
1796           if (!groupDone)
1797           {
1798             int groupId;
1799             SMESH_Group* aGroup = aMesh.AddGroup(SMDSAbs_Node, currentEnfVertex->grpName.c_str(), groupId);
1800             aGroup->SetName( currentEnfVertex->grpName.c_str() );
1801             SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1802             aGroupDS->SMDSGroup().Add(nodes[iv]);
1803             MESSAGE("Successfully created enforced vertex group " << currentEnfVertex->grpName);
1804             groupDone = true;
1805           }
1806           if (!groupDone)
1807             throw SALOME_Exception(LOCALIZED("An enforced vertex node was not added to a group"));
1808         }
1809         else
1810           MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
1811       }
1812     }
1813
1814
1815     // internal point are tagged to zero
1816     if(tag > 0 && tag <= pmap.Extent() ){
1817       meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1818       tags[iv] = false;
1819     } else {
1820       tags[iv] = true;
1821     }
1822   }
1823
1824   /* enumerate edges */
1825   for(int it=1;it<=ne;it++) {
1826
1827     mesh_get_edge_vertices(msh, it, vtx);
1828     mesh_get_edge_tag(msh, it, &tag);
1829     if (tags[vtx[0]]) {
1830       Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1831       tags[vtx[0]] = false;
1832     };
1833     if (tags[vtx[1]]) {
1834       Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1835       tags[vtx[1]] = false;
1836     };
1837     SMDS_MeshEdge* edg = helper.AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1838     meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1839   }
1840
1841   /* enumerate triangles */
1842   for(int it=1;it<=nt;it++) {
1843     mesh_get_triangle_vertices(msh, it, vtx);
1844     SMDS_MeshFace* tri = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1845     mesh_get_triangle_tag(msh, it, &tag);
1846     meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1847     if (tags[vtx[0]]) {
1848       meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1849       tags[vtx[0]] = false;
1850     };
1851     if (tags[vtx[1]]) {
1852       meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1853       tags[vtx[1]] = false;
1854     };
1855     if (tags[vtx[2]]) {
1856       meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1857       tags[vtx[2]] = false;
1858     };
1859   }
1860
1861   /* enumerate quadrangles */
1862   for(int it=1;it<=nq;it++) {
1863     mesh_get_quadrangle_vertices(msh, it, vtx);
1864     SMDS_MeshFace* quad = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1865     mesh_get_quadrangle_tag(msh, it, &tag);
1866     meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1867     if (tags[vtx[0]]) {
1868       meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1869       tags[vtx[0]] = false;
1870     };
1871     if (tags[vtx[1]]) {
1872       meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1873       tags[vtx[1]] = false;
1874     };
1875     if (tags[vtx[2]]) {
1876       meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1877       tags[vtx[2]] = false;
1878     };
1879     if (tags[vtx[3]]) {
1880       meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1881       tags[vtx[3]] = false;
1882     };
1883   }
1884
1885   // SetIsAlwaysComputed( true ) to sub-meshes of degenerated EDGEs
1886   TopLoc_Location loc; double f,l;
1887   for (int i = 1; i <= emap.Extent(); i++)
1888     if ( BRep_Tool::Curve(TopoDS::Edge( emap( i )), loc, f,l).IsNull() )
1889       if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
1890         sm->SetIsAlwaysComputed( true );
1891   for (int i = 1; i <= pmap.Extent(); i++)
1892       if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( pmap( i )))
1893         if ( !sm->IsMeshComputed() )
1894           sm->SetIsAlwaysComputed( true );
1895
1896   if ( needMerge )
1897   {
1898     set< SMESH_subMesh* >::iterator smIt = edgeSubmeshes.begin();
1899     for ( ; smIt != edgeSubmeshes.end(); ++smIt )
1900     {
1901       SMESH_subMesh* sm = *smIt;
1902       SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
1903                                                                    /*complexFirst=*/false);
1904       TIDSortedNodeSet nodesOnEdge;
1905       double mergeTolerance = 1e-7, tol;
1906       while ( subsmIt->more() )
1907       {
1908         // get nodes from an edge
1909         sm = subsmIt->next();
1910         if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() )
1911         {
1912           SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1913           while ( nIt->more() )
1914             nodesOnEdge.insert( nIt->next() );
1915         }
1916         // get max tolerance
1917         TopoDS_Shape subShape = sm->GetSubShape();
1918         if ( subShape.ShapeType() == TopAbs_EDGE )
1919           tol = BRep_Tool::Tolerance( TopoDS::Edge( subShape ));
1920         else
1921           tol = BRep_Tool::Tolerance( TopoDS::Vertex( subShape ));
1922         mergeTolerance = Max( tol, mergeTolerance );
1923       }
1924       // find nodes to merge
1925       SMESH_MeshEditor::TListOfListOfNodes nodeGroupsToMerge;
1926       SMESH_MeshEditor editor( &aMesh );
1927       editor.FindCoincidentNodes( nodesOnEdge, mergeTolerance*2, nodeGroupsToMerge );
1928       // merge
1929       editor.MergeNodes( nodeGroupsToMerge );
1930     }
1931   }
1932
1933   delete [] nodes;
1934
1935   /* release the mesh object */
1936   blsurf_data_regain_mesh(bls, msh);
1937
1938   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1939 #ifndef WNT
1940   if ( oldFEFlags > 0 )
1941     feenableexcept( oldFEFlags );
1942   feclearexcept( FE_ALL_EXCEPT );
1943 #endif
1944
1945   /*
1946   std::cout << "FacesWithSizeMap" << std::endl;
1947   FacesWithSizeMap.Statistics(std::cout);
1948   std::cout << "EdgesWithSizeMap" << std::endl;
1949   EdgesWithSizeMap.Statistics(std::cout);
1950   std::cout << "VerticesWithSizeMap" << std::endl;
1951   VerticesWithSizeMap.Statistics(std::cout);
1952   std::cout << "FacesWithEnforcedVertices" << std::endl;
1953   FacesWithEnforcedVertices.Statistics(std::cout);
1954   */
1955
1956   MESSAGE("END OF BLSURFPlugin_BLSURF::Compute()");
1957   return true;
1958 }
1959
1960 #ifdef WITH_SMESH_CANCEL_COMPUTE
1961 void BLSURFPlugin_BLSURF::CancelCompute()
1962 {
1963   _compute_canceled = true;
1964 }
1965 #endif
1966
1967 //=============================================================================
1968 /*!
1969  *  SetNodeOnEdge
1970  */
1971 //=============================================================================
1972
1973 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1974   const TopoDS_Edge edge = TopoDS::Edge(ed);
1975
1976   gp_Pnt pnt(node->X(), node->Y(), node->Z());
1977
1978   Standard_Real p0 = 0.0;
1979   Standard_Real p1 = 1.0;
1980   TopLoc_Location loc;
1981   Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1982
1983   if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1984   GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1985
1986   double pa = 0.;
1987   if ( proj.NbPoints() > 0 )
1988   {
1989     pa = (double)proj.LowerDistanceParameter();
1990     // Issue 0020656. Move node if it is too far from edge
1991     gp_Pnt curve_pnt = curve->Value( pa );
1992     double dist2 = pnt.SquareDistance( curve_pnt );
1993     double tol = BRep_Tool::Tolerance( edge );
1994     if ( 1e-14 < dist2 && dist2 <= 1000*tol ) // large enough and within tolerance
1995     {
1996       curve_pnt.Transform( loc );
1997       meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1998     }
1999   }
2000 //   GProp_GProps LProps;
2001 //   BRepGProp::LinearProperties(ed, LProps);
2002 //   double lg = (double)LProps.Mass();
2003
2004   meshDS->SetNodeOnEdge(node, edge, pa);
2005 }
2006
2007 //=============================================================================
2008 /*!
2009  *
2010  */
2011 //=============================================================================
2012
2013 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
2014 {
2015   return save;
2016 }
2017
2018 //=============================================================================
2019 /*!
2020  *
2021  */
2022 //=============================================================================
2023
2024 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
2025 {
2026   return load;
2027 }
2028
2029 //=============================================================================
2030 /*!
2031  *
2032  */
2033 //=============================================================================
2034
2035 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
2036 {
2037   return hyp.SaveTo( save );
2038 }
2039
2040 //=============================================================================
2041 /*!
2042  *
2043  */
2044 //=============================================================================
2045
2046 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
2047 {
2048   return hyp.LoadFrom( load );
2049 }
2050
2051 /* Curve definition function See cad_curv_t in file distene/cad.h for
2052  * more information.
2053  * NOTE : if when your CAD systems evaluates second
2054  * order derivatives it also computes first order derivatives and
2055  * function evaluation, you can optimize this example by making only
2056  * one CAD call and filling the necessary uv, dt, dtt arrays.
2057  */
2058 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
2059 {
2060   /* t is given. It contains the t (time) 1D parametric coordintaes
2061      of the point PreCAD/BLSurf is querying on the curve */
2062
2063   /* user_data identifies the edge PreCAD/BLSurf is querying
2064    * (see cad_edge_new later in this example) */
2065   const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
2066
2067   if (uv){
2068    /* BLSurf is querying the function evaluation */
2069     gp_Pnt2d P;
2070     P=pargeo->Value(t);
2071     uv[0]=P.X(); uv[1]=P.Y();
2072   }
2073
2074   if(dt) {
2075    /* query for the first order derivatives */
2076     gp_Vec2d V1;
2077     V1=pargeo->DN(t,1);
2078     dt[0]=V1.X(); dt[1]=V1.Y();
2079   }
2080
2081   if(dtt){
2082     /* query for the second order derivatives */
2083     gp_Vec2d V2;
2084     V2=pargeo->DN(t,2);
2085     dtt[0]=V2.X(); dtt[1]=V2.Y();
2086   }
2087
2088   return STATUS_OK;
2089 }
2090
2091 /* Surface definition function.
2092  * See cad_surf_t in file distene/cad.h for more information.
2093  * NOTE : if when your CAD systems evaluates second order derivatives it also
2094  * computes first order derivatives and function evaluation, you can optimize
2095  * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
2096  * arrays.
2097  */
2098 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
2099                   real *duu, real *duv, real *dvv, void *user_data)
2100 {
2101   /* uv[2] is given. It contains the u,v coordinates of the point
2102    * PreCAD/BLSurf is querying on the surface */
2103
2104   /* user_data identifies the face PreCAD/BLSurf is querying (see
2105    * cad_face_new later in this example)*/
2106   const Geom_Surface* geometry = (const Geom_Surface*) user_data;
2107
2108   if(xyz){
2109    gp_Pnt P;
2110    P=geometry->Value(uv[0],uv[1]);   // S.D0(U,V,P);
2111    xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
2112   }
2113
2114   if(du && dv){
2115     gp_Pnt P;
2116     gp_Vec D1U,D1V;
2117
2118     geometry->D1(uv[0],uv[1],P,D1U,D1V);
2119     du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
2120     dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
2121   }
2122
2123   if(duu && duv && dvv){
2124
2125     gp_Pnt P;
2126     gp_Vec D1U,D1V;
2127     gp_Vec D2U,D2V,D2UV;
2128
2129     geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
2130     duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
2131     duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
2132     dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
2133   }
2134
2135   return STATUS_OK;
2136 }
2137
2138
2139 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
2140 {
2141   if (face_id == 1) {
2142     if (my_u_min > uv[0]) {
2143       my_u_min = uv[0];
2144     }
2145     if (my_v_min > uv[1]) {
2146       my_v_min = uv[1];
2147     }
2148     if (my_u_max < uv[0]) {
2149       my_u_max = uv[0];
2150     }
2151     if (my_v_max < uv[1]) {
2152       my_v_max = uv[1];
2153     }
2154   }
2155   //MESSAGE("size_on_surface")
2156   if (FaceId2PythonSmp.count(face_id) != 0){
2157     //MESSAGE("A size map is used to calculate size on face : "<<face_id)
2158     PyObject * pyresult = NULL;
2159     PyObject* new_stderr = NULL;
2160     assert(Py_IsInitialized());
2161     PyGILState_STATE gstate;
2162     gstate = PyGILState_Ensure();
2163     pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
2164     double result;
2165     if ( pyresult == NULL){
2166       fflush(stderr);
2167       string err_description="";
2168       new_stderr = newPyStdOut(err_description);
2169       PySys_SetObject((char*)"stderr", new_stderr);
2170       PyErr_Print();
2171       PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
2172       Py_DECREF(new_stderr);
2173       MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
2174       result = *((double*)user_data);
2175       }
2176     else {
2177       result = PyFloat_AsDouble(pyresult);
2178       Py_DECREF(pyresult);
2179     }
2180    // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
2181     *size = result;
2182     PyGILState_Release(gstate);
2183   }
2184   else if (FaceIndex2ClassAttractor.count(face_id) !=0 && !FaceIndex2ClassAttractor[face_id]->Empty()){
2185 //    MESSAGE("attractor used on face :"<<face_id)
2186     // MESSAGE("List of attractor is not empty")
2187     // MESSAGE("Attractor empty : "<< FaceIndex2ClassAttractor[face_id]->Empty())
2188     double result = FaceIndex2ClassAttractor[face_id]->GetSize(uv[0],uv[1]);
2189     *size = result;
2190    // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
2191   }
2192   else {
2193     // MESSAGE("List of attractor is empty !!!")
2194     *size = *((double*)user_data);
2195   }
2196   return STATUS_OK;
2197 }
2198
2199 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
2200 {
2201   if (EdgeId2PythonSmp.count(edge_id) != 0){
2202     PyObject * pyresult = NULL;
2203     PyObject* new_stderr = NULL;
2204     assert(Py_IsInitialized());
2205     PyGILState_STATE gstate;
2206     gstate = PyGILState_Ensure();
2207     pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
2208     double result;
2209     if ( pyresult == NULL){
2210       fflush(stderr);
2211       string err_description="";
2212       new_stderr = newPyStdOut(err_description);
2213       PySys_SetObject((char*)"stderr", new_stderr);
2214       PyErr_Print();
2215       PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
2216       Py_DECREF(new_stderr);
2217       MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
2218       result = *((double*)user_data);
2219       }
2220     else {
2221       result = PyFloat_AsDouble(pyresult);
2222       Py_DECREF(pyresult);
2223     }
2224     *size = result;
2225     PyGILState_Release(gstate);
2226   }
2227   else {
2228     *size = *((double*)user_data);
2229   }
2230   return STATUS_OK;
2231 }
2232
2233 status_t size_on_vertex(integer point_id, real *size, void *user_data)
2234 {
2235   if (VertexId2PythonSmp.count(point_id) != 0){
2236     PyObject * pyresult = NULL;
2237     PyObject* new_stderr = NULL;
2238     assert(Py_IsInitialized());
2239     PyGILState_STATE gstate;
2240     gstate = PyGILState_Ensure();
2241     pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
2242     double result;
2243     if ( pyresult == NULL){
2244       fflush(stderr);
2245       string err_description="";
2246       new_stderr = newPyStdOut(err_description);
2247       PySys_SetObject((char*)"stderr", new_stderr);
2248       PyErr_Print();
2249       PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
2250       Py_DECREF(new_stderr);
2251       MESSAGE("Can't evaluate f()" << " error is " << err_description);
2252       result = *((double*)user_data);
2253       }
2254     else {
2255       result = PyFloat_AsDouble(pyresult);
2256       Py_DECREF(pyresult);
2257     }
2258     *size = result;
2259     PyGILState_Release(gstate);
2260   }
2261   else {
2262     *size = *((double*)user_data);
2263   }
2264  return STATUS_OK;
2265 }
2266
2267 /*
2268  * The following function will be called for PreCAD/BLSurf message
2269  * printing.  See context_set_message_callback (later in this
2270  * template) for how to set user_data.
2271  */
2272 status_t message_cb(message_t *msg, void *user_data)
2273 {
2274   integer errnumber = 0;
2275   char *desc;
2276   message_get_number(msg, &errnumber);
2277   message_get_description(msg, &desc);
2278   string err( desc );
2279   if ( errnumber < 0 || err.find("license") != string::npos ) {
2280     string * error = (string*)user_data;
2281 //   if ( !error->empty() )
2282 //     *error += "\n";
2283     // remove ^A from the tail
2284     int len = strlen( desc );
2285     while (len > 0 && desc[len-1] != '\n')
2286       len--;
2287     error->append( desc, len );
2288   }
2289   else {
2290       std::cout << desc << std::endl;
2291   }
2292   return STATUS_OK;
2293 }
2294
2295 /* This is the interrupt callback. PreCAD/BLSurf will call this
2296  * function regularily. See the file distene/interrupt.h
2297  */
2298 status_t interrupt_cb(integer *interrupt_status, void *user_data)
2299 {
2300   integer you_want_to_continue = 1;
2301 #ifdef WITH_SMESH_CANCEL_COMPUTE
2302   BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
2303   you_want_to_continue = !tmp->computeCanceled();
2304 #endif
2305
2306   if(you_want_to_continue)
2307   {
2308     *interrupt_status = INTERRUPT_CONTINUE;
2309     return STATUS_OK;
2310   }
2311   else /* you want to stop BLSurf */
2312   {
2313     *interrupt_status = INTERRUPT_STOP;
2314     return STATUS_ERROR;
2315   }
2316 }
2317
2318 //=============================================================================
2319 /*!
2320  *
2321  */
2322 //=============================================================================
2323 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh&         aMesh,
2324                                    const TopoDS_Shape& aShape,
2325                                    MapShapeNbElems&    aResMap)
2326 {
2327   int    _physicalMesh  = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
2328   double _phySize       = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
2329   //int    _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
2330   //double _angleMeshS    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
2331   double _angleMeshC    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
2332   bool   _quadAllowed   = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
2333   if(_hypothesis) {
2334     _physicalMesh  = (int) _hypothesis->GetPhysicalMesh();
2335     _phySize       = _hypothesis->GetPhySize();
2336     //_geometricMesh = (int) hyp->GetGeometricMesh();
2337     //_angleMeshS    = hyp->GetAngleMeshS();
2338     _angleMeshC    = _hypothesis->GetAngleMeshC();
2339     _quadAllowed   = _hypothesis->GetQuadAllowed();
2340   }
2341
2342   bool IsQuadratic = false;
2343
2344   // ----------------
2345   // evaluate 1D
2346   // ----------------
2347   TopTools_DataMapOfShapeInteger EdgesMap;
2348   double fullLen = 0.0;
2349   double fullNbSeg = 0;
2350   for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2351     TopoDS_Edge E = TopoDS::Edge( exp.Current() );
2352     if( EdgesMap.IsBound(E) )
2353       continue;
2354     SMESH_subMesh *sm = aMesh.GetSubMesh(E);
2355     double aLen = SMESH_Algo::EdgeLength(E);
2356     fullLen += aLen;
2357     int nb1d = 0;
2358     if(_physicalMesh==1) {
2359        nb1d = (int)( aLen/_phySize + 1 );
2360     }
2361     else {
2362       // use geometry
2363       double f,l;
2364       Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
2365       double fullAng = 0.0;
2366       double dp = (l-f)/200;
2367       gp_Pnt P1,P2,P3;
2368       C->D0(f,P1);
2369       C->D0(f+dp,P2);
2370       gp_Vec V1(P1,P2);
2371       for(int j=2; j<=200; j++) {
2372         C->D0(f+dp*j,P3);
2373         gp_Vec V2(P2,P3);
2374         fullAng += fabs(V1.Angle(V2));
2375         V1 = V2;
2376         P2 = P3;
2377       }
2378       nb1d = (int)( fullAng/_angleMeshC + 1 );
2379     }
2380     fullNbSeg += nb1d;
2381     std::vector<int> aVec(SMDSEntity_Last);
2382     for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2383     if( IsQuadratic > 0 ) {
2384       aVec[SMDSEntity_Node] = 2*nb1d - 1;
2385       aVec[SMDSEntity_Quad_Edge] = nb1d;
2386     }
2387     else {
2388       aVec[SMDSEntity_Node] = nb1d - 1;
2389       aVec[SMDSEntity_Edge] = nb1d;
2390     }
2391     aResMap.insert(std::make_pair(sm,aVec));
2392     EdgesMap.Bind(E,nb1d);
2393   }
2394   double ELen = fullLen/fullNbSeg;
2395   // ----------------
2396   // evaluate 2D
2397   // ----------------
2398   // try to evaluate as in MEFISTO
2399   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2400     TopoDS_Face F = TopoDS::Face( exp.Current() );
2401     SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2402     GProp_GProps G;
2403     BRepGProp::SurfaceProperties(F,G);
2404     double anArea = G.Mass();
2405     int nb1d = 0;
2406     std::vector<int> nb1dVec;
2407     for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
2408       int nbSeg = EdgesMap.Find(exp1.Current());
2409       nb1d += nbSeg;
2410       nb1dVec.push_back( nbSeg );
2411     }
2412     int nbQuad = 0;
2413     int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
2414     int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
2415     if ( _quadAllowed )
2416     {
2417       if ( nb1dVec.size() == 4 ) // quadrangle geom face
2418       {
2419         int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
2420         nbQuad = n1 * n2;
2421         nbNodes = (n1 + 1) * (n2 + 1);
2422         nbTria = 0;
2423       }
2424       else
2425       {
2426         nbTria = nbQuad = nbTria / 3 + 1;
2427       }
2428     }
2429     std::vector<int> aVec(SMDSEntity_Last,0);
2430     if( IsQuadratic ) {
2431       int nb1d_in = (nbTria*3 - nb1d) / 2;
2432       aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
2433       aVec[SMDSEntity_Quad_Triangle] = nbTria;
2434       aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
2435     }
2436     else {
2437       aVec[SMDSEntity_Node] = nbNodes;
2438       aVec[SMDSEntity_Triangle] = nbTria;
2439       aVec[SMDSEntity_Quadrangle] = nbQuad;
2440     }
2441     aResMap.insert(std::make_pair(sm,aVec));
2442   }
2443
2444   // ----------------
2445   // evaluate 3D
2446   // ----------------
2447   GProp_GProps G;
2448   BRepGProp::VolumeProperties(aShape,G);
2449   double aVolume = G.Mass();
2450   double tetrVol = 0.1179*ELen*ELen*ELen;
2451   int nbVols  = int(aVolume/tetrVol);
2452   int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
2453   std::vector<int> aVec(SMDSEntity_Last);
2454   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2455   if( IsQuadratic ) {
2456     aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
2457     aVec[SMDSEntity_Quad_Tetra] = nbVols;
2458   }
2459   else {
2460     aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
2461     aVec[SMDSEntity_Tetra] = nbVols;
2462   }
2463   SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2464   aResMap.insert(std::make_pair(sm,aVec));
2465
2466   return true;
2467 }
2468
2469 //=============================================================================
2470 /*!
2471  *  Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
2472  *  Following line was added:
2473  *        myExtrem.Perform(P);
2474  */
2475 //=============================================================================
2476 void  BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
2477                     const TopoDS_Face& face,
2478                     const gp_Pnt& P,
2479                     const Standard_Real Tol)
2480 {
2481   //-- Voir BRepExtrema_ExtPF.cxx
2482   BRepAdaptor_Surface Surf(face);
2483   Standard_Real U1, U2, V1, V2;
2484   BRepTools::UVBounds(face, U1, U2, V1, V2);
2485   Extrema_ExtPS myExtrem;
2486   myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
2487   myExtrem.Perform(P);
2488   //----------------------------------------------------------
2489   //-- On cherche le point le plus proche , PUIS
2490   //-- On le classifie.
2491   Standard_Integer nbv    = 0; // xpu
2492   Standard_Real MaxDist   =  RealLast();
2493   Standard_Integer indice = 0;
2494   if (myExtrem.IsDone()) {
2495     nbv = myExtrem.NbExt();
2496     for (Standard_Integer i = 1; i <= nbv; i++) {
2497 #if OCC_VERSION_LARGE > 0x06040000 // Porting to OCCT6.5.1
2498       Standard_Real d = myExtrem.SquareDistance(i);
2499 #else
2500       Standard_Real d = myExtrem.Value(i);
2501       d = Abs(d);
2502 #endif
2503       if (d <= MaxDist) {
2504         MaxDist = d;
2505         indice = i;
2506       }
2507     }
2508   }
2509   if (indice) {
2510     gp_Pnt2d Puv;
2511     Standard_Real U1,U2;
2512     myExtrem.Point(indice).Parameter(U1, U2);
2513     Puv.SetCoord(U1, U2);
2514     fc->Perform(face, Puv, Tol);
2515   }
2516   else {
2517     fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
2518     //-- le resultat est TopAbs_OUT;
2519   }
2520 }