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