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