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