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