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