1 // Copyright (C) 2007-2011 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
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)
27 #include "BLSURFPlugin_BLSURF.hxx"
28 #include "BLSURFPlugin_Hypothesis.hxx"
29 #include "BLSURFPlugin_Attractor.hxx"
32 #include <distene/api.h>
33 #include <distene/blsurf.h>
36 #include <structmember.h>
39 #include <Basics_Utils.hxx>
41 #include <SMESH_Gen.hxx>
42 #include <SMESH_Mesh.hxx>
43 #include <SMESH_ControlsDef.hxx>
44 #include <SMDSAbs_ElementType.hxx>
45 #include "SMESHDS_Group.hxx"
46 #include "SMESH_Group.hxx"
48 #include <utilities.h>
56 // OPENCASCADE includes
57 #include <BRep_Tool.hxx>
59 #include <TopExp_Explorer.hxx>
61 #include <NCollection_Map.hxx>
63 #include <Geom_Surface.hxx>
64 #include <Handle_Geom_Surface.hxx>
65 #include <Geom2d_Curve.hxx>
66 #include <Handle_Geom2d_Curve.hxx>
67 #include <Geom_Curve.hxx>
68 #include <Handle_Geom_Curve.hxx>
69 #include <Handle_AIS_InteractiveObject.hxx>
70 #include <TopoDS_Vertex.hxx>
71 #include <TopoDS_Edge.hxx>
72 #include <TopoDS_Wire.hxx>
73 #include <TopoDS_Face.hxx>
75 #include <gp_Pnt2d.hxx>
77 #include <TopTools_IndexedMapOfShape.hxx>
78 #include <TopoDS_Shape.hxx>
79 #include <BRep_Builder.hxx>
80 #include <BRepBuilderAPI_MakeVertex.hxx>
81 #include <BRepTools.hxx>
83 #include <TopTools_DataMapOfShapeInteger.hxx>
84 #include <GProp_GProps.hxx>
85 #include <BRepGProp.hxx>
91 #include <Standard_ErrorHandler.hxx>
92 #include <GeomAPI_ProjectPointOnCurve.hxx>
93 #include <GeomAPI_ProjectPointOnSurf.hxx>
96 // #include <BRepClass_FaceClassifier.hxx>
97 #include <TopTools_MapOfShape.hxx>
99 /* ==================================
100 * =========== PYTHON ==============
101 * ==================================*/
110 PyStdOut_dealloc(PyStdOut *self)
116 PyStdOut_write(PyStdOut *self, PyObject *args)
120 if (!PyArg_ParseTuple(args, "t#:write",&c, &l))
124 *(self->out)=*(self->out)+c;
130 static PyMethodDef PyStdOut_methods[] = {
131 {"write", (PyCFunction)PyStdOut_write, METH_VARARGS,
132 PyDoc_STR("write(string) -> None")},
133 {NULL, NULL} /* sentinel */
136 static PyMemberDef PyStdOut_memberlist[] = {
137 {(char*)"softspace", T_INT, offsetof(PyStdOut, softspace), 0,
138 (char*)"flag indicating that a space needs to be printed; used by print"},
139 {NULL} /* Sentinel */
142 static PyTypeObject PyStdOut_Type = {
143 /* The ob_type field must be initialized in the module init function
144 * to be portable to Windows without using C++. */
145 PyObject_HEAD_INIT(NULL)
148 sizeof(PyStdOut), /*tp_basicsize*/
151 (destructor)PyStdOut_dealloc, /*tp_dealloc*/
158 0, /*tp_as_sequence*/
163 PyObject_GenericGetAttr, /*tp_getattro*/
164 /* softspace is writable: we must supply tp_setattro */
165 PyObject_GenericSetAttr, /* tp_setattro */
167 Py_TPFLAGS_DEFAULT, /*tp_flags*/
171 0, /*tp_richcompare*/
172 0, /*tp_weaklistoffset*/
175 PyStdOut_methods, /*tp_methods*/
176 PyStdOut_memberlist, /*tp_members*/
190 PyObject * newPyStdOut( std::string& out )
193 self = PyObject_New(PyStdOut, &PyStdOut_Type);
198 return (PyObject*)self;
202 ////////////////////////END PYTHON///////////////////////////
204 //////////////////MY MAPS////////////////////////////////////////
205 TopTools_IndexedMapOfShape FacesWithSizeMap;
206 std::map<int,string> FaceId2SizeMap;
207 TopTools_IndexedMapOfShape EdgesWithSizeMap;
208 std::map<int,string> EdgeId2SizeMap;
209 TopTools_IndexedMapOfShape VerticesWithSizeMap;
210 std::map<int,string> VertexId2SizeMap;
212 std::map<int,PyObject*> FaceId2PythonSmp;
213 std::map<int,PyObject*> EdgeId2PythonSmp;
214 std::map<int,PyObject*> VertexId2PythonSmp;
216 std::map<int,std::vector<double> > FaceId2AttractorCoords;
217 std::map<int,BLSURFPlugin_Attractor*> FaceId2ClassAttractor;
218 std::map<int,BLSURFPlugin_Attractor*> FaceIndex2ClassAttractor;
220 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
221 std::map< int, BLSURFPlugin_Hypothesis::TEnfVertexCoordsList > FaceId2EnforcedVertexCoords;
222 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexCoords > EnfVertexCoords2ProjVertex;
223 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList > EnfVertexCoords2EnfVertexList;
225 bool HasSizeMapOnFace=false;
226 bool HasSizeMapOnEdge=false;
227 bool HasSizeMapOnVertex=false;
228 //bool HasAttractorOnFace=false;
230 //=============================================================================
234 //=============================================================================
236 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId, int studyId,
238 : SMESH_2D_Algo(hypId, studyId, gen)
240 MESSAGE("BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF");
243 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
244 _compatibleHypothesis.push_back("BLSURF_Parameters");
245 _requireDescretBoundary = false;
246 _onlyUnaryInput = false;
249 smeshGen_i = SMESH_Gen_i::GetSMESHGen();
250 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
251 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
253 MESSAGE("studyid = " << _studyId);
256 myStudy = aStudyMgr->GetStudyByID(_studyId);
258 MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
260 /* Initialize the Python interpreter */
261 assert(Py_IsInitialized());
262 PyGILState_STATE gstate;
263 gstate = PyGILState_Ensure();
266 main_mod = PyImport_AddModule("__main__");
269 main_dict = PyModule_GetDict(main_mod);
271 PyRun_SimpleString("from math import *");
272 PyGILState_Release(gstate);
274 FacesWithSizeMap.Clear();
275 FaceId2SizeMap.clear();
276 EdgesWithSizeMap.Clear();
277 EdgeId2SizeMap.clear();
278 VerticesWithSizeMap.Clear();
279 VertexId2SizeMap.clear();
280 FaceId2PythonSmp.clear();
281 EdgeId2PythonSmp.clear();
282 VertexId2PythonSmp.clear();
283 FaceId2AttractorCoords.clear();
284 FaceId2ClassAttractor.clear();
285 FaceIndex2ClassAttractor.clear();
286 FacesWithEnforcedVertices.Clear();
287 FaceId2EnforcedVertexCoords.clear();
288 EnfVertexCoords2ProjVertex.clear();
289 EnfVertexCoords2EnfVertexList.clear();
291 #ifdef WITH_SMESH_CANCEL_COMPUTE
292 _compute_canceled = false;
296 //=============================================================================
300 //=============================================================================
302 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
304 MESSAGE("BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF");
308 //=============================================================================
312 //=============================================================================
314 bool BLSURFPlugin_BLSURF::CheckHypothesis
316 const TopoDS_Shape& aShape,
317 SMESH_Hypothesis::Hypothesis_Status& aStatus)
321 list<const SMESHDS_Hypothesis*>::const_iterator itl;
322 const SMESHDS_Hypothesis* theHyp;
324 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
325 int nbHyp = hyps.size();
328 aStatus = SMESH_Hypothesis::HYP_OK;
329 return true; // can work with no hypothesis
333 theHyp = (*itl); // use only the first hypothesis
335 string hypName = theHyp->GetName();
337 if (hypName == "BLSURF_Parameters")
339 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
341 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
342 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
343 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
344 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
346 aStatus = SMESH_Hypothesis::HYP_OK;
349 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
351 return aStatus == SMESH_Hypothesis::HYP_OK;
354 //=============================================================================
356 * Pass parameters to BLSURF
358 //=============================================================================
360 inline std::string to_string(double d)
362 std::ostringstream o;
367 inline std::string to_string(int i)
369 std::ostringstream o;
374 double _smp_phy_size;
375 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
376 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
377 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
379 double my_u_min=1e6,my_v_min=1e6,my_u_max=-1e6,my_v_max=-1e6;
385 /////////////////////////////////////////////////////////
386 projectionPoint getProjectionPoint(const TopoDS_Face& face, const gp_Pnt& point)
388 projectionPoint myPoint;
389 Handle(Geom_Surface) surface = BRep_Tool::Surface(face);
390 GeomAPI_ProjectPointOnSurf projector( point, surface );
391 if ( !projector.IsDone() || projector.NbPoints()==0 )
392 throw "getProjectionPoint: Can't project";
394 Quantity_Parameter u,v;
395 projector.LowerDistanceParameters(u,v);
396 myPoint.uv = gp_XY(u,v);
397 gp_Pnt aPnt = projector.NearestPoint();
398 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
402 /////////////////////////////////////////////////////////
404 /////////////////////////////////////////////////////////
405 double getT(const TopoDS_Edge& edge, const gp_Pnt& point)
408 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, f,l);
409 GeomAPI_ProjectPointOnCurve projector( point, curve);
410 if ( projector.NbPoints() == 0 )
412 return projector.LowerDistanceParameter();
415 /////////////////////////////////////////////////////////
416 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
418 MESSAGE("BLSURFPlugin_BLSURF::entryToShape "<<entry );
419 GEOM::GEOM_Object_var aGeomObj;
420 TopoDS_Shape S = TopoDS_Shape();
421 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
422 SALOMEDS::GenericAttribute_var anAttr;
424 if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR")) {
425 SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
426 CORBA::String_var aVal = anIOR->Value();
427 CORBA::Object_var obj = myStudy->ConvertIORToObject(aVal);
428 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
430 if ( !aGeomObj->_is_nil() )
431 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
435 void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
437 BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
442 // Get the (u,v) values of the enforced vertex on the face
443 projectionPoint myPoint = getProjectionPoint(faceShape,aPnt);
445 MESSAGE("Enforced Vertex: " << aPnt.X() << ", " << aPnt.Y() << ", " << aPnt.Z());
446 MESSAGE("Projected Vertex: " << myPoint.xyz.X() << ", " << myPoint.xyz.Y() << ", " << myPoint.xyz.Z());
447 MESSAGE("Parametric coordinates: " << myPoint.uv.X() << ", " << myPoint.uv.Y() );
449 enf_coords.push_back(aPnt.X());
450 enf_coords.push_back(aPnt.Y());
451 enf_coords.push_back(aPnt.Z());
453 coords.push_back(myPoint.uv.X());
454 coords.push_back(myPoint.uv.Y());
455 coords.push_back(myPoint.xyz.X());
456 coords.push_back(myPoint.xyz.Y());
457 coords.push_back(myPoint.xyz.Z());
459 s_coords.push_back(myPoint.xyz.X());
460 s_coords.push_back(myPoint.xyz.Y());
461 s_coords.push_back(myPoint.xyz.Z());
463 // Save pair projected vertex / enf vertex
464 MESSAGE("Storing pair projected vertex / enf vertex:");
465 MESSAGE("("<< myPoint.xyz.X() << ", " << myPoint.xyz.Y() << ", " << myPoint.xyz.Z() <<") / (" << aPnt.X() << ", " << aPnt.Y() << ", " << aPnt.Z()<<")");
466 EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
467 MESSAGE("Group name is: \"" << enfVertex->grpName << "\"");
468 EnfVertexCoords2EnfVertexList[s_coords].insert(enfVertex);
471 if (! FacesWithEnforcedVertices.Contains(faceShape)) {
472 key = FacesWithEnforcedVertices.Add(faceShape);
475 key = FacesWithEnforcedVertices.FindIndex(faceShape);
478 // If a node is already created by an attractor, do not create enforced vertex
479 int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
480 bool sameAttractor = false;
481 if (attractorKey >= 0)
482 if (FaceId2AttractorCoords.count(attractorKey) > 0)
483 if (FaceId2AttractorCoords[attractorKey] == coords)
484 sameAttractor = true;
486 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
487 MESSAGE("Map of enf. vertex has key " << key)
488 MESSAGE("Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
490 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
492 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
493 MESSAGE("New Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
496 MESSAGE("Map of enf. vertex has not key " << key << ": creating it")
497 if (! sameAttractor) {
498 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList ens;
500 FaceId2EnforcedVertexCoords[key] = ens;
503 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
507 /////////////////////////////////////////////////////////
508 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
510 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
513 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
515 for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
516 enfVertex = *enfVertexListIt;
517 // Case of manual coords
518 if (enfVertex->coords.size() != 0) {
519 aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
520 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
523 // Case of geom vertex coords
524 if (enfVertex->geomEntry != "") {
525 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
526 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
527 if (GeomType == TopAbs_VERTEX){
528 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(GeomShape));
529 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
532 if (GeomType == TopAbs_COMPOUND){
533 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
534 if (it.Value().ShapeType() == TopAbs_VERTEX){
535 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
536 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
544 /////////////////////////////////////////////////////////
545 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction)
547 MESSAGE("Attractor function: "<< AttractorFunction);
548 double xa, ya, za; // Coordinates of attractor point
549 double a, b; // Attractor parameter
551 bool createNode=false; // To create a node on attractor projection
553 const char *sep = ";";
554 // atIt->second has the following pattern:
555 // ATTRACTOR(xa;ya;za;a;b;True|False;d)
557 // xa;ya;za : coordinates of attractor
558 // a : desired size on attractor
559 // b : distance of influence of attractor
560 // d : distance until which the size remains constant
562 // We search the parameters in the string
564 pos1 = AttractorFunction.find(sep);
565 if (pos1!=string::npos)
566 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
568 pos2 = AttractorFunction.find(sep, pos1+1);
569 if (pos2!=string::npos) {
570 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
574 pos2 = AttractorFunction.find(sep, pos1+1);
575 if (pos2!=string::npos) {
576 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
580 pos2 = AttractorFunction.find(sep, pos1+1);
581 if (pos2!=string::npos) {
582 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
586 pos2 = AttractorFunction.find(sep, pos1+1);
587 if (pos2!=string::npos) {
588 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
592 pos2 = AttractorFunction.find(sep, pos1+1);
593 if (pos2!=string::npos) {
594 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
595 MESSAGE("createNode: " << createNodeStr);
596 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
600 pos2 = AttractorFunction.find(")");
601 if (pos2!=string::npos) {
602 d = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
605 // Get the (u,v) values of the attractor on the face
606 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
607 gp_XY uvPoint = myPoint.uv;
608 gp_XYZ xyzPoint = myPoint.xyz;
609 Standard_Real u0 = uvPoint.X();
610 Standard_Real v0 = uvPoint.Y();
611 Standard_Real x0 = xyzPoint.X();
612 Standard_Real y0 = xyzPoint.Y();
613 Standard_Real z0 = xyzPoint.Z();
614 std::vector<double> coords;
615 coords.push_back(u0);
616 coords.push_back(v0);
617 coords.push_back(x0);
618 coords.push_back(y0);
619 coords.push_back(z0);
620 // We construct the python function
621 ostringstream attractorFunctionStream;
622 attractorFunctionStream << "def f(u,v): return ";
623 attractorFunctionStream << _smp_phy_size << "-(" << _smp_phy_size <<"-" << a << ")";
624 //attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
625 // rnc: make possible to keep the size constant until
626 // a defined distance. Distance is expressed as the positiv part
627 // of r-d where r is the distance to (u0,v0)
628 attractorFunctionStream << "*exp(-(0.5*(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"+abs(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"))/(" << b << "))**2)";
630 MESSAGE("Python function for attractor:" << std::endl << attractorFunctionStream.str());
633 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
634 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
637 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
639 FaceId2SizeMap[key] =attractorFunctionStream.str();
641 MESSAGE("Creating node on ("<<x0<<","<<y0<<","<<z0<<")");
642 FaceId2AttractorCoords[key] = coords;
644 // // Test for new attractors
645 // gp_Pnt myP(xyzPoint);
646 // TopoDS_Vertex myV = BRepBuilderAPI_MakeVertex(myP);
647 // BLSURFPlugin_Attractor myAttractor(TopoDS::Face(GeomShape),myV,200);
648 // myAttractor.SetParameters(a, _smp_phy_size, b, d);
649 // myAttractor.SetType(1);
650 // FaceId2ClassAttractor[key] = myAttractor;
651 // if(FaceId2ClassAttractor[key].GetFace().IsNull()){
652 // MESSAGE("face nulle ");
655 // MESSAGE("face OK");
657 // if (FaceId2ClassAttractor[key].GetAttractorShape().IsNull()){
658 // MESSAGE("pas de point");
661 // MESSAGE("point OK");
664 /////////////////////////////////////////////////////////
666 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
667 blsurf_session_t * bls,
670 int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
671 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
672 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
673 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
674 double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
675 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
676 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
677 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
678 bool _decimesh = BLSURFPlugin_Hypothesis::GetDefaultDecimesh();
679 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
682 MESSAGE("BLSURFPlugin_BLSURF::SetParameters");
683 _topology = (int) hyp->GetTopology();
684 _physicalMesh = (int) hyp->GetPhysicalMesh();
685 _phySize = hyp->GetPhySize();
686 _geometricMesh = (int) hyp->GetGeometricMesh();
687 _angleMeshS = hyp->GetAngleMeshS();
688 _angleMeshC = hyp->GetAngleMeshC();
689 _gradation = hyp->GetGradation();
690 _quadAllowed = hyp->GetQuadAllowed();
691 _decimesh = hyp->GetDecimesh();
692 _verb = hyp->GetVerbosity();
694 if ( hyp->GetPhyMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
695 blsurf_set_param(bls, "hphymin", to_string(hyp->GetPhyMin()).c_str());
696 if ( hyp->GetPhyMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
697 blsurf_set_param(bls, "hphymax", to_string(hyp->GetPhyMax()).c_str());
698 if ( hyp->GetGeoMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
699 blsurf_set_param(bls, "hgeomin", to_string(hyp->GetGeoMin()).c_str());
700 if ( hyp->GetGeoMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
701 blsurf_set_param(bls, "hgeomax", to_string(hyp->GetGeoMax()).c_str());
703 const BLSURFPlugin_Hypothesis::TOptionValues & opts = hyp->GetOptionValues();
704 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
705 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt )
706 if ( !opIt->second.empty() ) {
707 MESSAGE("blsurf_set_param(): " << opIt->first << " = " << opIt->second);
708 blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
712 //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
713 // GetDefaultPhySize() sometimes leads to computation failure
714 _phySize = mesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
715 MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
717 _smp_phy_size = _phySize;
718 blsurf_set_param(bls, "topo_points", _topology > 0 ? "1" : "0");
719 blsurf_set_param(bls, "topo_curves", _topology > 0 ? "1" : "0");
720 blsurf_set_param(bls, "topo_project", _topology > 0 ? "1" : "0");
721 blsurf_set_param(bls, "clean_boundary", _topology > 1 ? "1" : "0");
722 blsurf_set_param(bls, "close_boundary", _topology > 1 ? "1" : "0");
723 blsurf_set_param(bls, "hphy_flag", to_string(_physicalMesh).c_str());
724 // blsurf_set_param(bls, "hphy_flag", "2");
725 if ((to_string(_physicalMesh))=="2"){
726 TopoDS_Shape GeomShape;
727 TopoDS_Shape AttShape;
728 TopAbs_ShapeEnum GeomType;
730 // Standard Size Maps
732 MESSAGE("Setting a Size Map");
733 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
734 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
735 for ( ; smIt != sizeMaps.end(); ++smIt ) {
736 if ( !smIt->second.empty() ) {
737 MESSAGE("blsurf_set_sizeMap(): " << smIt->first << " = " << smIt->second);
738 GeomShape = entryToShape(smIt->first);
739 GeomType = GeomShape.ShapeType();
740 MESSAGE("Geomtype is " << GeomType);
743 if (GeomType == TopAbs_COMPOUND){
744 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
746 if (it.Value().ShapeType() == TopAbs_FACE){
747 HasSizeMapOnFace = true;
748 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
749 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
752 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
753 // MESSAGE("Face with key " << key << " already in map");
755 FaceId2SizeMap[key] = smIt->second;
758 if (it.Value().ShapeType() == TopAbs_EDGE){
759 HasSizeMapOnEdge = true;
760 HasSizeMapOnFace = true;
761 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
762 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
765 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
766 // MESSAGE("Edge with key " << key << " already in map");
768 EdgeId2SizeMap[key] = smIt->second;
771 if (it.Value().ShapeType() == TopAbs_VERTEX){
772 HasSizeMapOnVertex = true;
773 HasSizeMapOnEdge = true;
774 HasSizeMapOnFace = true;
775 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
776 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
779 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
780 MESSAGE("Group of vertices with key " << key << " already in map");
782 MESSAGE("Group of vertices with key " << key << " has a size map: " << smIt->second);
783 VertexId2SizeMap[key] = smIt->second;
788 if (GeomType == TopAbs_FACE){
789 HasSizeMapOnFace = true;
790 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
791 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
794 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
795 // MESSAGE("Face with key " << key << " already in map");
797 FaceId2SizeMap[key] = smIt->second;
800 if (GeomType == TopAbs_EDGE){
801 HasSizeMapOnEdge = true;
802 HasSizeMapOnFace = true;
803 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
804 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
807 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
808 // MESSAGE("Edge with key " << key << " already in map");
810 EdgeId2SizeMap[key] = smIt->second;
813 if (GeomType == TopAbs_VERTEX){
814 HasSizeMapOnVertex = true;
815 HasSizeMapOnEdge = true;
816 HasSizeMapOnFace = true;
817 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
818 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
821 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
822 MESSAGE("Vertex with key " << key << " already in map");
824 MESSAGE("Vertex with key " << key << " has a size map: " << smIt->second);
825 VertexId2SizeMap[key] = smIt->second;
833 // TODO appeler le constructeur des attracteurs directement ici
834 MESSAGE("Setting Attractors");
835 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
836 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
837 for ( ; atIt != attractors.end(); ++atIt ) {
838 if ( !atIt->second.empty() ) {
839 MESSAGE("blsurf_set_attractor(): " << atIt->first << " = " << atIt->second);
840 GeomShape = entryToShape(atIt->first);
841 GeomType = GeomShape.ShapeType();
843 if (GeomType == TopAbs_COMPOUND){
844 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
845 if (it.Value().ShapeType() == TopAbs_FACE){
846 HasSizeMapOnFace = true;
847 createAttractorOnFace(it.Value(), atIt->second);
852 if (GeomType == TopAbs_FACE){
853 HasSizeMapOnFace = true;
854 createAttractorOnFace(GeomShape, atIt->second);
857 if (GeomType == TopAbs_EDGE){
858 HasSizeMapOnEdge = true;
859 HasSizeMapOnFace = true;
860 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
862 if (GeomType == TopAbs_VERTEX){
863 HasSizeMapOnVertex = true;
864 HasSizeMapOnEdge = true;
865 HasSizeMapOnFace = true;
866 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
873 // temporary commented out for testing
875 // - Fill in the BLSURFPlugin_Hypothesis::TAttractorMap map in the hypothesis
876 // - Uncomment and complete this part to construct the attractors from the attractor shape and the passed parameters on each face of the map
877 // - To do this use the public methodss: SetParameters(several double parameters) and SetType(int type)
879 // - Construct the attractors with an empty dist. map in the hypothesis
880 // - build the map here for each face with an attractor set and only if the attractor shape as changed since the last call to _buildmap()
881 // -> define a bool _mapbuilt in the class that is set to false by default and set to true when calling _buildmap() OK
883 const BLSURFPlugin_Hypothesis::TAttractorMap class_attractors = BLSURFPlugin_Hypothesis::GetClassAttractorEntries(hyp);
885 BLSURFPlugin_Hypothesis::TAttractorMap::const_iterator AtIt = class_attractors.begin();
886 for ( ; AtIt != class_attractors.end(); ++AtIt ) {
887 if ( !AtIt->second->Empty() ) {
888 // MESSAGE("blsurf_set_attractor(): " << AtIt->first << " = " << AtIt->second);
889 GeomShape = entryToShape(AtIt->first);
890 AttShape = AtIt->second->GetAttractorShape();
891 GeomType = GeomShape.ShapeType();
893 // if (GeomType == TopAbs_COMPOUND){
894 // for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
895 // if (it.Value().ShapeType() == TopAbs_FACE){
896 // HasAttractorOnFace = true;
897 // myAttractor = BLSURFPluginAttractor(GeomShape, AttShape);
902 if (GeomType == TopAbs_FACE
903 && (AttShape.ShapeType() == TopAbs_VERTEX
904 || AttShape.ShapeType() == TopAbs_EDGE
905 || AttShape.ShapeType() == TopAbs_WIRE
906 || AttShape.ShapeType() == TopAbs_COMPOUND) ){
907 HasSizeMapOnFace = true;
909 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape)) ) {
910 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
913 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
914 // MESSAGE("Face with key " << key << " already in map");
917 FaceId2ClassAttractor[key] = AtIt->second;
920 MESSAGE("Wrong shape type !!")
931 MESSAGE("Setting Enforced Vertices");
932 const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
933 BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
934 for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
935 if ( !enfIt->second.empty() ) {
936 GeomShape = entryToShape(enfIt->first);
937 GeomType = GeomShape.ShapeType();
939 if (GeomType == TopAbs_COMPOUND){
940 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
941 if (it.Value().ShapeType() == TopAbs_FACE){
942 HasSizeMapOnFace = true;
943 createEnforcedVertexOnFace(it.Value(), enfIt->second);
948 if (GeomType == TopAbs_FACE){
949 HasSizeMapOnFace = true;
950 createEnforcedVertexOnFace(GeomShape, enfIt->second);
955 // if (HasSizeMapOnFace){
956 // In all size map cases (hphy_flag = 2), at least map on face must be defined
957 MESSAGE("Setting Size Map on FACES ");
958 blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
961 if (HasSizeMapOnEdge){
962 MESSAGE("Setting Size Map on EDGES ");
963 blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
965 if (HasSizeMapOnVertex){
966 MESSAGE("Setting Size Map on VERTICES ");
967 blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
970 blsurf_set_param(bls, "hphydef", to_string(_phySize).c_str());
971 blsurf_set_param(bls, "hgeo_flag", to_string(_geometricMesh).c_str());
972 blsurf_set_param(bls, "relax_size", _decimesh ? "0": to_string(_geometricMesh).c_str());
973 blsurf_set_param(bls, "angle_meshs", to_string(_angleMeshS).c_str());
974 blsurf_set_param(bls, "angle_meshc", to_string(_angleMeshC).c_str());
975 blsurf_set_param(bls, "gradation", to_string(_gradation).c_str());
976 blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
977 blsurf_set_param(bls, "element", _quadAllowed ? "q1.0" : "p1");
978 blsurf_set_param(bls, "verb", to_string(_verb).c_str());
981 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
982 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
983 real *duu, real *duv, real *dvv, void *user_data);
984 status_t message_cb(message_t *msg, void *user_data);
985 status_t interrupt_cb(integer *interrupt_status, void *user_data);
987 //=============================================================================
991 //=============================================================================
993 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
995 MESSAGE("BLSURFPlugin_BLSURF::Compute");
997 // if (aShape.ShapeType() == TopAbs_COMPOUND) {
998 // MESSAGE(" the shape is a COMPOUND");
1001 // MESSAGE(" the shape is UNKNOWN");
1004 // Fix problem with locales
1005 Kernel_Utils::Localizer aLocalizer;
1007 /* create a distene context (generic object) */
1008 status_t status = STATUS_ERROR;
1010 context_t *ctx = context_new();
1012 /* Set the message callback in the working context */
1013 context_set_message_callback(ctx, message_cb, &_comment);
1014 #ifdef WITH_SMESH_CANCEL_COMPUTE
1015 _compute_canceled = false;
1016 context_set_interrupt_callback(ctx, interrupt_cb, this);
1018 context_set_interrupt_callback(ctx, interrupt_cb, NULL);
1021 /* create the CAD object we will work on. It is associated to the context ctx. */
1022 cad_t *c = cad_new(ctx);
1024 blsurf_session_t *bls = blsurf_session_new(ctx);
1026 FacesWithSizeMap.Clear();
1027 FaceId2SizeMap.clear();
1028 FaceId2ClassAttractor.clear();
1029 FaceIndex2ClassAttractor.clear();
1030 EdgesWithSizeMap.Clear();
1031 EdgeId2SizeMap.clear();
1032 VerticesWithSizeMap.Clear();
1033 VertexId2SizeMap.clear();
1035 MESSAGE("BEGIN SetParameters");
1036 SetParameters(_hypothesis, bls, aMesh);
1037 MESSAGE("END SetParameters");
1039 /* Now fill the CAD object with data from your CAD
1040 * environement. This is the most complex part of a successfull
1044 // needed to prevent the opencascade memory managmement from freeing things
1045 vector<Handle(Geom2d_Curve)> curves;
1046 vector<Handle(Geom_Surface)> surfaces;
1051 TopTools_IndexedMapOfShape fmap;
1052 TopTools_IndexedMapOfShape emap;
1053 TopTools_IndexedMapOfShape pmap;
1056 FaceId2PythonSmp.clear();
1058 EdgeId2PythonSmp.clear();
1060 VertexId2PythonSmp.clear();
1062 assert(Py_IsInitialized());
1063 PyGILState_STATE gstate;
1064 gstate = PyGILState_Ensure();
1066 string theSizeMapStr;
1068 /****************************************************************************************
1070 *****************************************************************************************/
1072 string bad_end = "return";
1075 BLSURFPlugin_Attractor myAttractor;
1076 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
1077 TopoDS_Face f=TopoDS::Face(face_iter.Current());
1079 // make INTERNAL face oriented FORWARD (issue 0020993)
1080 if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
1081 f.Orientation(TopAbs_FORWARD);
1083 if (fmap.FindIndex(f) > 0)
1088 surfaces.push_back(BRep_Tool::Surface(f));
1090 /* create an object representing the face for blsurf */
1091 /* where face_id is an integer identifying the face.
1092 * surf_function is the function that defines the surface
1093 * (For this face, it will be called by blsurf with your_face_object_ptr
1094 * as last parameter.
1096 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
1098 /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
1099 cad_face_set_tag(fce, iface);
1101 /* Set face orientation (optional if you want a well oriented output mesh)*/
1102 if(f.Orientation() != TopAbs_FORWARD){
1103 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
1105 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
1108 if (HasSizeMapOnFace){
1109 // MESSAGE("A size map is defined on a face")
1110 // std::cout << "A size map is defined on a face" << std::endl;
1112 faceKey = FacesWithSizeMap.FindIndex(f);
1115 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
1116 MESSAGE("A size map is defined on face :"<<faceKey)
1117 theSizeMapStr = FaceId2SizeMap[faceKey];
1118 // check if function ends with "return"
1119 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1121 // Expr To Python function, verification is performed at validation in GUI
1122 PyObject * obj = NULL;
1123 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1125 PyObject * func = NULL;
1126 func = PyObject_GetAttrString(main_mod, "f");
1127 FaceId2PythonSmp[iface]=func;
1128 FaceId2SizeMap.erase(faceKey);
1131 // Specific size map = Attractor
1132 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
1134 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
1135 if (attractor_iter->first == faceKey) {
1136 MESSAGE("Face indice: " << iface);
1137 MESSAGE("Adding attractor");
1139 double xyzCoords[3] = {attractor_iter->second[2],
1140 attractor_iter->second[3],
1141 attractor_iter->second[4]};
1143 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1144 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1145 BRepClass_FaceClassifier scl(f,P,1e-7);
1146 // scl.Perform() is bugged. The function was rewritten
1148 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1149 TopAbs_State result = scl.State();
1150 MESSAGE("Position of point on face: "<<result);
1151 if ( result == TopAbs_OUT )
1152 MESSAGE("Point is out of face: node is not created");
1153 if ( result == TopAbs_UNKNOWN )
1154 MESSAGE("Point position on face is unknown: node is not created");
1155 if ( result == TopAbs_ON )
1156 MESSAGE("Point is on border of face: node is not created");
1157 if ( result == TopAbs_IN )
1159 // Point is inside face and not on border
1160 MESSAGE("Point is in face: node is created");
1161 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
1163 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << iatt);
1164 cad_point_t* point_p = cad_point_new(fce, iatt, uvCoords);
1165 cad_point_set_tag(point_p, iatt);
1167 FaceId2AttractorCoords.erase(faceKey);
1172 std::map<int,BLSURFPlugin_Attractor* >::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
1173 if (clAttractor_iter != FaceId2ClassAttractor.end()){
1174 MESSAGE("Face indice: " << iface);
1175 MESSAGE("Adding attractor");
1176 FaceIndex2ClassAttractor[iface]=clAttractor_iter->second;
1177 FaceId2ClassAttractor.erase(clAttractor_iter);
1181 // Enforced Vertices
1182 faceKey = FacesWithEnforcedVertices.FindIndex(f);
1183 std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1184 if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1185 MESSAGE("Some enforced vertices are defined");
1186 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl;
1187 MESSAGE("Face indice: " << iface);
1188 MESSAGE("Adding enforced vertices");
1189 evl = evmIt->second;
1190 MESSAGE("Number of vertices to add: "<< evl.size());
1191 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
1192 for (; evlIt != evl.end(); ++evlIt) {
1193 BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
1194 xyzCoords.push_back(evlIt->at(2));
1195 xyzCoords.push_back(evlIt->at(3));
1196 xyzCoords.push_back(evlIt->at(4));
1197 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1198 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1199 BRepClass_FaceClassifier scl(f,P,1e-7);
1200 // scl.Perform() is bugged. The function was rewritten
1202 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1203 TopAbs_State result = scl.State();
1204 MESSAGE("Position of point on face: "<<result);
1205 if ( result == TopAbs_OUT ) {
1206 MESSAGE("Point is out of face: node is not created");
1207 if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1208 EnfVertexCoords2ProjVertex.erase(xyzCoords);
1209 EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1212 if ( result == TopAbs_UNKNOWN ) {
1213 MESSAGE("Point position on face is unknown: node is not created");
1214 if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1215 EnfVertexCoords2ProjVertex.erase(xyzCoords);
1216 EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1219 if ( result == TopAbs_ON ) {
1220 MESSAGE("Point is on border of face: node is not created");
1221 if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1222 EnfVertexCoords2ProjVertex.erase(xyzCoords);
1223 EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1226 if ( result == TopAbs_IN )
1228 // Point is inside face and not on border
1229 MESSAGE("Point is in face: node is created");
1230 double uvCoords[2] = {evlIt->at(0),evlIt->at(1)};
1232 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1233 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1234 cad_point_set_tag(point_p, ienf);
1237 FaceId2EnforcedVertexCoords.erase(faceKey);
1240 // std::cout << "No enforced vertex defined" << std::endl;
1244 /****************************************************************************************
1246 now create the edges associated to this face
1247 *****************************************************************************************/
1249 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1250 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1251 int ic = emap.FindIndex(e);
1256 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1258 if (HasSizeMapOnEdge){
1259 edgeKey = EdgesWithSizeMap.FindIndex(e);
1260 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1261 MESSAGE("Sizemap defined on edge with index " << edgeKey);
1262 theSizeMapStr = EdgeId2SizeMap[edgeKey];
1263 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1265 // Expr To Python function, verification is performed at validation in GUI
1266 PyObject * obj = NULL;
1267 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1269 PyObject * func = NULL;
1270 func = PyObject_GetAttrString(main_mod, "f");
1271 EdgeId2PythonSmp[ic]=func;
1272 EdgeId2SizeMap.erase(edgeKey);
1276 /* attach the edge to the current blsurf face */
1277 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1279 /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1280 cad_edge_set_tag(edg, ic);
1282 /* by default, an edge does not necessalry appear in the resulting mesh,
1283 unless the following property is set :
1285 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1287 /* by default an edge is a boundary edge */
1288 if (e.Orientation() == TopAbs_INTERNAL)
1289 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1293 gp_Pnt2d e0 = curves.back()->Value(tmin);
1294 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1295 Standard_Real d1=0,d2=0;
1298 /****************************************************************************************
1300 *****************************************************************************************/
1302 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1303 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1307 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1310 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1312 *ip = pmap.FindIndex(v);
1316 if (HasSizeMapOnVertex){
1317 vertexKey = VerticesWithSizeMap.FindIndex(v);
1318 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1319 theSizeMapStr = VertexId2SizeMap[vertexKey];
1320 //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1321 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1323 // Expr To Python function, verification is performed at validation in GUI
1324 PyObject * obj = NULL;
1325 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1327 PyObject * func = NULL;
1328 func = PyObject_GetAttrString(main_mod, "f");
1329 VertexId2PythonSmp[*ip]=func;
1330 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
1335 // should not happen
1336 MESSAGE("An edge does not have 2 extremities.");
1339 // This defines the curves extremity connectivity
1340 cad_edge_set_extremities(edg, ip1, ip2);
1341 /* set the tag (color) to the same value as the extremity id : */
1342 cad_edge_set_extremities_tag(edg, ip1, ip2);
1345 cad_edge_set_extremities(edg, ip2, ip1);
1346 cad_edge_set_extremities_tag(edg, ip2, ip1);
1353 PyGILState_Release(gstate);
1355 blsurf_data_set_cad(bls, c);
1357 std::cout << std::endl;
1358 std::cout << "Beginning of Surface Mesh generation" << std::endl;
1359 std::cout << std::endl;
1361 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1363 feclearexcept( FE_ALL_EXCEPT );
1364 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1370 status = blsurf_compute_mesh(bls);
1373 catch ( std::exception& exc ) {
1374 _comment += exc.what();
1376 catch (Standard_Failure& ex) {
1377 _comment += ex.DynamicType()->Name();
1378 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1380 _comment += ex.GetMessageString();
1384 if ( _comment.empty() )
1385 _comment = "Exception in blsurf_compute_mesh()";
1387 if ( status != STATUS_OK) {
1388 // There was an error while meshing
1389 blsurf_session_delete(bls);
1391 context_delete(ctx);
1393 return error(_comment);
1396 std::cout << std::endl;
1397 std::cout << "End of Surface Mesh generation" << std::endl;
1398 std::cout << std::endl;
1401 blsurf_data_get_mesh(bls, &msh);
1403 blsurf_session_delete(bls);
1405 context_delete(ctx);
1407 return error(_comment);
1411 /* retrieve mesh data (see distene/mesh.h) */
1412 integer nv, ne, nt, nq, vtx[4], tag;
1415 mesh_get_vertex_count(msh, &nv);
1416 mesh_get_edge_count(msh, &ne);
1417 mesh_get_triangle_count(msh, &nt);
1418 mesh_get_quadrangle_count(msh, &nq);
1421 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1422 SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1423 bool* tags = new bool[nv+1];
1425 /* enumerated vertices */
1426 for(int iv=1;iv<=nv;iv++) {
1427 mesh_get_vertex_coordinates(msh, iv, xyz);
1428 mesh_get_vertex_tag(msh, iv, &tag);
1429 // Issue 0020656. Use vertex coordinates
1430 if ( tag > 0 && tag <= pmap.Extent() ) {
1431 TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
1432 double tol = BRep_Tool::Tolerance( v );
1433 gp_Pnt p = BRep_Tool::Pnt( v );
1434 if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
1435 xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
1437 tag = 0; // enforced or attracted vertex
1439 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1441 // Create group of enforced vertices if requested
1443 BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
1445 projVertex.push_back((double)xyz[0]);
1446 projVertex.push_back((double)xyz[1]);
1447 projVertex.push_back((double)xyz[2]);
1448 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
1449 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end()) {
1450 MESSAGE("Found enforced vertex @ " << xyz[0] << ", " << xyz[1] << ", " << xyz[2])
1451 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
1452 BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
1453 for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
1454 currentEnfVertex = (*enfListIt);
1455 if (currentEnfVertex->grpName != "") {
1456 bool groupDone = false;
1457 SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
1458 MESSAGE("currentEnfVertex->grpName: " << currentEnfVertex->grpName);
1459 MESSAGE("Parsing the groups of the mesh");
1460 while (grIt->more()) {
1461 SMESH_Group * group = grIt->next();
1462 if ( !group ) continue;
1463 MESSAGE("Group: " << group->GetName());
1464 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
1465 if ( !groupDS ) continue;
1466 MESSAGE("group->SMDSGroup().GetType(): " << (groupDS->GetType()));
1467 MESSAGE("group->SMDSGroup().GetType()==SMDSAbs_Node: " << (groupDS->GetType()==SMDSAbs_Node));
1468 MESSAGE("currentEnfVertex->grpName.compare(group->GetStoreName())==0: " << (currentEnfVertex->grpName.compare(group->GetName())==0));
1469 if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
1470 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
1471 aGroupDS->SMDSGroup().Add(nodes[iv]);
1472 MESSAGE("Node ID: " << nodes[iv]->GetID());
1473 // How can I inform the hypothesis ?
1474 // _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
1476 MESSAGE("Successfully added enforced vertex to existing group " << currentEnfVertex->grpName);
1483 SMESH_Group* aGroup = aMesh.AddGroup(SMDSAbs_Node, currentEnfVertex->grpName.c_str(), groupId);
1484 aGroup->SetName( currentEnfVertex->grpName.c_str() );
1485 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1486 aGroupDS->SMDSGroup().Add(nodes[iv]);
1487 MESSAGE("Successfully created enforced vertex group " << currentEnfVertex->grpName);
1491 throw SALOME_Exception(LOCALIZED("A enforced vertex node was not added to a group"));
1494 MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
1500 // internal point are tagged to zero
1501 if(tag > 0 && tag <= pmap.Extent() ){
1502 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1509 /* enumerate edges */
1510 for(int it=1;it<=ne;it++) {
1511 mesh_get_edge_vertices(msh, it, vtx);
1512 SMDS_MeshEdge* edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1513 mesh_get_edge_tag(msh, it, &tag);
1516 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1517 tags[vtx[0]] = false;
1520 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1521 tags[vtx[1]] = false;
1523 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1527 /* enumerate triangles */
1528 for(int it=1;it<=nt;it++) {
1529 mesh_get_triangle_vertices(msh, it, vtx);
1530 SMDS_MeshFace* tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1531 mesh_get_triangle_tag(msh, it, &tag);
1532 meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1534 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1535 tags[vtx[0]] = false;
1538 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1539 tags[vtx[1]] = false;
1542 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1543 tags[vtx[2]] = false;
1547 /* enumerate quadrangles */
1548 for(int it=1;it<=nq;it++) {
1549 mesh_get_quadrangle_vertices(msh, it, vtx);
1550 SMDS_MeshFace* quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1551 mesh_get_quadrangle_tag(msh, it, &tag);
1552 meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1554 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1555 tags[vtx[0]] = false;
1558 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1559 tags[vtx[1]] = false;
1562 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1563 tags[vtx[2]] = false;
1566 meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1567 tags[vtx[3]] = false;
1571 // SetIsAlwaysComputed( true ) to sub-meshes of degenerated EDGEs
1572 TopLoc_Location loc; double f,l;
1573 for (int i = 1; i <= emap.Extent(); i++)
1574 if ( BRep_Tool::Curve(TopoDS::Edge( emap( i )), loc, f,l).IsNull() )
1575 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
1576 sm->SetIsAlwaysComputed( true );
1580 /* release the mesh object */
1581 blsurf_data_regain_mesh(bls, msh);
1583 /* clean up everything */
1584 blsurf_session_delete(bls);
1587 context_delete(ctx);
1589 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1591 if ( oldFEFlags > 0 )
1592 feenableexcept( oldFEFlags );
1593 feclearexcept( FE_ALL_EXCEPT );
1597 std::cout << "FacesWithSizeMap" << std::endl;
1598 FacesWithSizeMap.Statistics(std::cout);
1599 std::cout << "EdgesWithSizeMap" << std::endl;
1600 EdgesWithSizeMap.Statistics(std::cout);
1601 std::cout << "VerticesWithSizeMap" << std::endl;
1602 VerticesWithSizeMap.Statistics(std::cout);
1603 std::cout << "FacesWithEnforcedVertices" << std::endl;
1604 FacesWithEnforcedVertices.Statistics(std::cout);
1610 #ifdef WITH_SMESH_CANCEL_COMPUTE
1611 void BLSURFPlugin_BLSURF::CancelCompute()
1613 _compute_canceled = true;
1617 //=============================================================================
1621 //=============================================================================
1623 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1624 const TopoDS_Edge edge = TopoDS::Edge(ed);
1626 gp_Pnt pnt(node->X(), node->Y(), node->Z());
1628 Standard_Real p0 = 0.0;
1629 Standard_Real p1 = 1.0;
1630 TopLoc_Location loc;
1631 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1633 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1634 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1637 if ( proj.NbPoints() > 0 )
1639 pa = (double)proj.LowerDistanceParameter();
1640 // Issue 0020656. Move node if it is too far from edge
1641 gp_Pnt curve_pnt = curve->Value( pa );
1642 double dist2 = pnt.SquareDistance( curve_pnt );
1643 double tol = BRep_Tool::Tolerance( edge );
1644 if ( 1e-14 < dist2 && dist2 <= 1000*tol ) // large enough and within tolerance
1646 curve_pnt.Transform( loc );
1647 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1650 // GProp_GProps LProps;
1651 // BRepGProp::LinearProperties(ed, LProps);
1652 // double lg = (double)LProps.Mass();
1654 meshDS->SetNodeOnEdge(node, edge, pa);
1657 //=============================================================================
1661 //=============================================================================
1663 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
1668 //=============================================================================
1672 //=============================================================================
1674 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
1679 //=============================================================================
1683 //=============================================================================
1685 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
1687 return hyp.SaveTo( save );
1690 //=============================================================================
1694 //=============================================================================
1696 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
1698 return hyp.LoadFrom( load );
1701 /* Curve definition function See cad_curv_t in file distene/cad.h for
1703 * NOTE : if when your CAD systems evaluates second
1704 * order derivatives it also computes first order derivatives and
1705 * function evaluation, you can optimize this example by making only
1706 * one CAD call and filling the necessary uv, dt, dtt arrays.
1708 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
1710 /* t is given. It contains the t (time) 1D parametric coordintaes
1711 of the point PreCAD/BLSurf is querying on the curve */
1713 /* user_data identifies the edge PreCAD/BLSurf is querying
1714 * (see cad_edge_new later in this example) */
1715 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
1718 /* BLSurf is querying the function evaluation */
1721 uv[0]=P.X(); uv[1]=P.Y();
1725 /* query for the first order derivatives */
1728 dt[0]=V1.X(); dt[1]=V1.Y();
1732 /* query for the second order derivatives */
1735 dtt[0]=V2.X(); dtt[1]=V2.Y();
1741 /* Surface definition function.
1742 * See cad_surf_t in file distene/cad.h for more information.
1743 * NOTE : if when your CAD systems evaluates second order derivatives it also
1744 * computes first order derivatives and function evaluation, you can optimize
1745 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
1748 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1749 real *duu, real *duv, real *dvv, void *user_data)
1751 /* uv[2] is given. It contains the u,v coordinates of the point
1752 * PreCAD/BLSurf is querying on the surface */
1754 /* user_data identifies the face PreCAD/BLSurf is querying (see
1755 * cad_face_new later in this example)*/
1756 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
1760 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
1761 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
1768 geometry->D1(uv[0],uv[1],P,D1U,D1V);
1769 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
1770 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
1773 if(duu && duv && dvv){
1777 gp_Vec D2U,D2V,D2UV;
1779 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
1780 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
1781 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
1782 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
1789 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
1792 if (my_u_min > uv[0]) {
1795 if (my_v_min > uv[1]) {
1798 if (my_u_max < uv[0]) {
1801 if (my_v_max < uv[1]) {
1805 //MESSAGE("size_on_surface")
1806 if (FaceId2PythonSmp.count(face_id) != 0){
1807 //MESSAGE("A size map is used to calculate size on face : "<<face_id)
1808 PyObject * pyresult = NULL;
1809 PyObject* new_stderr = NULL;
1810 assert(Py_IsInitialized());
1811 PyGILState_STATE gstate;
1812 gstate = PyGILState_Ensure();
1813 pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
1815 if ( pyresult == NULL){
1817 string err_description="";
1818 new_stderr = newPyStdOut(err_description);
1819 PySys_SetObject((char*)"stderr", new_stderr);
1821 PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
1822 Py_DECREF(new_stderr);
1823 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
1824 result = *((double*)user_data);
1827 result = PyFloat_AsDouble(pyresult);
1828 Py_DECREF(pyresult);
1830 // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1832 PyGILState_Release(gstate);
1834 else if (FaceIndex2ClassAttractor.count(face_id) !=0 && !FaceIndex2ClassAttractor[face_id]->Empty()){
1835 // MESSAGE("attractor used on face :"<<face_id)
1836 // MESSAGE("List of attractor is not empty")
1837 // MESSAGE("Attractor empty : "<< FaceIndex2ClassAttractor[face_id]->Empty())
1838 double result = FaceIndex2ClassAttractor[face_id]->GetSize(uv[0],uv[1]);
1840 // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1843 // MESSAGE("List of attractor is empty !!!")
1844 *size = *((double*)user_data);
1849 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
1851 if (EdgeId2PythonSmp.count(edge_id) != 0){
1852 PyObject * pyresult = NULL;
1853 PyObject* new_stderr = NULL;
1854 assert(Py_IsInitialized());
1855 PyGILState_STATE gstate;
1856 gstate = PyGILState_Ensure();
1857 pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
1859 if ( pyresult == NULL){
1861 string err_description="";
1862 new_stderr = newPyStdOut(err_description);
1863 PySys_SetObject((char*)"stderr", new_stderr);
1865 PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
1866 Py_DECREF(new_stderr);
1867 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
1868 result = *((double*)user_data);
1871 result = PyFloat_AsDouble(pyresult);
1872 Py_DECREF(pyresult);
1875 PyGILState_Release(gstate);
1878 *size = *((double*)user_data);
1883 status_t size_on_vertex(integer point_id, real *size, void *user_data)
1885 if (VertexId2PythonSmp.count(point_id) != 0){
1886 PyObject * pyresult = NULL;
1887 PyObject* new_stderr = NULL;
1888 assert(Py_IsInitialized());
1889 PyGILState_STATE gstate;
1890 gstate = PyGILState_Ensure();
1891 pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
1893 if ( pyresult == NULL){
1895 string err_description="";
1896 new_stderr = newPyStdOut(err_description);
1897 PySys_SetObject((char*)"stderr", new_stderr);
1899 PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
1900 Py_DECREF(new_stderr);
1901 MESSAGE("Can't evaluate f()" << " error is " << err_description);
1902 result = *((double*)user_data);
1905 result = PyFloat_AsDouble(pyresult);
1906 Py_DECREF(pyresult);
1909 PyGILState_Release(gstate);
1912 *size = *((double*)user_data);
1918 * The following function will be called for PreCAD/BLSurf message
1919 * printing. See context_set_message_callback (later in this
1920 * template) for how to set user_data.
1922 status_t message_cb(message_t *msg, void *user_data)
1924 integer errnumber = 0;
1926 message_get_number(msg, &errnumber);
1927 message_get_description(msg, &desc);
1928 if ( errnumber < 0 ) {
1929 string * error = (string*)user_data;
1930 // if ( !error->empty() )
1932 // remove ^A from the tail
1933 int len = strlen( desc );
1934 while (len > 0 && desc[len-1] != '\n')
1936 error->append( desc, len );
1939 std::cout << desc << std::endl;
1944 /* This is the interrupt callback. PreCAD/BLSurf will call this
1945 * function regularily. See the file distene/interrupt.h
1947 status_t interrupt_cb(integer *interrupt_status, void *user_data)
1949 integer you_want_to_continue = 1;
1950 #ifdef WITH_SMESH_CANCEL_COMPUTE
1951 BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
1952 you_want_to_continue = !tmp->computeCanceled();
1955 if(you_want_to_continue)
1957 *interrupt_status = INTERRUPT_CONTINUE;
1960 else /* you want to stop BLSurf */
1962 *interrupt_status = INTERRUPT_STOP;
1963 return STATUS_ERROR;
1967 //=============================================================================
1971 //=============================================================================
1972 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
1973 const TopoDS_Shape& aShape,
1974 MapShapeNbElems& aResMap)
1976 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
1977 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
1978 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
1979 //double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
1980 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
1981 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
1983 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
1984 _phySize = _hypothesis->GetPhySize();
1985 //_geometricMesh = (int) hyp->GetGeometricMesh();
1986 //_angleMeshS = hyp->GetAngleMeshS();
1987 _angleMeshC = _hypothesis->GetAngleMeshC();
1988 _quadAllowed = _hypothesis->GetQuadAllowed();
1991 bool IsQuadratic = false;
1996 TopTools_DataMapOfShapeInteger EdgesMap;
1997 double fullLen = 0.0;
1998 double fullNbSeg = 0;
1999 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2000 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
2001 if( EdgesMap.IsBound(E) )
2003 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
2004 double aLen = SMESH_Algo::EdgeLength(E);
2007 if(_physicalMesh==1) {
2008 nb1d = (int)( aLen/_phySize + 1 );
2013 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
2014 double fullAng = 0.0;
2015 double dp = (l-f)/200;
2020 for(int j=2; j<=200; j++) {
2023 fullAng += fabs(V1.Angle(V2));
2027 nb1d = (int)( fullAng/_angleMeshC + 1 );
2030 std::vector<int> aVec(SMDSEntity_Last);
2031 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2032 if( IsQuadratic > 0 ) {
2033 aVec[SMDSEntity_Node] = 2*nb1d - 1;
2034 aVec[SMDSEntity_Quad_Edge] = nb1d;
2037 aVec[SMDSEntity_Node] = nb1d - 1;
2038 aVec[SMDSEntity_Edge] = nb1d;
2040 aResMap.insert(std::make_pair(sm,aVec));
2041 EdgesMap.Bind(E,nb1d);
2043 double ELen = fullLen/fullNbSeg;
2047 // try to evaluate as in MEFISTO
2048 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2049 TopoDS_Face F = TopoDS::Face( exp.Current() );
2050 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2052 BRepGProp::SurfaceProperties(F,G);
2053 double anArea = G.Mass();
2055 std::vector<int> nb1dVec;
2056 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
2057 int nbSeg = EdgesMap.Find(exp1.Current());
2059 nb1dVec.push_back( nbSeg );
2062 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
2063 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
2066 if ( nb1dVec.size() == 4 ) // quadrangle geom face
2068 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
2070 nbNodes = (n1 + 1) * (n2 + 1);
2075 nbTria = nbQuad = nbTria / 3 + 1;
2078 std::vector<int> aVec(SMDSEntity_Last,0);
2080 int nb1d_in = (nbTria*3 - nb1d) / 2;
2081 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
2082 aVec[SMDSEntity_Quad_Triangle] = nbTria;
2083 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
2086 aVec[SMDSEntity_Node] = nbNodes;
2087 aVec[SMDSEntity_Triangle] = nbTria;
2088 aVec[SMDSEntity_Quadrangle] = nbQuad;
2090 aResMap.insert(std::make_pair(sm,aVec));
2097 BRepGProp::VolumeProperties(aShape,G);
2098 double aVolume = G.Mass();
2099 double tetrVol = 0.1179*ELen*ELen*ELen;
2100 int nbVols = int(aVolume/tetrVol);
2101 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
2102 std::vector<int> aVec(SMDSEntity_Last);
2103 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2105 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
2106 aVec[SMDSEntity_Quad_Tetra] = nbVols;
2109 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
2110 aVec[SMDSEntity_Tetra] = nbVols;
2112 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2113 aResMap.insert(std::make_pair(sm,aVec));
2118 //=============================================================================
2120 * Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
2121 * Following line was added:
2122 * myExtrem.Perform(P);
2124 //=============================================================================
2125 void BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
2126 const TopoDS_Face& face,
2128 const Standard_Real Tol)
2130 //-- Voir BRepExtrema_ExtPF.cxx
2131 BRepAdaptor_Surface Surf(face);
2132 Standard_Real U1, U2, V1, V2;
2133 BRepTools::UVBounds(face, U1, U2, V1, V2);
2134 Extrema_ExtPS myExtrem;
2135 myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
2136 myExtrem.Perform(P);
2137 //----------------------------------------------------------
2138 //-- On cherche le point le plus proche , PUIS
2139 //-- On le classifie.
2140 Standard_Integer nbv = 0; // xpu
2141 Standard_Real MaxDist = RealLast();
2142 Standard_Integer indice = 0;
2143 if(myExtrem.IsDone()) {
2144 nbv = myExtrem.NbExt();
2145 for (Standard_Integer i = 1; i <= nbv; i++) {
2146 Standard_Real d = myExtrem.Value(i);
2156 Standard_Real U1,U2;
2157 myExtrem.Point(indice).Parameter(U1, U2);
2158 Puv.SetCoord(U1, U2);
2159 fc->Perform(face, Puv, Tol);
2162 fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
2163 //-- le resultat est TopAbs_OUT;