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