]> SALOME platform Git repositories - plugins/blsurfplugin.git/blob - src/BLSURFPlugin/BLSURFPlugin_BLSURF.cxx
Salome HOME
0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
[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     if (fmap.FindIndex(f) > 0)
927       continue;
928
929     fmap.Add(f);
930     iface++;
931     surfaces.push_back(BRep_Tool::Surface(f));
932     
933     /* create an object representing the face for blsurf */
934     /* where face_id is an integer identifying the face.
935      * surf_function is the function that defines the surface
936      * (For this face, it will be called by blsurf with your_face_object_ptr
937      * as last parameter.
938      */
939     cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
940     
941     /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
942     cad_face_set_tag(fce, iface);
943       
944     /* Set face orientation (optional if you want a well oriented output mesh)*/
945     if(f.Orientation() != TopAbs_FORWARD){
946       cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
947     } else {
948       cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
949     }
950     
951     if (HasSizeMapOnFace){
952       std::cout << "A size map is defined on a face" << std::endl;
953       // Classic size map
954       faceKey = FacesWithSizeMap.FindIndex(f);
955       
956       if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
957         theSizeMapStr = FaceId2SizeMap[faceKey];
958         // check if function ends with "return"
959         if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
960           continue;
961         // Expr To Python function, verification is performed at validation in GUI
962         PyObject * obj = NULL;
963         obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
964         Py_DECREF(obj);
965         PyObject * func = NULL;
966         func = PyObject_GetAttrString(main_mod, "f");
967         FaceId2PythonSmp[iface]=func;
968         FaceId2SizeMap.erase(faceKey);
969       }
970       
971       // Specific size map = Attractor
972       std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
973       int iatt=0;
974       for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
975         if (attractor_iter->first == faceKey) {
976           MESSAGE("Face indice: " << iface);
977           MESSAGE("Adding attractor");
978           
979           double xyzCoords[3]  = {attractor_iter->second[2],
980                                   attractor_iter->second[3],
981                                   attractor_iter->second[4]};
982           
983           MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
984           gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
985           BRepClass_FaceClassifier scl(f,P,1e-7);
986           // scl.Perform() is bugged. The function was rewritten
987 //          scl.Perform();
988           BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
989           TopAbs_State result = scl.State();
990           MESSAGE("Position of point on face: "<<result);
991           if ( result == TopAbs_OUT )
992               MESSAGE("Point is out of face: node is not created");
993           if ( result == TopAbs_UNKNOWN )
994               MESSAGE("Point position on face is unknown: node is not created");
995           if ( result == TopAbs_ON )
996               MESSAGE("Point is on border of face: node is not created");
997           if ( result == TopAbs_IN )
998           {
999             // Point is inside face and not on border
1000             MESSAGE("Point is in face: node is created");
1001             double uvCoords[2]   = {attractor_iter->second[0],attractor_iter->second[1]};
1002             iatt++;
1003             MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << iatt);
1004             cad_point_t* point_p = cad_point_new(fce, iatt, uvCoords);
1005             cad_point_set_tag(point_p, iatt);
1006           }
1007           FaceId2AttractorCoords.erase(faceKey);
1008         }
1009       }
1010       
1011       // Enforced Vertices
1012       faceKey = FacesWithEnforcedVertices.FindIndex(f);
1013       std::map<int,std::set<std::vector<double> > >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1014       if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1015         std::cout << "Some enforced vertices are defined" << std::endl;
1016         int ienf = 0;
1017         std::set<std::vector<double> > evl;
1018 //         std::vector<double> ev;
1019         MESSAGE("Face indice: " << iface);
1020         MESSAGE("Adding enforced vertices");
1021         evl = evmIt->second;
1022         MESSAGE("Number of vertices to add: "<< evl.size());
1023         std::set< std::vector<double> >::const_iterator evlIt = evl.begin();
1024         for (; evlIt != evl.end(); ++evlIt) {
1025 //           ev = *evlIt;
1026 //         for (int i=0; i<evl.size() ; i++) {
1027 //           ev = evl[i];
1028           
1029 //           double xyzCoords[3]  = {ev[2], ev[3], ev[4]};
1030           double xyzCoords[3]  = {evlIt->at(2), evlIt->at(3), evlIt->at(4)};
1031           MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1032           gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1033           BRepClass_FaceClassifier scl(f,P,1e-7);
1034           // scl.Perform() is bugged. The function was rewritten
1035 //          scl.Perform();
1036           BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1037           TopAbs_State result = scl.State();
1038           MESSAGE("Position of point on face: "<<result);
1039           if ( result == TopAbs_OUT )
1040               MESSAGE("Point is out of face: node is not created");
1041           if ( result == TopAbs_UNKNOWN )
1042               MESSAGE("Point position on face is unknown: node is not created");
1043           if ( result == TopAbs_ON )
1044               MESSAGE("Point is on border of face: node is not created");
1045           if ( result == TopAbs_IN )
1046           {
1047             // Point is inside face and not on border
1048             MESSAGE("Point is in face: node is created");
1049 //             double uvCoords[2]   = {ev[0],ev[1]};
1050             double uvCoords[2]   = {evlIt->at(0),evlIt->at(1)};
1051             ienf++;
1052             MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1053             cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1054             cad_point_set_tag(point_p, ienf);
1055           }
1056         }
1057         FaceId2EnforcedVertexCoords.erase(faceKey);
1058       }
1059       else
1060         std::cout << "No enforced vertex defined" << std::endl;
1061     }
1062     
1063     
1064     /****************************************************************************************
1065                                     EDGES
1066                    now create the edges associated to this face
1067     *****************************************************************************************/
1068     int edgeKey = -1;
1069     for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1070       TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1071       int ic = emap.FindIndex(e);
1072       if (ic <= 0)
1073         ic = emap.Add(e);
1074
1075       double tmin,tmax;
1076       curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1077       
1078       if (HasSizeMapOnEdge){
1079         edgeKey = EdgesWithSizeMap.FindIndex(e);
1080         if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1081           MESSAGE("Sizemap defined on edge with index " << edgeKey);
1082           theSizeMapStr = EdgeId2SizeMap[edgeKey];
1083           if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1084             continue;
1085           // Expr To Python function, verification is performed at validation in GUI
1086           PyObject * obj = NULL;
1087           obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1088           Py_DECREF(obj);
1089           PyObject * func = NULL;
1090           func = PyObject_GetAttrString(main_mod, "f");
1091           EdgeId2PythonSmp[ic]=func;
1092           EdgeId2SizeMap.erase(edgeKey);
1093         }
1094       }
1095       
1096       /* attach the edge to the current blsurf face */
1097       cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1098       
1099       /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1100       cad_edge_set_tag(edg, ic);
1101       
1102       /* by default, an edge does not necessalry appear in the resulting mesh,
1103      unless the following property is set :
1104       */
1105       cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1106       
1107       /* by default an edge is a boundary edge */
1108       if (e.Orientation() == TopAbs_INTERNAL)
1109         cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1110
1111       int npts = 0;
1112       int ip1, ip2, *ip;
1113       gp_Pnt2d e0 = curves.back()->Value(tmin);
1114       gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1115       Standard_Real d1=0,d2=0;
1116       
1117       
1118       /****************************************************************************************
1119                                       VERTICES
1120       *****************************************************************************************/
1121       int vertexKey = -1;
1122       for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1123         TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1124         ++npts;
1125         if (npts == 1){
1126           ip = &ip1;
1127           d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1128         } else {
1129           ip = &ip2;
1130           d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1131         }
1132         *ip = pmap.FindIndex(v);
1133         if(*ip <= 0)
1134           *ip = pmap.Add(v);
1135         
1136         //vertexKey = VerticesWithSizeMap.FindIndex(v);
1137         if (HasSizeMapOnVertex){
1138           vertexKey = VerticesWithSizeMap.FindIndex(v);
1139           if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1140             theSizeMapStr = VertexId2SizeMap[vertexKey];
1141             //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1142             if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1143               continue;
1144             // Expr To Python function, verification is performed at validation in GUI
1145             PyObject * obj = NULL;
1146             obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1147             Py_DECREF(obj);
1148             PyObject * func = NULL;
1149             func = PyObject_GetAttrString(main_mod, "f");
1150             VertexId2PythonSmp[*ip]=func;
1151             VertexId2SizeMap.erase(vertexKey);   // do not erase if using a vector
1152           }
1153         }
1154       }
1155       if (npts != 2) {
1156         // should not happen
1157         MESSAGE("An edge does not have 2 extremities.");
1158       } else {
1159         if (d1 < d2) {
1160           // This defines the curves extremity connectivity
1161           cad_edge_set_extremities(edg, ip1, ip2);
1162           /* set the tag (color) to the same value as the extremity id : */
1163           cad_edge_set_extremities_tag(edg, ip1, ip2);
1164         }
1165         else {
1166           cad_edge_set_extremities(edg, ip2, ip1);
1167           cad_edge_set_extremities_tag(edg, ip2, ip1);
1168         }
1169       }
1170     } // for edge
1171   } //for face
1172
1173
1174   PyGILState_Release(gstate);
1175
1176   blsurf_data_set_cad(bls, c);
1177
1178   std::cout << std::endl;
1179   std::cout << "Beginning of Surface Mesh generation" << std::endl;
1180   std::cout << std::endl;
1181
1182   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1183 #ifndef WNT
1184   feclearexcept( FE_ALL_EXCEPT );
1185   int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1186 #endif
1187
1188   try {
1189     OCC_CATCH_SIGNALS;
1190
1191     status = blsurf_compute_mesh(bls);
1192
1193   }
1194   catch ( std::exception& exc ) {
1195     _comment += exc.what();
1196   }
1197   catch (Standard_Failure& ex) {
1198     _comment += ex.DynamicType()->Name();
1199     if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1200       _comment += ": ";
1201       _comment += ex.GetMessageString();
1202     }
1203   }
1204   catch (...) {
1205     if ( _comment.empty() )
1206       _comment = "Exception in blsurf_compute_mesh()";
1207   }
1208   if ( status != STATUS_OK) {
1209     // There was an error while meshing
1210     blsurf_session_delete(bls);
1211     cad_delete(c);
1212     context_delete(ctx);
1213
1214     return error(_comment);
1215   }
1216
1217   std::cout << std::endl;
1218   std::cout << "End of Surface Mesh generation" << std::endl;
1219   std::cout << std::endl;
1220
1221   mesh_t *msh = NULL;
1222   blsurf_data_get_mesh(bls, &msh);
1223   if(!msh){
1224     blsurf_session_delete(bls);
1225     cad_delete(c);
1226     context_delete(ctx);
1227
1228     return error(_comment);
1229     //return false;
1230   }
1231
1232   /* retrieve mesh data (see distene/mesh.h) */
1233   integer nv, ne, nt, nq, vtx[4], tag;
1234   real xyz[3];
1235
1236   mesh_get_vertex_count(msh, &nv);
1237   mesh_get_edge_count(msh, &ne);
1238   mesh_get_triangle_count(msh, &nt);
1239   mesh_get_quadrangle_count(msh, &nq);
1240
1241
1242   SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1243   SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1244   bool* tags = new bool[nv+1];
1245
1246   /* enumerated vertices */
1247   for(int iv=1;iv<=nv;iv++) {
1248     mesh_get_vertex_coordinates(msh, iv, xyz);
1249     mesh_get_vertex_tag(msh, iv, &tag);
1250     // Issue 0020656. Use vertex coordinates
1251     if ( tag > 0 && tag <= pmap.Extent() ) {
1252       TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
1253       double tol = BRep_Tool::Tolerance( v );
1254       gp_Pnt p = BRep_Tool::Pnt( v );
1255       if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
1256         xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
1257       else
1258         tag = 0; // enforced or attracted vertex
1259     }
1260     nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1261     // internal point are tagged to zero
1262     if(tag > 0 && tag <= pmap.Extent() ){
1263       meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1264       tags[iv] = false;
1265     } else {
1266       tags[iv] = true;
1267     }
1268   }
1269
1270   /* enumerate edges */
1271   for(int it=1;it<=ne;it++) {
1272     mesh_get_edge_vertices(msh, it, vtx);
1273     SMDS_MeshEdge* edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1274     mesh_get_edge_tag(msh, it, &tag);
1275
1276     if (tags[vtx[0]]) {
1277       Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1278       tags[vtx[0]] = false;
1279     };
1280     if (tags[vtx[1]]) {
1281       Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1282       tags[vtx[1]] = false;
1283     };
1284     meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1285
1286   }
1287
1288   /* enumerate triangles */
1289   for(int it=1;it<=nt;it++) {
1290     mesh_get_triangle_vertices(msh, it, vtx);
1291     SMDS_MeshFace* tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1292     mesh_get_triangle_tag(msh, it, &tag);
1293     meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1294     if (tags[vtx[0]]) {
1295       meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1296       tags[vtx[0]] = false;
1297     };
1298     if (tags[vtx[1]]) {
1299       meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1300       tags[vtx[1]] = false;
1301     };
1302     if (tags[vtx[2]]) {
1303       meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1304       tags[vtx[2]] = false;
1305     };
1306   }
1307
1308   /* enumerate quadrangles */
1309   for(int it=1;it<=nq;it++) {
1310     mesh_get_quadrangle_vertices(msh, it, vtx);
1311     SMDS_MeshFace* quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1312     mesh_get_quadrangle_tag(msh, it, &tag);
1313     meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1314     if (tags[vtx[0]]) {
1315       meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1316       tags[vtx[0]] = false;
1317     };
1318     if (tags[vtx[1]]) {
1319       meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1320       tags[vtx[1]] = false;
1321     };
1322     if (tags[vtx[2]]) {
1323       meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1324       tags[vtx[2]] = false;
1325     };
1326     if (tags[vtx[3]]) {
1327       meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1328       tags[vtx[3]] = false;
1329     };
1330   }
1331
1332   delete nodes;
1333
1334   /* release the mesh object */
1335   blsurf_data_regain_mesh(bls, msh);
1336
1337   /* clean up everything */
1338   blsurf_session_delete(bls);
1339   cad_delete(c);
1340
1341   context_delete(ctx);
1342
1343   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1344 #ifndef WNT
1345   if ( oldFEFlags > 0 )
1346     feenableexcept( oldFEFlags );
1347   feclearexcept( FE_ALL_EXCEPT );
1348 #endif
1349
1350   /*  
1351   std::cout << "FacesWithSizeMap" << std::endl;
1352   FacesWithSizeMap.Statistics(std::cout);
1353   std::cout << "EdgesWithSizeMap" << std::endl;
1354   EdgesWithSizeMap.Statistics(std::cout);
1355   std::cout << "VerticesWithSizeMap" << std::endl;
1356   VerticesWithSizeMap.Statistics(std::cout);
1357   std::cout << "FacesWithEnforcedVertices" << std::endl;
1358   FacesWithEnforcedVertices.Statistics(std::cout);
1359   */
1360   
1361   return true;
1362 }
1363
1364 //=============================================================================
1365 /*!
1366  *  SetNodeOnEdge
1367  */
1368 //=============================================================================
1369
1370 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1371   const TopoDS_Edge edge = TopoDS::Edge(ed);
1372
1373   gp_Pnt pnt(node->X(), node->Y(), node->Z());
1374
1375   Standard_Real p0 = 0.0;
1376   Standard_Real p1 = 1.0;
1377   TopLoc_Location loc;
1378   Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1379
1380   if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1381   GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1382
1383   double pa = 0.;
1384   if ( proj.NbPoints() > 0 )
1385   {
1386     pa = (double)proj.LowerDistanceParameter();
1387     // Issue 0020656. Move node if it is too far from edge
1388     gp_Pnt curve_pnt = curve->Value( pa );
1389     double dist2 = pnt.SquareDistance( curve_pnt );
1390     double tol = BRep_Tool::Tolerance( edge );
1391     if ( 1e-12 < dist2 && dist2 <= 2*tol*tol ) // large enough and within tolerance
1392     {
1393       curve_pnt.Transform( loc );
1394       meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1395     }
1396   }
1397 //   GProp_GProps LProps;
1398 //   BRepGProp::LinearProperties(ed, LProps);
1399 //   double lg = (double)LProps.Mass();
1400
1401   meshDS->SetNodeOnEdge(node, edge, pa);
1402 }
1403
1404 //=============================================================================
1405 /*!
1406  *
1407  */
1408 //=============================================================================
1409
1410 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
1411 {
1412   return save;
1413 }
1414
1415 //=============================================================================
1416 /*!
1417  *
1418  */
1419 //=============================================================================
1420
1421 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
1422 {
1423   return load;
1424 }
1425
1426 //=============================================================================
1427 /*!
1428  *
1429  */
1430 //=============================================================================
1431
1432 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
1433 {
1434   return hyp.SaveTo( save );
1435 }
1436
1437 //=============================================================================
1438 /*!
1439  *
1440  */
1441 //=============================================================================
1442
1443 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
1444 {
1445   return hyp.LoadFrom( load );
1446 }
1447
1448 /* Curve definition function See cad_curv_t in file distene/cad.h for
1449  * more information.
1450  * NOTE : if when your CAD systems evaluates second
1451  * order derivatives it also computes first order derivatives and
1452  * function evaluation, you can optimize this example by making only
1453  * one CAD call and filling the necessary uv, dt, dtt arrays.
1454  */
1455 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
1456 {
1457   /* t is given. It contains the t (time) 1D parametric coordintaes
1458      of the point PreCAD/BLSurf is querying on the curve */
1459   
1460   /* user_data identifies the edge PreCAD/BLSurf is querying
1461    * (see cad_edge_new later in this example) */
1462   const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
1463
1464   if (uv){
1465    /* BLSurf is querying the function evaluation */
1466     gp_Pnt2d P;
1467     P=pargeo->Value(t);
1468     uv[0]=P.X(); uv[1]=P.Y();
1469   }
1470
1471   if(dt) {
1472    /* query for the first order derivatives */
1473     gp_Vec2d V1;
1474     V1=pargeo->DN(t,1);
1475     dt[0]=V1.X(); dt[1]=V1.Y();
1476   }
1477
1478   if(dtt){
1479     /* query for the second order derivatives */
1480     gp_Vec2d V2;
1481     V2=pargeo->DN(t,2);
1482     dtt[0]=V2.X(); dtt[1]=V2.Y();
1483   }
1484
1485   return STATUS_OK;
1486 }
1487
1488 /* Surface definition function.
1489  * See cad_surf_t in file distene/cad.h for more information.
1490  * NOTE : if when your CAD systems evaluates second order derivatives it also
1491  * computes first order derivatives and function evaluation, you can optimize 
1492  * this example by making only one CAD call and filling the necessary xyz, du, dv, etc.. 
1493  * arrays.
1494  */
1495 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1496                   real *duu, real *duv, real *dvv, void *user_data)
1497 {
1498   /* uv[2] is given. It contains the u,v coordinates of the point
1499    * PreCAD/BLSurf is querying on the surface */
1500   
1501   /* user_data identifies the face PreCAD/BLSurf is querying (see
1502    * cad_face_new later in this example)*/
1503   const Geom_Surface* geometry = (const Geom_Surface*) user_data;
1504
1505   if(xyz){
1506    gp_Pnt P;
1507    P=geometry->Value(uv[0],uv[1]);   // S.D0(U,V,P);
1508    xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
1509   }
1510
1511   if(du && dv){
1512     gp_Pnt P;
1513     gp_Vec D1U,D1V;
1514
1515     geometry->D1(uv[0],uv[1],P,D1U,D1V);
1516     du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
1517     dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
1518   }
1519
1520   if(duu && duv && dvv){
1521
1522     gp_Pnt P;
1523     gp_Vec D1U,D1V;
1524     gp_Vec D2U,D2V,D2UV;
1525
1526     geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
1527     duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
1528     duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
1529     dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
1530   }
1531
1532   return STATUS_OK;
1533 }
1534
1535
1536 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
1537 {
1538   if (face_id == 1) {
1539     if (my_u_min > uv[0]) {
1540       my_u_min = uv[0];
1541     }
1542     if (my_v_min > uv[1]) {
1543       my_v_min = uv[1];
1544     }
1545     if (my_u_max < uv[0]) {
1546       my_u_max = uv[0];
1547     }
1548     if (my_v_max < uv[1]) {
1549       my_v_max = uv[1];
1550     }
1551   }
1552
1553   if (FaceId2PythonSmp.count(face_id) != 0){
1554     PyObject * pyresult = NULL;
1555     PyObject* new_stderr = NULL;
1556     assert(Py_IsInitialized());
1557     PyGILState_STATE gstate;
1558     gstate = PyGILState_Ensure();
1559     pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],"(f,f)",uv[0],uv[1]);
1560     double result;
1561     if ( pyresult == NULL){
1562       fflush(stderr);
1563       string err_description="";
1564       new_stderr = newPyStdOut(err_description);
1565       PySys_SetObject("stderr", new_stderr);
1566       PyErr_Print();
1567       PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1568       Py_DECREF(new_stderr);
1569       MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
1570       result = *((double*)user_data);
1571       }
1572     else {
1573       result = PyFloat_AsDouble(pyresult);
1574       Py_DECREF(pyresult);
1575     }
1576     *size = result;
1577     //MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1578     PyGILState_Release(gstate);
1579   }
1580   else {
1581     *size = *((double*)user_data);
1582   }
1583   return STATUS_OK;
1584 }
1585
1586 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
1587 {
1588   if (EdgeId2PythonSmp.count(edge_id) != 0){
1589     PyObject * pyresult = NULL;
1590     PyObject* new_stderr = NULL;
1591     assert(Py_IsInitialized());
1592     PyGILState_STATE gstate;
1593     gstate = PyGILState_Ensure();
1594     pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],"(f)",t);
1595     double result;
1596     if ( pyresult == NULL){
1597       fflush(stderr);
1598       string err_description="";
1599       new_stderr = newPyStdOut(err_description);
1600       PySys_SetObject("stderr", new_stderr);
1601       PyErr_Print();
1602       PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1603       Py_DECREF(new_stderr);
1604       MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
1605       result = *((double*)user_data);
1606       }
1607     else {
1608       result = PyFloat_AsDouble(pyresult);
1609       Py_DECREF(pyresult);
1610     }
1611     *size = result;
1612     PyGILState_Release(gstate);
1613   }
1614   else {
1615     *size = *((double*)user_data);
1616   }
1617   return STATUS_OK;
1618 }
1619
1620 status_t size_on_vertex(integer point_id, real *size, void *user_data)
1621 {
1622   if (VertexId2PythonSmp.count(point_id) != 0){
1623     PyObject * pyresult = NULL;
1624     PyObject* new_stderr = NULL;
1625     assert(Py_IsInitialized());
1626     PyGILState_STATE gstate;
1627     gstate = PyGILState_Ensure();
1628     pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],"");
1629     double result;
1630     if ( pyresult == NULL){
1631       fflush(stderr);
1632       string err_description="";
1633       new_stderr = newPyStdOut(err_description);
1634       PySys_SetObject("stderr", new_stderr);
1635       PyErr_Print();
1636       PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1637       Py_DECREF(new_stderr);
1638       MESSAGE("Can't evaluate f()" << " error is " << err_description);
1639       result = *((double*)user_data);
1640       }
1641     else {
1642       result = PyFloat_AsDouble(pyresult);
1643       Py_DECREF(pyresult);
1644     }
1645     *size = result;
1646     PyGILState_Release(gstate);
1647   }
1648   else {
1649     *size = *((double*)user_data);
1650   }
1651  return STATUS_OK;
1652 }
1653
1654 /*
1655  * The following function will be called for PreCAD/BLSurf message
1656  * printing.  See context_set_message_callback (later in this
1657  * template) for how to set user_data.
1658  */
1659 status_t message_cb(message_t *msg, void *user_data)
1660 {
1661   integer errnumber = 0;
1662   char *desc;
1663   message_get_number(msg, &errnumber);
1664   message_get_description(msg, &desc);
1665   if ( errnumber < 0 ) {
1666     string * error = (string*)user_data;
1667 //   if ( !error->empty() )
1668 //     *error += "\n";
1669     // remove ^A from the tail
1670     int len = strlen( desc );
1671     while (len > 0 && desc[len-1] != '\n')
1672       len--;
1673     error->append( desc, len );
1674   }
1675   else {
1676       std::cout << desc << std::endl;
1677   }
1678   return STATUS_OK;
1679 }
1680
1681 /* This is the interrupt callback. PreCAD/BLSurf will call this
1682  * function regularily. See the file distene/interrupt.h
1683  */
1684 status_t interrupt_cb(integer *interrupt_status, void *user_data)
1685 {
1686   integer you_want_to_continue = 1;
1687
1688   if(you_want_to_continue)
1689     *interrupt_status = INTERRUPT_CONTINUE;
1690   else /* you want to stop BLSurf */
1691     *interrupt_status = INTERRUPT_STOP;
1692   
1693   return STATUS_OK;
1694 }
1695
1696 //=============================================================================
1697 /*!
1698  *  
1699  */
1700 //=============================================================================
1701 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh&         aMesh,
1702                                    const TopoDS_Shape& aShape,
1703                                    MapShapeNbElems&    aResMap)
1704 {
1705   int    _physicalMesh  = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
1706   double _phySize       = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
1707   //int    _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
1708   //double _angleMeshS    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
1709   double _angleMeshC    = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
1710   bool   _quadAllowed   = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
1711   if(_hypothesis) {
1712     _physicalMesh  = (int) _hypothesis->GetPhysicalMesh();
1713     _phySize       = _hypothesis->GetPhySize();
1714     //_geometricMesh = (int) hyp->GetGeometricMesh();
1715     //_angleMeshS    = hyp->GetAngleMeshS();
1716     _angleMeshC    = _hypothesis->GetAngleMeshC();
1717     _quadAllowed   = _hypothesis->GetQuadAllowed();
1718   }
1719
1720   bool IsQuadratic = false;
1721
1722   // ----------------
1723   // evaluate 1D 
1724   // ----------------
1725   TopTools_DataMapOfShapeInteger EdgesMap;
1726   double fullLen = 0.0;
1727   double fullNbSeg = 0;
1728   for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1729     TopoDS_Edge E = TopoDS::Edge( exp.Current() );
1730     if( EdgesMap.IsBound(E) )
1731       continue;
1732     SMESH_subMesh *sm = aMesh.GetSubMesh(E);
1733     double aLen = SMESH_Algo::EdgeLength(E);
1734     fullLen += aLen;
1735     int nb1d = 0;
1736     if(_physicalMesh==1) {
1737        nb1d = (int)( aLen/_phySize + 1 );
1738     }
1739     else {
1740       // use geometry
1741       double f,l;
1742       Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
1743       double fullAng = 0.0;
1744       double dp = (l-f)/200;
1745       gp_Pnt P1,P2,P3;
1746       C->D0(f,P1);
1747       C->D0(f+dp,P2);
1748       gp_Vec V1(P1,P2);
1749       for(int j=2; j<=200; j++) {
1750         C->D0(f+dp*j,P3);
1751         gp_Vec V2(P2,P3);
1752         fullAng += fabs(V1.Angle(V2));
1753         V1 = V2;
1754         P2 = P3;
1755       }
1756       nb1d = (int)( fullAng/_angleMeshC + 1 );
1757     }
1758     fullNbSeg += nb1d;
1759     std::vector<int> aVec(SMDSEntity_Last);
1760     for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1761     if( IsQuadratic > 0 ) {
1762       aVec[SMDSEntity_Node] = 2*nb1d - 1;
1763       aVec[SMDSEntity_Quad_Edge] = nb1d;
1764     }
1765     else {
1766       aVec[SMDSEntity_Node] = nb1d - 1;
1767       aVec[SMDSEntity_Edge] = nb1d;
1768     }
1769     aResMap.insert(std::make_pair(sm,aVec));
1770     EdgesMap.Bind(E,nb1d);
1771   }
1772   double ELen = fullLen/fullNbSeg;
1773   // ----------------
1774   // evaluate 2D 
1775   // ----------------
1776   // try to evaluate as in MEFISTO
1777   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1778     TopoDS_Face F = TopoDS::Face( exp.Current() );
1779     SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1780     GProp_GProps G;
1781     BRepGProp::SurfaceProperties(F,G);
1782     double anArea = G.Mass();
1783     int nb1d = 0;
1784     std::vector<int> nb1dVec;
1785     for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
1786       int nbSeg = EdgesMap.Find(exp1.Current());
1787       nb1d += nbSeg;
1788       nb1dVec.push_back( nbSeg );
1789     }
1790     int nbQuad = 0;
1791     int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
1792     int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
1793     if ( _quadAllowed )
1794     {
1795       if ( nb1dVec.size() == 4 ) // quadrangle geom face
1796       {
1797         int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
1798         nbQuad = n1 * n2;
1799         nbNodes = (n1 + 1) * (n2 + 1);
1800         nbTria = 0;
1801       }
1802       else
1803       {
1804         nbTria = nbQuad = nbTria / 3 + 1;
1805       }
1806     }
1807     std::vector<int> aVec(SMDSEntity_Last,0);
1808     if( IsQuadratic ) {
1809       int nb1d_in = (nbTria*3 - nb1d) / 2;
1810       aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
1811       aVec[SMDSEntity_Quad_Triangle] = nbTria;
1812       aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
1813     }
1814     else {
1815       aVec[SMDSEntity_Node] = nbNodes;
1816       aVec[SMDSEntity_Triangle] = nbTria;
1817       aVec[SMDSEntity_Quadrangle] = nbQuad;
1818     }
1819     aResMap.insert(std::make_pair(sm,aVec));
1820   }
1821
1822   // ----------------
1823   // evaluate 3D
1824   // ----------------
1825   GProp_GProps G;
1826   BRepGProp::VolumeProperties(aShape,G);
1827   double aVolume = G.Mass();
1828   double tetrVol = 0.1179*ELen*ELen*ELen;
1829   int nbVols  = int(aVolume/tetrVol);
1830   int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
1831   std::vector<int> aVec(SMDSEntity_Last);
1832   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1833   if( IsQuadratic ) {
1834     aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
1835     aVec[SMDSEntity_Quad_Tetra] = nbVols;
1836   }
1837   else {
1838     aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
1839     aVec[SMDSEntity_Tetra] = nbVols;
1840   }
1841   SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1842   aResMap.insert(std::make_pair(sm,aVec));
1843
1844   return true;
1845 }
1846
1847 //=============================================================================
1848 /*!
1849  *  Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
1850  *  Following line was added:
1851  *        myExtrem.Perform(P);
1852  */
1853 //=============================================================================
1854 void  BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
1855                     const TopoDS_Face& face, 
1856                     const gp_Pnt& P, 
1857                     const Standard_Real Tol)
1858 {
1859   //-- Voir BRepExtrema_ExtPF.cxx 
1860   BRepAdaptor_Surface Surf(face);
1861   Standard_Real U1, U2, V1, V2;
1862   BRepTools::UVBounds(face, U1, U2, V1, V2);
1863   Extrema_ExtPS myExtrem;
1864   myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
1865   myExtrem.Perform(P);
1866   //----------------------------------------------------------
1867   //-- On cherche le point le plus proche , PUIS 
1868   //-- On le classifie. 
1869   Standard_Integer nbv    = 0; // xpu
1870   Standard_Real MaxDist   =  RealLast();
1871   Standard_Integer indice = 0;
1872   if(myExtrem.IsDone()) {
1873     nbv = myExtrem.NbExt();
1874     for (Standard_Integer i = 1; i <= nbv; i++) {
1875       Standard_Real d = myExtrem.Value(i);
1876       d = Abs(d);
1877       if(d <= MaxDist) { 
1878     MaxDist = d;
1879     indice = i;
1880       }
1881     }
1882   }
1883   if(indice) { 
1884     gp_Pnt2d Puv;
1885     Standard_Real U1,U2;
1886     myExtrem.Point(indice).Parameter(U1, U2);
1887     Puv.SetCoord(U1, U2);
1888     fc->Perform(face, Puv, Tol);
1889   }
1890   else { 
1891     fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
1892     //-- le resultat est TopAbs_OUT;
1893   }
1894 }
1895