]> SALOME platform Git repositories - plugins/blsurfplugin.git/blob - src/BLSURFPlugin/BLSURFPlugin_BLSURF.cxx
Salome HOME
Merge from BR_ATTRACTOR_GEO
[plugins/blsurfplugin.git] / src / BLSURFPlugin / BLSURFPlugin_BLSURF.cxx
1 //  Copyright (C) 2007-2010  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 && (AttShape.ShapeType() == TopAbs_VERTEX || AttShape.ShapeType() == TopAbs_EDGE)){
903           HasSizeMapOnFace = true;
904           if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape)) ) {
905                 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
906               }
907               else {
908                 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
909 //                 MESSAGE("Face with key " << key << " already in map");
910               }
911               FaceId2ClassAttractor[key] = AtIt->second;
912         }
913         else{
914           MESSAGE("Wrong shape type !!")
915         }
916
917       }
918     }
919
920     
921
922     //
923     // Enforced Vertices
924     //
925     MESSAGE("Setting Enforced Vertices");
926     const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
927     BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
928     for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
929       if ( !enfIt->second.empty() ) {
930         GeomShape = entryToShape(enfIt->first);
931         GeomType  = GeomShape.ShapeType();
932         // Group Management
933         if (GeomType == TopAbs_COMPOUND){
934           for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
935             if (it.Value().ShapeType() == TopAbs_FACE){
936               HasSizeMapOnFace = true;
937               createEnforcedVertexOnFace(it.Value(), enfIt->second);
938             }
939           }
940         }
941             
942         if (GeomType == TopAbs_FACE){
943           HasSizeMapOnFace = true;
944           createEnforcedVertexOnFace(GeomShape, enfIt->second);
945         }
946       }
947     }
948
949 //    if (HasSizeMapOnFace){
950     // In all size map cases (hphy_flag = 2), at least map on face must be defined
951     MESSAGE("Setting Size Map on FACES ");
952     blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
953 //    }
954
955     if (HasSizeMapOnEdge){
956       MESSAGE("Setting Size Map on EDGES ");
957       blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
958     }
959     if (HasSizeMapOnVertex){
960       MESSAGE("Setting Size Map on VERTICES ");
961       blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
962     }
963   }
964   blsurf_set_param(bls, "hphydef",           to_string(_phySize).c_str());
965   blsurf_set_param(bls, "hgeo_flag",         to_string(_geometricMesh).c_str());
966   blsurf_set_param(bls, "relax_size",        _decimesh ? "0": to_string(_geometricMesh).c_str());
967   blsurf_set_param(bls, "angle_meshs",       to_string(_angleMeshS).c_str());
968   blsurf_set_param(bls, "angle_meshc",       to_string(_angleMeshC).c_str());
969   blsurf_set_param(bls, "gradation",         to_string(_gradation).c_str());
970   blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
971   blsurf_set_param(bls, "element",           _quadAllowed ? "q1.0" : "p1");
972   blsurf_set_param(bls, "verb",              to_string(_verb).c_str());
973 }
974
975 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
976 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
977                   real *duu, real *duv, real *dvv, void *user_data);
978 status_t message_cb(message_t *msg, void *user_data);
979 status_t interrupt_cb(integer *interrupt_status, void *user_data);
980
981 //=============================================================================
982 /*!
983  *
984  */
985 //=============================================================================
986
987 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
988
989   MESSAGE("BLSURFPlugin_BLSURF::Compute");
990
991 //   if (aShape.ShapeType() == TopAbs_COMPOUND) {
992 //     MESSAGE("  the shape is a COMPOUND");
993 //   }
994 //   else {
995 //     MESSAGE("  the shape is UNKNOWN");
996 //   };
997
998   // Fix problem with locales
999   Kernel_Utils::Localizer aLocalizer;
1000
1001   /* create a distene context (generic object) */
1002   status_t status = STATUS_ERROR;
1003   
1004   context_t *ctx =  context_new();
1005   
1006   /* Set the message callback in the working context */
1007   context_set_message_callback(ctx, message_cb, &_comment);
1008 #ifdef WITH_SMESH_CANCEL_COMPUTE
1009   _compute_canceled = false;
1010   context_set_interrupt_callback(ctx, interrupt_cb, this);
1011 #else
1012   context_set_interrupt_callback(ctx, interrupt_cb, NULL);
1013 #endif
1014
1015   /* create the CAD object we will work on. It is associated to the context ctx. */
1016   cad_t *c = cad_new(ctx);
1017
1018   blsurf_session_t *bls = blsurf_session_new(ctx);
1019
1020   FacesWithSizeMap.Clear();
1021   FaceId2SizeMap.clear();
1022   FaceId2ClassAttractor.clear();
1023   FaceIndex2ClassAttractor.clear();
1024   EdgesWithSizeMap.Clear();
1025   EdgeId2SizeMap.clear();
1026   VerticesWithSizeMap.Clear();
1027   VertexId2SizeMap.clear();
1028
1029   MESSAGE("BEGIN SetParameters");
1030   SetParameters(_hypothesis, bls, aMesh);
1031   MESSAGE("END SetParameters");
1032
1033   /* Now fill the CAD object with data from your CAD
1034    * environement. This is the most complex part of a successfull
1035    * integration.
1036    */
1037
1038   // needed to prevent the opencascade memory managmement from freeing things
1039   vector<Handle(Geom2d_Curve)> curves;
1040   vector<Handle(Geom_Surface)> surfaces;
1041
1042   surfaces.resize(0);
1043   curves.resize(0);
1044   
1045   TopTools_IndexedMapOfShape fmap;
1046   TopTools_IndexedMapOfShape emap;
1047   TopTools_IndexedMapOfShape pmap;
1048   
1049   fmap.Clear();
1050   FaceId2PythonSmp.clear();
1051   emap.Clear();
1052   EdgeId2PythonSmp.clear();
1053   pmap.Clear();
1054   VertexId2PythonSmp.clear();
1055
1056   assert(Py_IsInitialized());
1057   PyGILState_STATE gstate;
1058   gstate = PyGILState_Ensure();
1059
1060   string theSizeMapStr;
1061   
1062   /****************************************************************************************
1063                                   FACES
1064   *****************************************************************************************/
1065   int iface = 0;
1066   string bad_end = "return";
1067   int faceKey = -1;
1068   int ienf = 0;
1069   BLSURFPlugin_Attractor myAttractor;
1070   for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
1071     TopoDS_Face f=TopoDS::Face(face_iter.Current());
1072
1073     // make INTERNAL face oriented FORWARD (issue 0020993)
1074     if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
1075       f.Orientation(TopAbs_FORWARD);
1076     
1077     if (fmap.FindIndex(f) > 0)
1078       continue;
1079
1080     fmap.Add(f);
1081     iface++;
1082     surfaces.push_back(BRep_Tool::Surface(f));
1083     
1084     /* create an object representing the face for blsurf */
1085     /* where face_id is an integer identifying the face.
1086      * surf_function is the function that defines the surface
1087      * (For this face, it will be called by blsurf with your_face_object_ptr
1088      * as last parameter.
1089      */
1090     cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
1091     
1092     /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
1093     cad_face_set_tag(fce, iface);
1094       
1095     /* Set face orientation (optional if you want a well oriented output mesh)*/
1096     if(f.Orientation() != TopAbs_FORWARD){
1097       cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
1098     } else {
1099       cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
1100     }
1101     
1102     if (HasSizeMapOnFace){
1103 //       MESSAGE("A size map is defined on a face")
1104 //       std::cout << "A size map is defined on a face" << std::endl;
1105       // Classic size map
1106       faceKey = FacesWithSizeMap.FindIndex(f);
1107       
1108       
1109       if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
1110         MESSAGE("A size map is defined on face :"<<faceKey)
1111         theSizeMapStr = FaceId2SizeMap[faceKey];
1112         // check if function ends with "return"
1113         if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1114           continue;
1115         // Expr To Python function, verification is performed at validation in GUI
1116         PyObject * obj = NULL;
1117         obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1118         Py_DECREF(obj);
1119         PyObject * func = NULL;
1120         func = PyObject_GetAttrString(main_mod, "f");
1121         FaceId2PythonSmp[iface]=func;
1122         FaceId2SizeMap.erase(faceKey);
1123       }
1124       
1125       // Specific size map = Attractor
1126       std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
1127       int iatt=0;
1128       for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
1129         if (attractor_iter->first == faceKey) {
1130           MESSAGE("Face indice: " << iface);
1131           MESSAGE("Adding attractor");
1132           
1133           double xyzCoords[3]  = {attractor_iter->second[2],
1134                                   attractor_iter->second[3],
1135                                   attractor_iter->second[4]};
1136           
1137           MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1138           gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1139           BRepClass_FaceClassifier scl(f,P,1e-7);
1140           // scl.Perform() is bugged. The function was rewritten
1141 //          scl.Perform();
1142           BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1143           TopAbs_State result = scl.State();
1144           MESSAGE("Position of point on face: "<<result);
1145           if ( result == TopAbs_OUT )
1146               MESSAGE("Point is out of face: node is not created");
1147           if ( result == TopAbs_UNKNOWN )
1148               MESSAGE("Point position on face is unknown: node is not created");
1149           if ( result == TopAbs_ON )
1150               MESSAGE("Point is on border of face: node is not created");
1151           if ( result == TopAbs_IN )
1152           {
1153             // Point is inside face and not on border
1154             MESSAGE("Point is in face: node is created");
1155             double uvCoords[2]   = {attractor_iter->second[0],attractor_iter->second[1]};
1156             iatt++;
1157             MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << iatt);
1158             cad_point_t* point_p = cad_point_new(fce, iatt, uvCoords);
1159             cad_point_set_tag(point_p, iatt);
1160           }
1161           FaceId2AttractorCoords.erase(faceKey);
1162         }
1163       }
1164       
1165       // Class Attractors
1166       std::map<int,BLSURFPlugin_Attractor* >::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
1167       if (clAttractor_iter != FaceId2ClassAttractor.end()){
1168           MESSAGE("Face indice: " << iface);
1169           MESSAGE("Adding attractor");
1170           FaceIndex2ClassAttractor[iface]=clAttractor_iter->second;
1171           FaceId2ClassAttractor.erase(clAttractor_iter);
1172         }
1173       }     
1174       
1175       // Enforced Vertices
1176       faceKey = FacesWithEnforcedVertices.FindIndex(f);
1177       std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1178       if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1179         MESSAGE("Some enforced vertices are defined");
1180         BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl;
1181         MESSAGE("Face indice: " << iface);
1182         MESSAGE("Adding enforced vertices");
1183         evl = evmIt->second;
1184         MESSAGE("Number of vertices to add: "<< evl.size());
1185         BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
1186         for (; evlIt != evl.end(); ++evlIt) {
1187           BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
1188           xyzCoords.push_back(evlIt->at(2));
1189           xyzCoords.push_back(evlIt->at(3));
1190           xyzCoords.push_back(evlIt->at(4));
1191           MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1192           gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1193           BRepClass_FaceClassifier scl(f,P,1e-7);
1194           // scl.Perform() is bugged. The function was rewritten
1195 //          scl.Perform();
1196           BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1197           TopAbs_State result = scl.State();
1198           MESSAGE("Position of point on face: "<<result);
1199           if ( result == TopAbs_OUT ) {
1200             MESSAGE("Point is out of face: node is not created");
1201             if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1202               EnfVertexCoords2ProjVertex.erase(xyzCoords);
1203               EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1204             }
1205           }
1206           if ( result == TopAbs_UNKNOWN ) {
1207             MESSAGE("Point position on face is unknown: node is not created");
1208             if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1209               EnfVertexCoords2ProjVertex.erase(xyzCoords);
1210               EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1211             }
1212           }
1213           if ( result == TopAbs_ON ) {
1214             MESSAGE("Point is on border of face: node is not created");
1215             if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1216               EnfVertexCoords2ProjVertex.erase(xyzCoords);
1217               EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1218             }
1219           }
1220           if ( result == TopAbs_IN )
1221           {
1222             // Point is inside face and not on border
1223             MESSAGE("Point is in face: node is created");
1224             double uvCoords[2]   = {evlIt->at(0),evlIt->at(1)};
1225             ienf++;
1226             MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1227             cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1228             cad_point_set_tag(point_p, ienf);
1229           }
1230         }
1231         FaceId2EnforcedVertexCoords.erase(faceKey);
1232       }
1233 //       else
1234 //         std::cout << "No enforced vertex defined" << std::endl;
1235 //     }
1236     
1237     
1238     /****************************************************************************************
1239                                     EDGES
1240                    now create the edges associated to this face
1241     *****************************************************************************************/
1242     int edgeKey = -1;
1243     for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1244       TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1245       int ic = emap.FindIndex(e);
1246       if (ic <= 0)
1247         ic = emap.Add(e);
1248
1249       double tmin,tmax;
1250       curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1251       
1252       if (HasSizeMapOnEdge){
1253         edgeKey = EdgesWithSizeMap.FindIndex(e);
1254         if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1255           MESSAGE("Sizemap defined on edge with index " << edgeKey);
1256           theSizeMapStr = EdgeId2SizeMap[edgeKey];
1257           if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1258             continue;
1259           // Expr To Python function, verification is performed at validation in GUI
1260           PyObject * obj = NULL;
1261           obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1262           Py_DECREF(obj);
1263           PyObject * func = NULL;
1264           func = PyObject_GetAttrString(main_mod, "f");
1265           EdgeId2PythonSmp[ic]=func;
1266           EdgeId2SizeMap.erase(edgeKey);
1267         }
1268       }
1269       
1270       /* attach the edge to the current blsurf face */
1271       cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1272       
1273       /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1274       cad_edge_set_tag(edg, ic);
1275       
1276       /* by default, an edge does not necessalry appear in the resulting mesh,
1277      unless the following property is set :
1278       */
1279       cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1280       
1281       /* by default an edge is a boundary edge */
1282       if (e.Orientation() == TopAbs_INTERNAL)
1283         cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1284
1285       int npts = 0;
1286       int ip1, ip2, *ip;
1287       gp_Pnt2d e0 = curves.back()->Value(tmin);
1288       gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1289       Standard_Real d1=0,d2=0;
1290       
1291       
1292       /****************************************************************************************
1293                                       VERTICES
1294       *****************************************************************************************/
1295       int vertexKey = -1;
1296       for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1297         TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1298         ++npts;
1299         if (npts == 1){
1300           ip = &ip1;
1301           d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1302         } else {
1303           ip = &ip2;
1304           d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1305         }
1306         *ip = pmap.FindIndex(v);
1307         if(*ip <= 0)
1308           *ip = pmap.Add(v);
1309         
1310         if (HasSizeMapOnVertex){
1311           vertexKey = VerticesWithSizeMap.FindIndex(v);
1312           if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1313             theSizeMapStr = VertexId2SizeMap[vertexKey];
1314             //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1315             if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1316               continue;
1317             // Expr To Python function, verification is performed at validation in GUI
1318             PyObject * obj = NULL;
1319             obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1320             Py_DECREF(obj);
1321             PyObject * func = NULL;
1322             func = PyObject_GetAttrString(main_mod, "f");
1323             VertexId2PythonSmp[*ip]=func;
1324             VertexId2SizeMap.erase(vertexKey);   // do not erase if using a vector
1325           }
1326         }
1327       }
1328       if (npts != 2) {
1329         // should not happen
1330         MESSAGE("An edge does not have 2 extremities.");
1331       } else {
1332         if (d1 < d2) {
1333           // This defines the curves extremity connectivity
1334           cad_edge_set_extremities(edg, ip1, ip2);
1335           /* set the tag (color) to the same value as the extremity id : */
1336           cad_edge_set_extremities_tag(edg, ip1, ip2);
1337         }
1338         else {
1339           cad_edge_set_extremities(edg, ip2, ip1);
1340           cad_edge_set_extremities_tag(edg, ip2, ip1);
1341         }
1342       }
1343     } // for edge
1344   } //for face
1345
1346
1347   PyGILState_Release(gstate);
1348
1349   blsurf_data_set_cad(bls, c);
1350
1351   std::cout << std::endl;
1352   std::cout << "Beginning of Surface Mesh generation" << std::endl;
1353   std::cout << std::endl;
1354
1355   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1356 #ifndef WNT
1357   feclearexcept( FE_ALL_EXCEPT );
1358   int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1359 #endif
1360
1361   try {
1362     OCC_CATCH_SIGNALS;
1363
1364     status = blsurf_compute_mesh(bls);
1365
1366   }
1367   catch ( std::exception& exc ) {
1368     _comment += exc.what();
1369   }
1370   catch (Standard_Failure& ex) {
1371     _comment += ex.DynamicType()->Name();
1372     if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1373       _comment += ": ";
1374       _comment += ex.GetMessageString();
1375     }
1376   }
1377   catch (...) {
1378     if ( _comment.empty() )
1379       _comment = "Exception in blsurf_compute_mesh()";
1380   }
1381   if ( status != STATUS_OK) {
1382     // There was an error while meshing
1383     blsurf_session_delete(bls);
1384     cad_delete(c);
1385     context_delete(ctx);
1386
1387     return error(_comment);
1388   }
1389
1390   std::cout << std::endl;
1391   std::cout << "End of Surface Mesh generation" << std::endl;
1392   std::cout << std::endl;
1393
1394   mesh_t *msh = NULL;
1395   blsurf_data_get_mesh(bls, &msh);
1396   if(!msh){
1397     blsurf_session_delete(bls);
1398     cad_delete(c);
1399     context_delete(ctx);
1400
1401     return error(_comment);
1402     //return false;
1403   }
1404
1405   /* retrieve mesh data (see distene/mesh.h) */
1406   integer nv, ne, nt, nq, vtx[4], tag;
1407   real xyz[3];
1408
1409   mesh_get_vertex_count(msh, &nv);
1410   mesh_get_edge_count(msh, &ne);
1411   mesh_get_triangle_count(msh, &nt);
1412   mesh_get_quadrangle_count(msh, &nq);
1413
1414
1415   SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1416   SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1417   bool* tags = new bool[nv+1];
1418
1419   /* enumerated vertices */
1420   for(int iv=1;iv<=nv;iv++) {
1421     mesh_get_vertex_coordinates(msh, iv, xyz);
1422     mesh_get_vertex_tag(msh, iv, &tag);
1423     // Issue 0020656. Use vertex coordinates
1424     if ( tag > 0 && tag <= pmap.Extent() ) {
1425       TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
1426       double tol = BRep_Tool::Tolerance( v );
1427       gp_Pnt p = BRep_Tool::Pnt( v );
1428       if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
1429         xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
1430       else
1431         tag = 0; // enforced or attracted vertex
1432     }
1433     nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1434
1435     // Create group of enforced vertices if requested
1436     if(_hypothesis) {
1437       BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
1438       projVertex.clear();
1439       projVertex.push_back((double)xyz[0]);
1440       projVertex.push_back((double)xyz[1]);
1441       projVertex.push_back((double)xyz[2]);
1442       std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
1443       if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end()) {
1444         MESSAGE("Found enforced vertex @ " << xyz[0] << ", " << xyz[1] << ", " << xyz[2])
1445         BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
1446         BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
1447         for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
1448           currentEnfVertex = (*enfListIt);
1449           if (currentEnfVertex->grpName != "") {
1450             bool groupDone = false;
1451             SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
1452             MESSAGE("currentEnfVertex->grpName: " << currentEnfVertex->grpName);
1453             MESSAGE("Parsing the groups of the mesh");
1454             while (grIt->more()) {
1455               SMESH_Group * group = grIt->next();
1456               if ( !group ) continue;
1457               MESSAGE("Group: " << group->GetName());
1458               SMESHDS_GroupBase* groupDS = group->GetGroupDS();
1459               if ( !groupDS ) continue;
1460               MESSAGE("group->SMDSGroup().GetType(): " << (groupDS->GetType()));
1461               MESSAGE("group->SMDSGroup().GetType()==SMDSAbs_Node: " << (groupDS->GetType()==SMDSAbs_Node));
1462               MESSAGE("currentEnfVertex->grpName.compare(group->GetStoreName())==0: " << (currentEnfVertex->grpName.compare(group->GetName())==0));
1463               if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
1464                 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
1465                 aGroupDS->SMDSGroup().Add(nodes[iv]);
1466                 MESSAGE("Node ID: " << nodes[iv]->GetID());
1467                 // How can I inform the hypothesis ?
1468 //                 _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
1469                 groupDone = true;
1470                 MESSAGE("Successfully added enforced vertex to existing group " << currentEnfVertex->grpName);
1471                 break;
1472               }
1473             }
1474             if (!groupDone)
1475             {
1476               int groupId;
1477               SMESH_Group* aGroup = aMesh.AddGroup(SMDSAbs_Node, currentEnfVertex->grpName.c_str(), groupId);
1478               aGroup->SetName( currentEnfVertex->grpName.c_str() );
1479               SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1480               aGroupDS->SMDSGroup().Add(nodes[iv]);
1481               MESSAGE("Successfully created enforced vertex group " << currentEnfVertex->grpName);
1482               groupDone = true;
1483             }
1484             if (!groupDone)
1485               throw SALOME_Exception(LOCALIZED("A enforced vertex node was not added to a group"));
1486           }
1487           else
1488             MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
1489         }
1490       }
1491     }
1492
1493
1494     // internal point are tagged to zero
1495     if(tag > 0 && tag <= pmap.Extent() ){
1496       meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1497       tags[iv] = false;
1498     } else {
1499       tags[iv] = true;
1500     }
1501   }
1502
1503   /* enumerate edges */
1504   for(int it=1;it<=ne;it++) {
1505     mesh_get_edge_vertices(msh, it, vtx);
1506     SMDS_MeshEdge* edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1507     mesh_get_edge_tag(msh, it, &tag);
1508
1509     if (tags[vtx[0]]) {
1510       Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1511       tags[vtx[0]] = false;
1512     };
1513     if (tags[vtx[1]]) {
1514       Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1515       tags[vtx[1]] = false;
1516     };
1517     meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1518
1519   }
1520
1521   /* enumerate triangles */
1522   for(int it=1;it<=nt;it++) {
1523     mesh_get_triangle_vertices(msh, it, vtx);
1524     SMDS_MeshFace* tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1525     mesh_get_triangle_tag(msh, it, &tag);
1526     meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1527     if (tags[vtx[0]]) {
1528       meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1529       tags[vtx[0]] = false;
1530     };
1531     if (tags[vtx[1]]) {
1532       meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1533       tags[vtx[1]] = false;
1534     };
1535     if (tags[vtx[2]]) {
1536       meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1537       tags[vtx[2]] = false;
1538     };
1539   }
1540
1541   /* enumerate quadrangles */
1542   for(int it=1;it<=nq;it++) {
1543     mesh_get_quadrangle_vertices(msh, it, vtx);
1544     SMDS_MeshFace* quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1545     mesh_get_quadrangle_tag(msh, it, &tag);
1546     meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1547     if (tags[vtx[0]]) {
1548       meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1549       tags[vtx[0]] = false;
1550     };
1551     if (tags[vtx[1]]) {
1552       meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1553       tags[vtx[1]] = false;
1554     };
1555     if (tags[vtx[2]]) {
1556       meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1557       tags[vtx[2]] = false;
1558     };
1559     if (tags[vtx[3]]) {
1560       meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1561       tags[vtx[3]] = false;
1562     };
1563   }
1564
1565   // SetIsAlwaysComputed( true ) to sub-meshes of degenerated EDGEs
1566   TopLoc_Location loc; double f,l;
1567   for (int i = 1; i <= emap.Extent(); i++)
1568     if ( BRep_Tool::Curve(TopoDS::Edge( emap( i )), loc, f,l).IsNull() )
1569       if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
1570         sm->SetIsAlwaysComputed( true );
1571
1572   delete [] nodes;
1573
1574   /* release the mesh object */
1575   blsurf_data_regain_mesh(bls, msh);
1576
1577   /* clean up everything */
1578   blsurf_session_delete(bls);
1579   cad_delete(c);
1580
1581   context_delete(ctx);
1582
1583   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1584 #ifndef WNT
1585   if ( oldFEFlags > 0 )
1586     feenableexcept( oldFEFlags );
1587   feclearexcept( FE_ALL_EXCEPT );
1588 #endif
1589
1590   /*  
1591   std::cout << "FacesWithSizeMap" << std::endl;
1592   FacesWithSizeMap.Statistics(std::cout);
1593   std::cout << "EdgesWithSizeMap" << std::endl;
1594   EdgesWithSizeMap.Statistics(std::cout);
1595   std::cout << "VerticesWithSizeMap" << std::endl;
1596   VerticesWithSizeMap.Statistics(std::cout);
1597   std::cout << "FacesWithEnforcedVertices" << std::endl;
1598   FacesWithEnforcedVertices.Statistics(std::cout);
1599   */
1600   
1601   return true;
1602 }
1603
1604 #ifdef WITH_SMESH_CANCEL_COMPUTE
1605 void BLSURFPlugin_BLSURF::CancelCompute()
1606 {
1607   _compute_canceled = true;
1608 }
1609 #endif
1610
1611 //=============================================================================
1612 /*!
1613  *  SetNodeOnEdge
1614  */
1615 //=============================================================================
1616
1617 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1618   const TopoDS_Edge edge = TopoDS::Edge(ed);
1619
1620   gp_Pnt pnt(node->X(), node->Y(), node->Z());
1621
1622   Standard_Real p0 = 0.0;
1623   Standard_Real p1 = 1.0;
1624   TopLoc_Location loc;
1625   Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1626
1627   if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1628   GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1629
1630   double pa = 0.;
1631   if ( proj.NbPoints() > 0 )
1632   {
1633     pa = (double)proj.LowerDistanceParameter();
1634     // Issue 0020656. Move node if it is too far from edge
1635     gp_Pnt curve_pnt = curve->Value( pa );
1636     double dist2 = pnt.SquareDistance( curve_pnt );
1637     double tol = BRep_Tool::Tolerance( edge );
1638     if ( 1e-12 < dist2 && dist2 <= 2*tol*tol ) // large enough and within tolerance
1639     {
1640       curve_pnt.Transform( loc );
1641       meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1642     }
1643   }
1644 //   GProp_GProps LProps;
1645 //   BRepGProp::LinearProperties(ed, LProps);
1646 //   double lg = (double)LProps.Mass();
1647
1648   meshDS->SetNodeOnEdge(node, edge, pa);
1649 }
1650
1651 //=============================================================================
1652 /*!
1653  *
1654  */
1655 //=============================================================================
1656
1657 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
1658 {
1659   return save;
1660 }
1661
1662 //=============================================================================
1663 /*!
1664  *
1665  */
1666 //=============================================================================
1667
1668 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
1669 {
1670   return load;
1671 }
1672
1673 //=============================================================================
1674 /*!
1675  *
1676  */
1677 //=============================================================================
1678
1679 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
1680 {
1681   return hyp.SaveTo( save );
1682 }
1683
1684 //=============================================================================
1685 /*!
1686  *
1687  */
1688 //=============================================================================
1689
1690 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
1691 {
1692   return hyp.LoadFrom( load );
1693 }
1694
1695 /* Curve definition function See cad_curv_t in file distene/cad.h for
1696  * more information.
1697  * NOTE : if when your CAD systems evaluates second
1698  * order derivatives it also computes first order derivatives and
1699  * function evaluation, you can optimize this example by making only
1700  * one CAD call and filling the necessary uv, dt, dtt arrays.
1701  */
1702 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
1703 {
1704   /* t is given. It contains the t (time) 1D parametric coordintaes
1705      of the point PreCAD/BLSurf is querying on the curve */
1706   
1707   /* user_data identifies the edge PreCAD/BLSurf is querying
1708    * (see cad_edge_new later in this example) */
1709   const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
1710
1711   if (uv){
1712    /* BLSurf is querying the function evaluation */
1713     gp_Pnt2d P;
1714     P=pargeo->Value(t);
1715     uv[0]=P.X(); uv[1]=P.Y();
1716   }
1717
1718   if(dt) {
1719    /* query for the first order derivatives */
1720     gp_Vec2d V1;
1721     V1=pargeo->DN(t,1);
1722     dt[0]=V1.X(); dt[1]=V1.Y();
1723   }
1724
1725   if(dtt){
1726     /* query for the second order derivatives */
1727     gp_Vec2d V2;
1728     V2=pargeo->DN(t,2);
1729     dtt[0]=V2.X(); dtt[1]=V2.Y();
1730   }
1731
1732   return STATUS_OK;
1733 }
1734
1735 /* Surface definition function.
1736  * See cad_surf_t in file distene/cad.h for more information.
1737  * NOTE : if when your CAD systems evaluates second order derivatives it also
1738  * computes first order derivatives and function evaluation, you can optimize 
1739  * this example by making only one CAD call and filling the necessary xyz, du, dv, etc.. 
1740  * arrays.
1741  */
1742 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1743                   real *duu, real *duv, real *dvv, void *user_data)
1744 {
1745   /* uv[2] is given. It contains the u,v coordinates of the point
1746    * PreCAD/BLSurf is querying on the surface */
1747   
1748   /* user_data identifies the face PreCAD/BLSurf is querying (see
1749    * cad_face_new later in this example)*/
1750   const Geom_Surface* geometry = (const Geom_Surface*) user_data;
1751
1752   if(xyz){
1753    gp_Pnt P;
1754    P=geometry->Value(uv[0],uv[1]);   // S.D0(U,V,P);
1755    xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
1756   }
1757
1758   if(du && dv){
1759     gp_Pnt P;
1760     gp_Vec D1U,D1V;
1761
1762     geometry->D1(uv[0],uv[1],P,D1U,D1V);
1763     du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
1764     dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
1765   }
1766
1767   if(duu && duv && dvv){
1768
1769     gp_Pnt P;
1770     gp_Vec D1U,D1V;
1771     gp_Vec D2U,D2V,D2UV;
1772
1773     geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
1774     duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
1775     duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
1776     dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
1777   }
1778
1779   return STATUS_OK;
1780 }
1781
1782
1783 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
1784 {
1785   if (face_id == 1) {
1786     if (my_u_min > uv[0]) {
1787       my_u_min = uv[0];
1788     }
1789     if (my_v_min > uv[1]) {
1790       my_v_min = uv[1];
1791     }
1792     if (my_u_max < uv[0]) {
1793       my_u_max = uv[0];
1794     }
1795     if (my_v_max < uv[1]) {
1796       my_v_max = uv[1];
1797     }
1798   }
1799   //MESSAGE("size_on_surface")
1800   if (FaceId2PythonSmp.count(face_id) != 0){
1801     //MESSAGE("A size map is used to calculate size on face : "<<face_id)
1802     PyObject * pyresult = NULL;
1803     PyObject* new_stderr = NULL;
1804     assert(Py_IsInitialized());
1805     PyGILState_STATE gstate;
1806     gstate = PyGILState_Ensure();
1807     pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
1808     double result;
1809     if ( pyresult == NULL){
1810       fflush(stderr);
1811       string err_description="";
1812       new_stderr = newPyStdOut(err_description);
1813       PySys_SetObject((char*)"stderr", new_stderr);
1814       PyErr_Print();
1815       PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
1816       Py_DECREF(new_stderr);
1817       MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
1818       result = *((double*)user_data);
1819       }
1820     else {
1821       result = PyFloat_AsDouble(pyresult);
1822       Py_DECREF(pyresult);
1823     }
1824    // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1825     *size = result;
1826     PyGILState_Release(gstate);
1827   }
1828   else if (FaceIndex2ClassAttractor.count(face_id) !=0 && !FaceIndex2ClassAttractor[face_id]->Empty()){
1829 //    MESSAGE("attractor used on face :"<<face_id)
1830     // MESSAGE("List of attractor is not empty")
1831     // MESSAGE("Attractor empty : "<< FaceIndex2ClassAttractor[face_id]->Empty())
1832     double result = FaceIndex2ClassAttractor[face_id]->GetSize(uv[0],uv[1]);
1833     *size = result;
1834    // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1835   }
1836   else {
1837     // MESSAGE("List of attractor is empty !!!")
1838     *size = *((double*)user_data);
1839   }
1840   return STATUS_OK;
1841 }
1842
1843 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
1844 {
1845   if (EdgeId2PythonSmp.count(edge_id) != 0){
1846     PyObject * pyresult = NULL;
1847     PyObject* new_stderr = NULL;
1848     assert(Py_IsInitialized());
1849     PyGILState_STATE gstate;
1850     gstate = PyGILState_Ensure();
1851     pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
1852     double result;
1853     if ( pyresult == NULL){
1854       fflush(stderr);
1855       string err_description="";
1856       new_stderr = newPyStdOut(err_description);
1857       PySys_SetObject((char*)"stderr", new_stderr);
1858       PyErr_Print();
1859       PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
1860       Py_DECREF(new_stderr);
1861       MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
1862       result = *((double*)user_data);
1863       }
1864     else {
1865       result = PyFloat_AsDouble(pyresult);
1866       Py_DECREF(pyresult);
1867     }
1868     *size = result;
1869     PyGILState_Release(gstate);
1870   }
1871   else {
1872     *size = *((double*)user_data);
1873   }
1874   return STATUS_OK;
1875 }
1876
1877 status_t size_on_vertex(integer point_id, real *size, void *user_data)
1878 {
1879   if (VertexId2PythonSmp.count(point_id) != 0){
1880     PyObject * pyresult = NULL;
1881     PyObject* new_stderr = NULL;
1882     assert(Py_IsInitialized());
1883     PyGILState_STATE gstate;
1884     gstate = PyGILState_Ensure();
1885     pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
1886     double result;
1887     if ( pyresult == NULL){
1888       fflush(stderr);
1889       string err_description="";
1890       new_stderr = newPyStdOut(err_description);
1891       PySys_SetObject((char*)"stderr", new_stderr);
1892       PyErr_Print();
1893       PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
1894       Py_DECREF(new_stderr);
1895       MESSAGE("Can't evaluate f()" << " error is " << err_description);
1896       result = *((double*)user_data);
1897       }
1898     else {
1899       result = PyFloat_AsDouble(pyresult);
1900       Py_DECREF(pyresult);
1901     }
1902     *size = result;
1903     PyGILState_Release(gstate);
1904   }
1905   else {
1906     *size = *((double*)user_data);
1907   }
1908  return STATUS_OK;
1909 }
1910
1911 /*
1912  * The following function will be called for PreCAD/BLSurf message
1913  * printing.  See context_set_message_callback (later in this
1914  * template) for how to set user_data.
1915  */
1916 status_t message_cb(message_t *msg, void *user_data)
1917 {
1918   integer errnumber = 0;
1919   char *desc;
1920   message_get_number(msg, &errnumber);
1921   message_get_description(msg, &desc);
1922   if ( errnumber < 0 ) {
1923     string * error = (string*)user_data;
1924 //   if ( !error->empty() )
1925 //     *error += "\n";
1926     // remove ^A from the tail
1927     int len = strlen( desc );
1928     while (len > 0 && desc[len-1] != '\n')
1929       len--;
1930     error->append( desc, len );
1931   }
1932   else {
1933       std::cout << desc << std::endl;
1934   }
1935   return STATUS_OK;
1936 }
1937
1938 /* This is the interrupt callback. PreCAD/BLSurf will call this
1939  * function regularily. See the file distene/interrupt.h
1940  */
1941 status_t interrupt_cb(integer *interrupt_status, void *user_data)
1942 {
1943   integer you_want_to_continue = 1;
1944 #ifdef WITH_SMESH_CANCEL_COMPUTE
1945   BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
1946   you_want_to_continue = !tmp->computeCanceled();
1947 #endif
1948
1949   if(you_want_to_continue)
1950   {
1951     *interrupt_status = INTERRUPT_CONTINUE;
1952     return STATUS_OK;
1953   }
1954   else /* you want to stop BLSurf */
1955   {
1956     *interrupt_status = INTERRUPT_STOP;
1957     return STATUS_ERROR;
1958   }
1959 }
1960
1961 //=============================================================================
1962 /*!
1963  *  
1964  */
1965 //=============================================================================
1966 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh&         aMesh,
1967                                    const TopoDS_Shape& aShape,
1968                                    MapShapeNbElems&    aResMap)
1969 {
1970   int    _physicalMesh  = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
1971   double _phySize       = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
1972   //int    _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
1973   //double _angleMeshS    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
1974   double _angleMeshC    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
1975   bool   _quadAllowed   = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
1976   if(_hypothesis) {
1977     _physicalMesh  = (int) _hypothesis->GetPhysicalMesh();
1978     _phySize       = _hypothesis->GetPhySize();
1979     //_geometricMesh = (int) hyp->GetGeometricMesh();
1980     //_angleMeshS    = hyp->GetAngleMeshS();
1981     _angleMeshC    = _hypothesis->GetAngleMeshC();
1982     _quadAllowed   = _hypothesis->GetQuadAllowed();
1983   }
1984
1985   bool IsQuadratic = false;
1986
1987   // ----------------
1988   // evaluate 1D 
1989   // ----------------
1990   TopTools_DataMapOfShapeInteger EdgesMap;
1991   double fullLen = 0.0;
1992   double fullNbSeg = 0;
1993   for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1994     TopoDS_Edge E = TopoDS::Edge( exp.Current() );
1995     if( EdgesMap.IsBound(E) )
1996       continue;
1997     SMESH_subMesh *sm = aMesh.GetSubMesh(E);
1998     double aLen = SMESH_Algo::EdgeLength(E);
1999     fullLen += aLen;
2000     int nb1d = 0;
2001     if(_physicalMesh==1) {
2002        nb1d = (int)( aLen/_phySize + 1 );
2003     }
2004     else {
2005       // use geometry
2006       double f,l;
2007       Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
2008       double fullAng = 0.0;
2009       double dp = (l-f)/200;
2010       gp_Pnt P1,P2,P3;
2011       C->D0(f,P1);
2012       C->D0(f+dp,P2);
2013       gp_Vec V1(P1,P2);
2014       for(int j=2; j<=200; j++) {
2015         C->D0(f+dp*j,P3);
2016         gp_Vec V2(P2,P3);
2017         fullAng += fabs(V1.Angle(V2));
2018         V1 = V2;
2019         P2 = P3;
2020       }
2021       nb1d = (int)( fullAng/_angleMeshC + 1 );
2022     }
2023     fullNbSeg += nb1d;
2024     std::vector<int> aVec(SMDSEntity_Last);
2025     for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2026     if( IsQuadratic > 0 ) {
2027       aVec[SMDSEntity_Node] = 2*nb1d - 1;
2028       aVec[SMDSEntity_Quad_Edge] = nb1d;
2029     }
2030     else {
2031       aVec[SMDSEntity_Node] = nb1d - 1;
2032       aVec[SMDSEntity_Edge] = nb1d;
2033     }
2034     aResMap.insert(std::make_pair(sm,aVec));
2035     EdgesMap.Bind(E,nb1d);
2036   }
2037   double ELen = fullLen/fullNbSeg;
2038   // ----------------
2039   // evaluate 2D 
2040   // ----------------
2041   // try to evaluate as in MEFISTO
2042   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2043     TopoDS_Face F = TopoDS::Face( exp.Current() );
2044     SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2045     GProp_GProps G;
2046     BRepGProp::SurfaceProperties(F,G);
2047     double anArea = G.Mass();
2048     int nb1d = 0;
2049     std::vector<int> nb1dVec;
2050     for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
2051       int nbSeg = EdgesMap.Find(exp1.Current());
2052       nb1d += nbSeg;
2053       nb1dVec.push_back( nbSeg );
2054     }
2055     int nbQuad = 0;
2056     int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
2057     int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
2058     if ( _quadAllowed )
2059     {
2060       if ( nb1dVec.size() == 4 ) // quadrangle geom face
2061       {
2062         int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
2063         nbQuad = n1 * n2;
2064         nbNodes = (n1 + 1) * (n2 + 1);
2065         nbTria = 0;
2066       }
2067       else
2068       {
2069         nbTria = nbQuad = nbTria / 3 + 1;
2070       }
2071     }
2072     std::vector<int> aVec(SMDSEntity_Last,0);
2073     if( IsQuadratic ) {
2074       int nb1d_in = (nbTria*3 - nb1d) / 2;
2075       aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
2076       aVec[SMDSEntity_Quad_Triangle] = nbTria;
2077       aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
2078     }
2079     else {
2080       aVec[SMDSEntity_Node] = nbNodes;
2081       aVec[SMDSEntity_Triangle] = nbTria;
2082       aVec[SMDSEntity_Quadrangle] = nbQuad;
2083     }
2084     aResMap.insert(std::make_pair(sm,aVec));
2085   }
2086
2087   // ----------------
2088   // evaluate 3D
2089   // ----------------
2090   GProp_GProps G;
2091   BRepGProp::VolumeProperties(aShape,G);
2092   double aVolume = G.Mass();
2093   double tetrVol = 0.1179*ELen*ELen*ELen;
2094   int nbVols  = int(aVolume/tetrVol);
2095   int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
2096   std::vector<int> aVec(SMDSEntity_Last);
2097   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2098   if( IsQuadratic ) {
2099     aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
2100     aVec[SMDSEntity_Quad_Tetra] = nbVols;
2101   }
2102   else {
2103     aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
2104     aVec[SMDSEntity_Tetra] = nbVols;
2105   }
2106   SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2107   aResMap.insert(std::make_pair(sm,aVec));
2108
2109   return true;
2110 }
2111
2112 //=============================================================================
2113 /*!
2114  *  Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
2115  *  Following line was added:
2116  *        myExtrem.Perform(P);
2117  */
2118 //=============================================================================
2119 void  BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
2120                     const TopoDS_Face& face, 
2121                     const gp_Pnt& P, 
2122                     const Standard_Real Tol)
2123 {
2124   //-- Voir BRepExtrema_ExtPF.cxx 
2125   BRepAdaptor_Surface Surf(face);
2126   Standard_Real U1, U2, V1, V2;
2127   BRepTools::UVBounds(face, U1, U2, V1, V2);
2128   Extrema_ExtPS myExtrem;
2129   myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
2130   myExtrem.Perform(P);
2131   //----------------------------------------------------------
2132   //-- On cherche le point le plus proche , PUIS 
2133   //-- On le classifie. 
2134   Standard_Integer nbv    = 0; // xpu
2135   Standard_Real MaxDist   =  RealLast();
2136   Standard_Integer indice = 0;
2137   if(myExtrem.IsDone()) {
2138     nbv = myExtrem.NbExt();
2139     for (Standard_Integer i = 1; i <= nbv; i++) {
2140       Standard_Real d = myExtrem.Value(i);
2141       d = Abs(d);
2142       if(d <= MaxDist) { 
2143     MaxDist = d;
2144     indice = i;
2145       }
2146     }
2147   }
2148   if(indice) { 
2149     gp_Pnt2d Puv;
2150     Standard_Real U1,U2;
2151     myExtrem.Point(indice).Parameter(U1, U2);
2152     Puv.SetCoord(U1, U2);
2153     fc->Perform(face, Puv, Tol);
2154   }
2155   else { 
2156     fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
2157     //-- le resultat est TopAbs_OUT;
2158   }
2159 }
2160