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