1 // Copyright (C) 2007-2010 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"
31 #include <distene/api.h>
32 #include <distene/blsurf.h>
35 #include <structmember.h>
38 #include <Basics_Utils.hxx>
40 #include <SMESH_Gen.hxx>
41 #include <SMESH_Mesh.hxx>
42 #include <SMESH_ControlsDef.hxx>
44 #include <utilities.h>
52 // OPENCASCADE includes
53 #include <BRep_Tool.hxx>
55 #include <TopExp_Explorer.hxx>
57 #include <NCollection_Map.hxx>
59 #include <Geom_Surface.hxx>
60 #include <Handle_Geom_Surface.hxx>
61 #include <Geom2d_Curve.hxx>
62 #include <Handle_Geom2d_Curve.hxx>
63 #include <Geom_Curve.hxx>
64 #include <Handle_Geom_Curve.hxx>
65 #include <Handle_AIS_InteractiveObject.hxx>
66 #include <TopoDS_Vertex.hxx>
67 #include <TopoDS_Edge.hxx>
68 #include <TopoDS_Wire.hxx>
69 #include <TopoDS_Face.hxx>
71 #include <gp_Pnt2d.hxx>
72 #include <TopTools_IndexedMapOfShape.hxx>
73 #include <TopoDS_Shape.hxx>
74 #include <BRep_Builder.hxx>
75 #include <BRepTools.hxx>
77 #include <TopTools_DataMapOfShapeInteger.hxx>
78 #include <GProp_GProps.hxx>
79 #include <BRepGProp.hxx>
85 #include <Standard_ErrorHandler.hxx>
86 #include <GeomAPI_ProjectPointOnCurve.hxx>
87 #include <GeomAPI_ProjectPointOnSurf.hxx>
90 // #include <BRepClass_FaceClassifier.hxx>
91 #include <TopTools_MapOfShape.hxx>
93 /* ==================================
94 * =========== PYTHON ==============
95 * ==================================*/
104 PyStdOut_dealloc(PyStdOut *self)
110 PyStdOut_write(PyStdOut *self, PyObject *args)
114 if (!PyArg_ParseTuple(args, "t#:write",&c, &l))
118 *(self->out)=*(self->out)+c;
124 static PyMethodDef PyStdOut_methods[] = {
125 {"write", (PyCFunction)PyStdOut_write, METH_VARARGS,
126 PyDoc_STR("write(string) -> None")},
127 {NULL, NULL} /* sentinel */
130 static PyMemberDef PyStdOut_memberlist[] = {
131 {"softspace", T_INT, offsetof(PyStdOut, softspace), 0,
132 "flag indicating that a space needs to be printed; used by print"},
133 {NULL} /* Sentinel */
136 static PyTypeObject PyStdOut_Type = {
137 /* The ob_type field must be initialized in the module init function
138 * to be portable to Windows without using C++. */
139 PyObject_HEAD_INIT(NULL)
142 sizeof(PyStdOut), /*tp_basicsize*/
145 (destructor)PyStdOut_dealloc, /*tp_dealloc*/
152 0, /*tp_as_sequence*/
157 PyObject_GenericGetAttr, /*tp_getattro*/
158 /* softspace is writable: we must supply tp_setattro */
159 PyObject_GenericSetAttr, /* tp_setattro */
161 Py_TPFLAGS_DEFAULT, /*tp_flags*/
165 0, /*tp_richcompare*/
166 0, /*tp_weaklistoffset*/
169 PyStdOut_methods, /*tp_methods*/
170 PyStdOut_memberlist, /*tp_members*/
184 PyObject * newPyStdOut( std::string& out )
187 self = PyObject_New(PyStdOut, &PyStdOut_Type);
192 return (PyObject*)self;
196 ////////////////////////END PYTHON///////////////////////////
198 //////////////////MY MAPS////////////////////////////////////////
199 TopTools_IndexedMapOfShape FacesWithSizeMap;
200 std::map<int,string> FaceId2SizeMap;
201 TopTools_IndexedMapOfShape EdgesWithSizeMap;
202 std::map<int,string> EdgeId2SizeMap;
203 TopTools_IndexedMapOfShape VerticesWithSizeMap;
204 std::map<int,string> VertexId2SizeMap;
206 std::map<int,PyObject*> FaceId2PythonSmp;
207 std::map<int,PyObject*> EdgeId2PythonSmp;
208 std::map<int,PyObject*> VertexId2PythonSmp;
210 std::map<int,std::vector<double> > FaceId2AttractorCoords;
212 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
213 std::map< int, std::set< std::vector<double> > > FaceId2EnforcedVertexCoords;
215 bool HasSizeMapOnFace=false;
216 bool HasSizeMapOnEdge=false;
217 bool HasSizeMapOnVertex=false;
219 //=============================================================================
223 //=============================================================================
225 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId, int studyId,
227 : SMESH_2D_Algo(hypId, studyId, gen)
229 MESSAGE("BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF");
232 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
233 _compatibleHypothesis.push_back("BLSURF_Parameters");
234 _requireDescretBoundary = false;
235 _onlyUnaryInput = false;
238 smeshGen_i = SMESH_Gen_i::GetSMESHGen();
239 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
240 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
242 MESSAGE("studyid = " << _studyId);
245 myStudy = aStudyMgr->GetStudyByID(_studyId);
246 MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
248 /* Initialize the Python interpreter */
249 assert(Py_IsInitialized());
250 PyGILState_STATE gstate;
251 gstate = PyGILState_Ensure();
254 main_mod = PyImport_AddModule("__main__");
257 main_dict = PyModule_GetDict(main_mod);
259 PyRun_SimpleString("from math import *");
260 PyGILState_Release(gstate);
262 FacesWithSizeMap.Clear();
263 FaceId2SizeMap.clear();
264 EdgesWithSizeMap.Clear();
265 EdgeId2SizeMap.clear();
266 VerticesWithSizeMap.Clear();
267 VertexId2SizeMap.clear();
268 FaceId2PythonSmp.clear();
269 EdgeId2PythonSmp.clear();
270 VertexId2PythonSmp.clear();
271 FaceId2AttractorCoords.clear();
272 FacesWithEnforcedVertices.Clear();
273 FaceId2EnforcedVertexCoords.clear();
277 //=============================================================================
281 //=============================================================================
283 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
285 MESSAGE("BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF");
289 //=============================================================================
293 //=============================================================================
295 bool BLSURFPlugin_BLSURF::CheckHypothesis
297 const TopoDS_Shape& aShape,
298 SMESH_Hypothesis::Hypothesis_Status& aStatus)
302 list<const SMESHDS_Hypothesis*>::const_iterator itl;
303 const SMESHDS_Hypothesis* theHyp;
305 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
306 int nbHyp = hyps.size();
309 aStatus = SMESH_Hypothesis::HYP_OK;
310 return true; // can work with no hypothesis
314 theHyp = (*itl); // use only the first hypothesis
316 string hypName = theHyp->GetName();
318 if (hypName == "BLSURF_Parameters")
320 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
322 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
323 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
324 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
325 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
327 aStatus = SMESH_Hypothesis::HYP_OK;
330 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
332 return aStatus == SMESH_Hypothesis::HYP_OK;
335 //=============================================================================
337 * Pass parameters to BLSURF
339 //=============================================================================
341 inline std::string to_string(double d)
343 std::ostringstream o;
348 inline std::string to_string(int i)
350 std::ostringstream o;
355 double _smp_phy_size;
356 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
357 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
358 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
360 double my_u_min=1e6,my_v_min=1e6,my_u_max=-1e6,my_v_max=-1e6;
366 /////////////////////////////////////////////////////////
367 projectionPoint getProjectionPoint(const TopoDS_Face& face, const gp_XYZ& point)
369 projectionPoint myPoint;
370 Handle(Geom_Surface) surface = BRep_Tool::Surface(face);
371 GeomAPI_ProjectPointOnSurf projector( point, surface );
372 if ( !projector.IsDone() || projector.NbPoints()==0 )
373 throw "getProjectionPoint: Can't project";
375 Quantity_Parameter u,v;
376 projector.LowerDistanceParameters(u,v);
377 myPoint.uv = gp_XY(u,v);
378 gp_Pnt aPnt = projector.NearestPoint();
379 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
383 /////////////////////////////////////////////////////////
385 /////////////////////////////////////////////////////////
386 double getT(const TopoDS_Edge& edge, const gp_XYZ& point)
389 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, f,l);
390 GeomAPI_ProjectPointOnCurve projector( point, curve);
391 if ( projector.NbPoints() == 0 )
393 return projector.LowerDistanceParameter();
396 /////////////////////////////////////////////////////////
397 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
399 MESSAGE("BLSURFPlugin_BLSURF::entryToShape "<<entry );
400 GEOM::GEOM_Object_var aGeomObj;
401 TopoDS_Shape S = TopoDS_Shape();
402 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
403 SALOMEDS::GenericAttribute_var anAttr;
405 if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR")) {
406 SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
407 CORBA::String_var aVal = anIOR->Value();
408 CORBA::Object_var obj = myStudy->ConvertIORToObject(aVal);
409 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
411 if ( !aGeomObj->_is_nil() )
412 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
416 /////////////////////////////////////////////////////////
417 void createEnforcedVertexOnFace(TopoDS_Shape GeomShape, BLSURFPlugin_Hypothesis::TEnforcedVertexList enforcedVertexList)
420 std::vector<double> coords;
421 BLSURFPlugin_Hypothesis::TEnforcedVertex enforcedVertex;
422 BLSURFPlugin_Hypothesis::TEnforcedVertexList::const_iterator evlIt = enforcedVertexList.begin();
424 for( ; evlIt != enforcedVertexList.end() ; ++evlIt ) {
426 enforcedVertex = *evlIt;
427 xe = enforcedVertex[0];
428 ye = enforcedVertex[1];
429 ze = enforcedVertex[2];
430 MESSAGE("Enforced Vertex: " << xe << ", " << ye << ", " << ze);
431 // Get the (u,v) values of the enforced vertex on the face
432 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_XYZ(xe,ye,ze));
433 gp_XY uvPoint = myPoint.uv;
434 gp_XYZ xyzPoint = myPoint.xyz;
435 Standard_Real u0 = uvPoint.X();
436 Standard_Real v0 = uvPoint.Y();
437 Standard_Real x0 = xyzPoint.X();
438 Standard_Real y0 = xyzPoint.Y();
439 Standard_Real z0 = xyzPoint.Z();
440 MESSAGE("Projected Vertex: " << x0 << ", " << y0 << ", " << z0);
441 MESSAGE("Parametric coordinates: " << u0 << ", " << v0 );
442 coords.push_back(u0);
443 coords.push_back(v0);
444 coords.push_back(x0);
445 coords.push_back(y0);
446 coords.push_back(z0);
449 if (! FacesWithEnforcedVertices.Contains(TopoDS::Face(GeomShape))) {
450 key = FacesWithEnforcedVertices.Add(TopoDS::Face(GeomShape));
453 key = FacesWithEnforcedVertices.FindIndex(TopoDS::Face(GeomShape));
456 // If a node is already created by an attractor, do not create enforced vertex
457 int attractorKey = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
458 bool sameAttractor = false;
459 if (attractorKey >= 0)
460 if (FaceId2AttractorCoords.count(attractorKey) > 0)
461 if (FaceId2AttractorCoords[attractorKey] == coords)
462 sameAttractor = true;
464 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
465 MESSAGE("Map of enf. vertex has key " << key)
466 MESSAGE("Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
468 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
470 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
471 MESSAGE("New Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
474 MESSAGE("Map of enf. vertex has not key " << key << ": creating it")
475 if (! sameAttractor) {
476 std::set< std::vector<double> > ens;
478 FaceId2EnforcedVertexCoords[key] = ens;
481 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
486 /////////////////////////////////////////////////////////
487 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction)
489 MESSAGE("Attractor function: "<< AttractorFunction);
490 double xa, ya, za; // Coordinates of attractor point
491 double a, b; // Attractor parameter
492 bool createNode=false; // To create a node on attractor projection
494 const char *sep = ";";
495 // atIt->second has the following pattern:
496 // ATTRACTOR(xa;ya;za;a;b)
498 // xa;ya;za : coordinates of attractor
499 // a : desired size on attractor
500 // b : distance of influence of attractor
502 // We search the parameters in the string
504 pos1 = AttractorFunction.find(sep);
505 if (pos1!=string::npos)
506 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
508 pos2 = AttractorFunction.find(sep, pos1+1);
509 if (pos2!=string::npos) {
510 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
514 pos2 = AttractorFunction.find(sep, pos1+1);
515 if (pos2!=string::npos) {
516 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
520 pos2 = AttractorFunction.find(sep, pos1+1);
521 if (pos2!=string::npos) {
522 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
526 pos2 = AttractorFunction.find(sep, pos1+1);
527 if (pos2!=string::npos) {
528 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
532 pos2 = AttractorFunction.find(")");
533 if (pos2!=string::npos) {
534 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
535 MESSAGE("createNode: " << createNodeStr);
536 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
539 // Get the (u,v) values of the attractor on the face
540 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_XYZ(xa,ya,za));
541 gp_XY uvPoint = myPoint.uv;
542 gp_XYZ xyzPoint = myPoint.xyz;
543 Standard_Real u0 = uvPoint.X();
544 Standard_Real v0 = uvPoint.Y();
545 Standard_Real x0 = xyzPoint.X();
546 Standard_Real y0 = xyzPoint.Y();
547 Standard_Real z0 = xyzPoint.Z();
548 std::vector<double> coords;
549 coords.push_back(u0);
550 coords.push_back(v0);
551 coords.push_back(x0);
552 coords.push_back(y0);
553 coords.push_back(z0);
554 // We construct the python function
555 ostringstream attractorFunctionStream;
556 attractorFunctionStream << "def f(u,v): return ";
557 attractorFunctionStream << _smp_phy_size << "-(" << _smp_phy_size <<"-" << a << ")";
558 attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
560 MESSAGE("Python function for attractor:" << std::endl << attractorFunctionStream.str());
563 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
564 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
567 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
569 FaceId2SizeMap[key] =attractorFunctionStream.str();
571 MESSAGE("Creating node on ("<<x0<<","<<y0<<","<<z0<<")");
572 FaceId2AttractorCoords[key] = coords;
576 /////////////////////////////////////////////////////////
578 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp, blsurf_session_t *bls)
580 int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
581 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
582 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
583 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
584 double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
585 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
586 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
587 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
588 bool _decimesh = BLSURFPlugin_Hypothesis::GetDefaultDecimesh();
589 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
592 MESSAGE("BLSURFPlugin_BLSURF::SetParameters");
593 _topology = (int) hyp->GetTopology();
594 _physicalMesh = (int) hyp->GetPhysicalMesh();
595 _phySize = hyp->GetPhySize();
596 _geometricMesh = (int) hyp->GetGeometricMesh();
597 _angleMeshS = hyp->GetAngleMeshS();
598 _angleMeshC = hyp->GetAngleMeshC();
599 _gradation = hyp->GetGradation();
600 _quadAllowed = hyp->GetQuadAllowed();
601 _decimesh = hyp->GetDecimesh();
602 _verb = hyp->GetVerbosity();
604 if ( hyp->GetPhyMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
605 blsurf_set_param(bls, "hphymin", to_string(hyp->GetPhyMin()).c_str());
606 if ( hyp->GetPhyMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
607 blsurf_set_param(bls, "hphymax", to_string(hyp->GetPhyMax()).c_str());
608 if ( hyp->GetGeoMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
609 blsurf_set_param(bls, "hgeomin", to_string(hyp->GetGeoMin()).c_str());
610 if ( hyp->GetGeoMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
611 blsurf_set_param(bls, "hgeomax", to_string(hyp->GetGeoMax()).c_str());
613 const BLSURFPlugin_Hypothesis::TOptionValues & opts = hyp->GetOptionValues();
614 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
615 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt )
616 if ( !opIt->second.empty() ) {
617 MESSAGE("blsurf_set_param(): " << opIt->first << " = " << opIt->second);
618 blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
622 MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
624 _smp_phy_size = _phySize;
625 blsurf_set_param(bls, "topo_points", _topology > 0 ? "1" : "0");
626 blsurf_set_param(bls, "topo_curves", _topology > 0 ? "1" : "0");
627 blsurf_set_param(bls, "topo_project", _topology > 0 ? "1" : "0");
628 blsurf_set_param(bls, "clean_boundary", _topology > 1 ? "1" : "0");
629 blsurf_set_param(bls, "close_boundary", _topology > 1 ? "1" : "0");
630 blsurf_set_param(bls, "hphy_flag", to_string(_physicalMesh).c_str());
631 // blsurf_set_param(bls, "hphy_flag", "2");
632 if ((to_string(_physicalMesh))=="2"){
633 TopoDS_Shape GeomShape;
634 TopAbs_ShapeEnum GeomType;
636 // Standard Size Maps
638 MESSAGE("Setting a Size Map");
639 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
640 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
641 for ( ; smIt != sizeMaps.end(); ++smIt ) {
642 if ( !smIt->second.empty() ) {
643 MESSAGE("blsurf_set_sizeMap(): " << smIt->first << " = " << smIt->second);
644 GeomShape = entryToShape(smIt->first);
645 GeomType = GeomShape.ShapeType();
646 MESSAGE("Geomtype is " << GeomType);
649 if (GeomType == TopAbs_COMPOUND){
650 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
652 if (it.Value().ShapeType() == TopAbs_FACE){
653 HasSizeMapOnFace = true;
654 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
655 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
658 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
659 // MESSAGE("Face with key " << key << " already in map");
661 FaceId2SizeMap[key] = smIt->second;
664 if (it.Value().ShapeType() == TopAbs_EDGE){
665 HasSizeMapOnEdge = true;
666 HasSizeMapOnFace = true;
667 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
668 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
671 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
672 // MESSAGE("Edge with key " << key << " already in map");
674 EdgeId2SizeMap[key] = smIt->second;
677 if (it.Value().ShapeType() == TopAbs_VERTEX){
678 HasSizeMapOnVertex = true;
679 HasSizeMapOnEdge = true;
680 HasSizeMapOnFace = true;
681 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
682 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
685 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
686 MESSAGE("Group of vertices with key " << key << " already in map");
688 MESSAGE("Group of vertices with key " << key << " has a size map: " << smIt->second);
689 VertexId2SizeMap[key] = smIt->second;
694 if (GeomType == TopAbs_FACE){
695 HasSizeMapOnFace = true;
696 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
697 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
700 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
701 // MESSAGE("Face with key " << key << " already in map");
703 FaceId2SizeMap[key] = smIt->second;
706 if (GeomType == TopAbs_EDGE){
707 HasSizeMapOnEdge = true;
708 HasSizeMapOnFace = true;
709 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
710 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
713 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
714 // MESSAGE("Edge with key " << key << " already in map");
716 EdgeId2SizeMap[key] = smIt->second;
719 if (GeomType == TopAbs_VERTEX){
720 HasSizeMapOnVertex = true;
721 HasSizeMapOnEdge = true;
722 HasSizeMapOnFace = true;
723 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
724 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
727 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
728 MESSAGE("Vertex with key " << key << " already in map");
730 MESSAGE("Vertex with key " << key << " has a size map: " << smIt->second);
731 VertexId2SizeMap[key] = smIt->second;
739 MESSAGE("Setting Attractors");
740 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
741 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
742 for ( ; atIt != attractors.end(); ++atIt ) {
743 if ( !atIt->second.empty() ) {
744 MESSAGE("blsurf_set_attractor(): " << atIt->first << " = " << atIt->second);
745 GeomShape = entryToShape(atIt->first);
746 GeomType = GeomShape.ShapeType();
748 if (GeomType == TopAbs_COMPOUND){
749 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
750 if (it.Value().ShapeType() == TopAbs_FACE){
751 HasSizeMapOnFace = true;
752 createAttractorOnFace(it.Value(), atIt->second);
757 if (GeomType == TopAbs_FACE){
758 HasSizeMapOnFace = true;
759 createAttractorOnFace(GeomShape, atIt->second);
762 if (GeomType == TopAbs_EDGE){
763 HasSizeMapOnEdge = true;
764 HasSizeMapOnFace = true;
765 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
767 if (GeomType == TopAbs_VERTEX){
768 HasSizeMapOnVertex = true;
769 HasSizeMapOnEdge = true;
770 HasSizeMapOnFace = true;
771 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
781 MESSAGE("Setting Enforced Vertices");
782 const BLSURFPlugin_Hypothesis::TEnforcedVertexMap enforcedVertexMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVertices(hyp);
783 BLSURFPlugin_Hypothesis::TEnforcedVertexMap::const_iterator enfIt = enforcedVertexMap.begin();
784 for ( ; enfIt != enforcedVertexMap.end(); ++enfIt ) {
785 if ( !enfIt->second.empty() ) {
786 GeomShape = entryToShape(enfIt->first);
787 GeomType = GeomShape.ShapeType();
789 if (GeomType == TopAbs_COMPOUND){
790 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
791 if (it.Value().ShapeType() == TopAbs_FACE){
792 HasSizeMapOnFace = true;
793 createEnforcedVertexOnFace(it.Value(), enfIt->second);
798 if (GeomType == TopAbs_FACE){
799 HasSizeMapOnFace = true;
800 createEnforcedVertexOnFace(GeomShape, enfIt->second);
805 // if (HasSizeMapOnFace){
806 // In all size map cases (hphy_flag = 2), at least map on face must be defined
807 MESSAGE("Setting Size Map on FACES ");
808 blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
811 if (HasSizeMapOnEdge){
812 MESSAGE("Setting Size Map on EDGES ");
813 blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
815 if (HasSizeMapOnVertex){
816 MESSAGE("Setting Size Map on VERTICES ");
817 blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
820 blsurf_set_param(bls, "hphydef", to_string(_phySize).c_str());
821 blsurf_set_param(bls, "hgeo_flag", to_string(_geometricMesh).c_str());
822 blsurf_set_param(bls, "relax_size", _decimesh ? "0": to_string(_geometricMesh).c_str());
823 blsurf_set_param(bls, "angle_meshs", to_string(_angleMeshS).c_str());
824 blsurf_set_param(bls, "angle_meshc", to_string(_angleMeshC).c_str());
825 blsurf_set_param(bls, "gradation", to_string(_gradation).c_str());
826 blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
827 blsurf_set_param(bls, "element", _quadAllowed ? "q1.0" : "p1");
828 blsurf_set_param(bls, "verb", to_string(_verb).c_str());
831 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
832 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
833 real *duu, real *duv, real *dvv, void *user_data);
834 status_t message_cb(message_t *msg, void *user_data);
835 status_t interrupt_cb(integer *interrupt_status, void *user_data);
837 //=============================================================================
841 //=============================================================================
843 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
845 MESSAGE("BLSURFPlugin_BLSURF::Compute");
847 // if (aShape.ShapeType() == TopAbs_COMPOUND) {
848 // MESSAGE(" the shape is a COMPOUND");
851 // MESSAGE(" the shape is UNKNOWN");
854 // Fix problem with locales
855 Kernel_Utils::Localizer loc;
857 /* create a distene context (generic object) */
858 status_t status = STATUS_ERROR;
860 context_t *ctx = context_new();
862 /* Set the message callback in the working context */
863 context_set_message_callback(ctx, message_cb, &_comment);
864 context_set_interrupt_callback(ctx, interrupt_cb, NULL);
866 /* create the CAD object we will work on. It is associated to the context ctx. */
867 cad_t *c = cad_new(ctx);
869 blsurf_session_t *bls = blsurf_session_new(ctx);
871 FacesWithSizeMap.Clear();
872 FaceId2SizeMap.clear();
873 EdgesWithSizeMap.Clear();
874 EdgeId2SizeMap.clear();
875 VerticesWithSizeMap.Clear();
876 VertexId2SizeMap.clear();
878 MESSAGE("BEGIN SetParameters");
879 SetParameters(_hypothesis, bls);
880 MESSAGE("END SetParameters");
882 /* Now fill the CAD object with data from your CAD
883 * environement. This is the most complex part of a successfull
887 // needed to prevent the opencascade memory managmement from freeing things
888 vector<Handle(Geom2d_Curve)> curves;
889 vector<Handle(Geom_Surface)> surfaces;
894 TopTools_IndexedMapOfShape fmap;
895 TopTools_IndexedMapOfShape emap;
896 TopTools_IndexedMapOfShape pmap;
899 FaceId2PythonSmp.clear();
901 EdgeId2PythonSmp.clear();
903 VertexId2PythonSmp.clear();
905 assert(Py_IsInitialized());
906 PyGILState_STATE gstate;
907 gstate = PyGILState_Ensure();
909 string theSizeMapStr;
911 /****************************************************************************************
913 *****************************************************************************************/
915 string bad_end = "return";
917 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
918 TopoDS_Face f=TopoDS::Face(face_iter.Current());
920 if (fmap.FindIndex(f) > 0)
925 surfaces.push_back(BRep_Tool::Surface(f));
927 /* create an object representing the face for blsurf */
928 /* where face_id is an integer identifying the face.
929 * surf_function is the function that defines the surface
930 * (For this face, it will be called by blsurf with your_face_object_ptr
933 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
935 /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
936 cad_face_set_tag(fce, iface);
938 /* Set face orientation (optional if you want a well oriented output mesh)*/
939 if(f.Orientation() != TopAbs_FORWARD){
940 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
942 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
945 if (HasSizeMapOnFace){
946 std::cout << "A size map is defined on a face" << std::endl;
948 faceKey = FacesWithSizeMap.FindIndex(f);
950 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
951 theSizeMapStr = FaceId2SizeMap[faceKey];
952 // check if function ends with "return"
953 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
955 // Expr To Python function, verification is performed at validation in GUI
956 PyObject * obj = NULL;
957 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
959 PyObject * func = NULL;
960 func = PyObject_GetAttrString(main_mod, "f");
961 FaceId2PythonSmp[iface]=func;
962 FaceId2SizeMap.erase(faceKey);
965 // Specific size map = Attractor
966 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
968 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
969 if (attractor_iter->first == faceKey) {
970 MESSAGE("Face indice: " << iface);
971 MESSAGE("Adding attractor");
973 double xyzCoords[3] = {attractor_iter->second[2],
974 attractor_iter->second[3],
975 attractor_iter->second[4]};
977 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
978 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
979 BRepClass_FaceClassifier scl(f,P,1e-7);
980 // scl.Perform() is bugged. The function was rewritten
982 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
983 TopAbs_State result = scl.State();
984 MESSAGE("Position of point on face: "<<result);
985 if ( result == TopAbs_OUT )
986 MESSAGE("Point is out of face: node is not created");
987 if ( result == TopAbs_UNKNOWN )
988 MESSAGE("Point position on face is unknown: node is not created");
989 if ( result == TopAbs_ON )
990 MESSAGE("Point is on border of face: node is not created");
991 if ( result == TopAbs_IN )
993 // Point is inside face and not on border
994 MESSAGE("Point is in face: node is created");
995 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
997 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << iatt);
998 cad_point_t* point_p = cad_point_new(fce, iatt, uvCoords);
999 cad_point_set_tag(point_p, iatt);
1001 FaceId2AttractorCoords.erase(faceKey);
1005 // Enforced Vertices
1006 faceKey = FacesWithEnforcedVertices.FindIndex(f);
1007 std::map<int,std::set<std::vector<double> > >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1008 if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1009 std::cout << "Some enforced vertices are defined" << std::endl;
1011 std::set<std::vector<double> > evl;
1012 // std::vector<double> ev;
1013 MESSAGE("Face indice: " << iface);
1014 MESSAGE("Adding enforced vertices");
1015 evl = evmIt->second;
1016 MESSAGE("Number of vertices to add: "<< evl.size());
1017 std::set< std::vector<double> >::const_iterator evlIt = evl.begin();
1018 for (; evlIt != evl.end(); ++evlIt) {
1020 // for (int i=0; i<evl.size() ; i++) {
1023 // double xyzCoords[3] = {ev[2], ev[3], ev[4]};
1024 double xyzCoords[3] = {evlIt->at(2), evlIt->at(3), evlIt->at(4)};
1025 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1026 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1027 BRepClass_FaceClassifier scl(f,P,1e-7);
1028 // scl.Perform() is bugged. The function was rewritten
1030 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1031 TopAbs_State result = scl.State();
1032 MESSAGE("Position of point on face: "<<result);
1033 if ( result == TopAbs_OUT )
1034 MESSAGE("Point is out of face: node is not created");
1035 if ( result == TopAbs_UNKNOWN )
1036 MESSAGE("Point position on face is unknown: node is not created");
1037 if ( result == TopAbs_ON )
1038 MESSAGE("Point is on border of face: node is not created");
1039 if ( result == TopAbs_IN )
1041 // Point is inside face and not on border
1042 MESSAGE("Point is in face: node is created");
1043 // double uvCoords[2] = {ev[0],ev[1]};
1044 double uvCoords[2] = {evlIt->at(0),evlIt->at(1)};
1046 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1047 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1048 cad_point_set_tag(point_p, ienf);
1051 FaceId2EnforcedVertexCoords.erase(faceKey);
1054 std::cout << "No enforced vertex defined" << std::endl;
1058 /****************************************************************************************
1060 now create the edges associated to this face
1061 *****************************************************************************************/
1063 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1064 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1065 int ic = emap.FindIndex(e);
1070 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1072 if (HasSizeMapOnEdge){
1073 edgeKey = EdgesWithSizeMap.FindIndex(e);
1074 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1075 MESSAGE("Sizemap defined on edge with index " << edgeKey);
1076 theSizeMapStr = EdgeId2SizeMap[edgeKey];
1077 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1079 // Expr To Python function, verification is performed at validation in GUI
1080 PyObject * obj = NULL;
1081 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1083 PyObject * func = NULL;
1084 func = PyObject_GetAttrString(main_mod, "f");
1085 EdgeId2PythonSmp[ic]=func;
1086 EdgeId2SizeMap.erase(edgeKey);
1090 /* attach the edge to the current blsurf face */
1091 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1093 /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1094 cad_edge_set_tag(edg, ic);
1096 /* by default, an edge does not necessalry appear in the resulting mesh,
1097 unless the following property is set :
1099 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1101 /* by default an edge is a boundary edge */
1102 if (e.Orientation() == TopAbs_INTERNAL)
1103 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1107 gp_Pnt2d e0 = curves.back()->Value(tmin);
1108 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1109 Standard_Real d1=0,d2=0;
1112 /****************************************************************************************
1114 *****************************************************************************************/
1116 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1117 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1121 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1124 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1126 *ip = pmap.FindIndex(v);
1130 //vertexKey = VerticesWithSizeMap.FindIndex(v);
1131 if (HasSizeMapOnVertex){
1132 vertexKey = VerticesWithSizeMap.FindIndex(v);
1133 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1134 theSizeMapStr = VertexId2SizeMap[vertexKey];
1135 //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1136 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1138 // Expr To Python function, verification is performed at validation in GUI
1139 PyObject * obj = NULL;
1140 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1142 PyObject * func = NULL;
1143 func = PyObject_GetAttrString(main_mod, "f");
1144 VertexId2PythonSmp[*ip]=func;
1145 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
1150 // should not happen
1151 MESSAGE("An edge does not have 2 extremities.");
1154 // This defines the curves extremity connectivity
1155 cad_edge_set_extremities(edg, ip1, ip2);
1156 /* set the tag (color) to the same value as the extremity id : */
1157 cad_edge_set_extremities_tag(edg, ip1, ip2);
1160 cad_edge_set_extremities(edg, ip2, ip1);
1161 cad_edge_set_extremities_tag(edg, ip2, ip1);
1168 PyGILState_Release(gstate);
1170 blsurf_data_set_cad(bls, c);
1172 std::cout << std::endl;
1173 std::cout << "Beginning of Surface Mesh generation" << std::endl;
1174 std::cout << std::endl;
1176 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1178 feclearexcept( FE_ALL_EXCEPT );
1179 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1185 status = blsurf_compute_mesh(bls);
1188 catch ( std::exception& exc ) {
1189 _comment += exc.what();
1191 catch (Standard_Failure& ex) {
1192 _comment += ex.DynamicType()->Name();
1193 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1195 _comment += ex.GetMessageString();
1199 if ( _comment.empty() )
1200 _comment = "Exception in blsurf_compute_mesh()";
1202 if ( status != STATUS_OK) {
1203 // There was an error while meshing
1204 blsurf_session_delete(bls);
1206 context_delete(ctx);
1208 return error(_comment);
1211 std::cout << std::endl;
1212 std::cout << "End of Surface Mesh generation" << std::endl;
1213 std::cout << std::endl;
1216 blsurf_data_get_mesh(bls, &msh);
1218 blsurf_session_delete(bls);
1220 context_delete(ctx);
1222 return error(_comment);
1226 /* retrieve mesh data (see distene/mesh.h) */
1227 integer nv, ne, nt, nq, vtx[4], tag;
1230 mesh_get_vertex_count(msh, &nv);
1231 mesh_get_edge_count(msh, &ne);
1232 mesh_get_triangle_count(msh, &nt);
1233 mesh_get_quadrangle_count(msh, &nq);
1236 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1237 SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1238 bool* tags = new bool[nv+1];
1240 /* enumerated vertices */
1241 for(int iv=1;iv<=nv;iv++) {
1242 mesh_get_vertex_coordinates(msh, iv, xyz);
1243 mesh_get_vertex_tag(msh, iv, &tag);
1244 // Issue 0020656. Use vertex coordinates
1245 if ( tag > 0 && tag <= pmap.Extent() ) {
1246 TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
1247 double tol = BRep_Tool::Tolerance( v );
1248 gp_Pnt p = BRep_Tool::Pnt( v );
1249 if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
1250 xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
1252 tag = 0; // enforced or attracted vertex
1254 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1255 // internal point are tagged to zero
1256 if(tag > 0 && tag <= pmap.Extent() ){
1257 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1264 /* enumerate edges */
1265 for(int it=1;it<=ne;it++) {
1266 mesh_get_edge_vertices(msh, it, vtx);
1267 SMDS_MeshEdge* edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1268 mesh_get_edge_tag(msh, it, &tag);
1271 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1272 tags[vtx[0]] = false;
1275 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1276 tags[vtx[1]] = false;
1278 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1282 /* enumerate triangles */
1283 for(int it=1;it<=nt;it++) {
1284 mesh_get_triangle_vertices(msh, it, vtx);
1285 SMDS_MeshFace* tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1286 mesh_get_triangle_tag(msh, it, &tag);
1287 meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1289 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1290 tags[vtx[0]] = false;
1293 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1294 tags[vtx[1]] = false;
1297 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1298 tags[vtx[2]] = false;
1302 /* enumerate quadrangles */
1303 for(int it=1;it<=nq;it++) {
1304 mesh_get_quadrangle_vertices(msh, it, vtx);
1305 SMDS_MeshFace* quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1306 mesh_get_quadrangle_tag(msh, it, &tag);
1307 meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1309 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1310 tags[vtx[0]] = false;
1313 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1314 tags[vtx[1]] = false;
1317 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1318 tags[vtx[2]] = false;
1321 meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1322 tags[vtx[3]] = false;
1328 /* release the mesh object */
1329 blsurf_data_regain_mesh(bls, msh);
1331 /* clean up everything */
1332 blsurf_session_delete(bls);
1335 context_delete(ctx);
1337 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1339 if ( oldFEFlags > 0 )
1340 feenableexcept( oldFEFlags );
1341 feclearexcept( FE_ALL_EXCEPT );
1345 std::cout << "FacesWithSizeMap" << std::endl;
1346 FacesWithSizeMap.Statistics(std::cout);
1347 std::cout << "EdgesWithSizeMap" << std::endl;
1348 EdgesWithSizeMap.Statistics(std::cout);
1349 std::cout << "VerticesWithSizeMap" << std::endl;
1350 VerticesWithSizeMap.Statistics(std::cout);
1351 std::cout << "FacesWithEnforcedVertices" << std::endl;
1352 FacesWithEnforcedVertices.Statistics(std::cout);
1358 //=============================================================================
1362 //=============================================================================
1364 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1365 const TopoDS_Edge edge = TopoDS::Edge(ed);
1367 gp_Pnt pnt(node->X(), node->Y(), node->Z());
1369 Standard_Real p0 = 0.0;
1370 Standard_Real p1 = 1.0;
1371 TopLoc_Location loc;
1372 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1374 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1375 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1378 if ( proj.NbPoints() > 0 )
1380 pa = (double)proj.LowerDistanceParameter();
1381 // Issue 0020656. Move node if it is too far from edge
1382 gp_Pnt curve_pnt = curve->Value( pa );
1383 double dist2 = pnt.SquareDistance( curve_pnt );
1384 double tol = BRep_Tool::Tolerance( edge );
1385 if ( 1e-12 < dist2 && dist2 <= 2*tol*tol ) // large enough and within tolerance
1387 curve_pnt.Transform( loc );
1388 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1391 // GProp_GProps LProps;
1392 // BRepGProp::LinearProperties(ed, LProps);
1393 // double lg = (double)LProps.Mass();
1395 meshDS->SetNodeOnEdge(node, edge, pa);
1398 //=============================================================================
1402 //=============================================================================
1404 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
1409 //=============================================================================
1413 //=============================================================================
1415 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
1420 //=============================================================================
1424 //=============================================================================
1426 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
1428 return hyp.SaveTo( save );
1431 //=============================================================================
1435 //=============================================================================
1437 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
1439 return hyp.LoadFrom( load );
1442 /* Curve definition function See cad_curv_t in file distene/cad.h for
1444 * NOTE : if when your CAD systems evaluates second
1445 * order derivatives it also computes first order derivatives and
1446 * function evaluation, you can optimize this example by making only
1447 * one CAD call and filling the necessary uv, dt, dtt arrays.
1449 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
1451 /* t is given. It contains the t (time) 1D parametric coordintaes
1452 of the point PreCAD/BLSurf is querying on the curve */
1454 /* user_data identifies the edge PreCAD/BLSurf is querying
1455 * (see cad_edge_new later in this example) */
1456 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
1459 /* BLSurf is querying the function evaluation */
1462 uv[0]=P.X(); uv[1]=P.Y();
1466 /* query for the first order derivatives */
1469 dt[0]=V1.X(); dt[1]=V1.Y();
1473 /* query for the second order derivatives */
1476 dtt[0]=V2.X(); dtt[1]=V2.Y();
1482 /* Surface definition function.
1483 * See cad_surf_t in file distene/cad.h for more information.
1484 * NOTE : if when your CAD systems evaluates second order derivatives it also
1485 * computes first order derivatives and function evaluation, you can optimize
1486 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
1489 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1490 real *duu, real *duv, real *dvv, void *user_data)
1492 /* uv[2] is given. It contains the u,v coordinates of the point
1493 * PreCAD/BLSurf is querying on the surface */
1495 /* user_data identifies the face PreCAD/BLSurf is querying (see
1496 * cad_face_new later in this example)*/
1497 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
1501 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
1502 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
1509 geometry->D1(uv[0],uv[1],P,D1U,D1V);
1510 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
1511 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
1514 if(duu && duv && dvv){
1518 gp_Vec D2U,D2V,D2UV;
1520 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
1521 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
1522 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
1523 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
1530 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
1533 if (my_u_min > uv[0]) {
1536 if (my_v_min > uv[1]) {
1539 if (my_u_max < uv[0]) {
1542 if (my_v_max < uv[1]) {
1547 if (FaceId2PythonSmp.count(face_id) != 0){
1548 PyObject * pyresult = NULL;
1549 PyObject* new_stderr = NULL;
1550 assert(Py_IsInitialized());
1551 PyGILState_STATE gstate;
1552 gstate = PyGILState_Ensure();
1553 pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],"(f,f)",uv[0],uv[1]);
1555 if ( pyresult == NULL){
1557 string err_description="";
1558 new_stderr = newPyStdOut(err_description);
1559 PySys_SetObject("stderr", new_stderr);
1561 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1562 Py_DECREF(new_stderr);
1563 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
1564 result = *((double*)user_data);
1567 result = PyFloat_AsDouble(pyresult);
1568 Py_DECREF(pyresult);
1571 //MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1572 PyGILState_Release(gstate);
1575 *size = *((double*)user_data);
1580 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
1582 if (EdgeId2PythonSmp.count(edge_id) != 0){
1583 PyObject * pyresult = NULL;
1584 PyObject* new_stderr = NULL;
1585 assert(Py_IsInitialized());
1586 PyGILState_STATE gstate;
1587 gstate = PyGILState_Ensure();
1588 pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],"(f)",t);
1590 if ( pyresult == NULL){
1592 string err_description="";
1593 new_stderr = newPyStdOut(err_description);
1594 PySys_SetObject("stderr", new_stderr);
1596 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1597 Py_DECREF(new_stderr);
1598 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
1599 result = *((double*)user_data);
1602 result = PyFloat_AsDouble(pyresult);
1603 Py_DECREF(pyresult);
1606 PyGILState_Release(gstate);
1609 *size = *((double*)user_data);
1614 status_t size_on_vertex(integer point_id, real *size, void *user_data)
1616 if (VertexId2PythonSmp.count(point_id) != 0){
1617 PyObject * pyresult = NULL;
1618 PyObject* new_stderr = NULL;
1619 assert(Py_IsInitialized());
1620 PyGILState_STATE gstate;
1621 gstate = PyGILState_Ensure();
1622 pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],"");
1624 if ( pyresult == NULL){
1626 string err_description="";
1627 new_stderr = newPyStdOut(err_description);
1628 PySys_SetObject("stderr", new_stderr);
1630 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1631 Py_DECREF(new_stderr);
1632 MESSAGE("Can't evaluate f()" << " error is " << err_description);
1633 result = *((double*)user_data);
1636 result = PyFloat_AsDouble(pyresult);
1637 Py_DECREF(pyresult);
1640 PyGILState_Release(gstate);
1643 *size = *((double*)user_data);
1649 * The following function will be called for PreCAD/BLSurf message
1650 * printing. See context_set_message_callback (later in this
1651 * template) for how to set user_data.
1653 status_t message_cb(message_t *msg, void *user_data)
1655 integer errnumber = 0;
1657 message_get_number(msg, &errnumber);
1658 message_get_description(msg, &desc);
1659 if ( errnumber < 0 ) {
1660 string * error = (string*)user_data;
1661 // if ( !error->empty() )
1663 // remove ^A from the tail
1664 int len = strlen( desc );
1665 while (len > 0 && desc[len-1] != '\n')
1667 error->append( desc, len );
1670 std::cout << desc << std::endl;
1675 /* This is the interrupt callback. PreCAD/BLSurf will call this
1676 * function regularily. See the file distene/interrupt.h
1678 status_t interrupt_cb(integer *interrupt_status, void *user_data)
1680 integer you_want_to_continue = 1;
1682 if(you_want_to_continue)
1683 *interrupt_status = INTERRUPT_CONTINUE;
1684 else /* you want to stop BLSurf */
1685 *interrupt_status = INTERRUPT_STOP;
1690 //=============================================================================
1694 //=============================================================================
1695 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
1696 const TopoDS_Shape& aShape,
1697 MapShapeNbElems& aResMap)
1699 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
1700 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
1701 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
1702 //double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
1703 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
1704 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
1706 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
1707 _phySize = _hypothesis->GetPhySize();
1708 //_geometricMesh = (int) hyp->GetGeometricMesh();
1709 //_angleMeshS = hyp->GetAngleMeshS();
1710 _angleMeshC = _hypothesis->GetAngleMeshC();
1711 _quadAllowed = _hypothesis->GetQuadAllowed();
1714 bool IsQuadratic = false;
1719 TopTools_DataMapOfShapeInteger EdgesMap;
1720 double fullLen = 0.0;
1721 double fullNbSeg = 0;
1722 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1723 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
1724 if( EdgesMap.IsBound(E) )
1726 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
1727 double aLen = SMESH_Algo::EdgeLength(E);
1730 if(_physicalMesh==1) {
1731 nb1d = (int)( aLen/_phySize + 1 );
1736 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
1737 double fullAng = 0.0;
1738 double dp = (l-f)/200;
1743 for(int j=2; j<=200; j++) {
1746 fullAng += fabs(V1.Angle(V2));
1750 nb1d = (int)( fullAng/_angleMeshC + 1 );
1753 std::vector<int> aVec(SMDSEntity_Last);
1754 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1755 if( IsQuadratic > 0 ) {
1756 aVec[SMDSEntity_Node] = 2*nb1d - 1;
1757 aVec[SMDSEntity_Quad_Edge] = nb1d;
1760 aVec[SMDSEntity_Node] = nb1d - 1;
1761 aVec[SMDSEntity_Edge] = nb1d;
1763 aResMap.insert(std::make_pair(sm,aVec));
1764 EdgesMap.Bind(E,nb1d);
1766 double ELen = fullLen/fullNbSeg;
1770 // try to evaluate as in MEFISTO
1771 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1772 TopoDS_Face F = TopoDS::Face( exp.Current() );
1773 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1775 BRepGProp::SurfaceProperties(F,G);
1776 double anArea = G.Mass();
1778 std::vector<int> nb1dVec;
1779 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
1780 int nbSeg = EdgesMap.Find(exp1.Current());
1782 nb1dVec.push_back( nbSeg );
1785 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
1786 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
1789 if ( nb1dVec.size() == 4 ) // quadrangle geom face
1791 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
1793 nbNodes = (n1 + 1) * (n2 + 1);
1798 nbTria = nbQuad = nbTria / 3 + 1;
1801 std::vector<int> aVec(SMDSEntity_Last,0);
1803 int nb1d_in = (nbTria*3 - nb1d) / 2;
1804 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
1805 aVec[SMDSEntity_Quad_Triangle] = nbTria;
1806 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
1809 aVec[SMDSEntity_Node] = nbNodes;
1810 aVec[SMDSEntity_Triangle] = nbTria;
1811 aVec[SMDSEntity_Quadrangle] = nbQuad;
1813 aResMap.insert(std::make_pair(sm,aVec));
1820 BRepGProp::VolumeProperties(aShape,G);
1821 double aVolume = G.Mass();
1822 double tetrVol = 0.1179*ELen*ELen*ELen;
1823 int nbVols = int(aVolume/tetrVol);
1824 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
1825 std::vector<int> aVec(SMDSEntity_Last);
1826 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1828 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
1829 aVec[SMDSEntity_Quad_Tetra] = nbVols;
1832 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
1833 aVec[SMDSEntity_Tetra] = nbVols;
1835 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1836 aResMap.insert(std::make_pair(sm,aVec));
1841 //=============================================================================
1843 * Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
1844 * Following line was added:
1845 * myExtrem.Perform(P);
1847 //=============================================================================
1848 void BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
1849 const TopoDS_Face& face,
1851 const Standard_Real Tol)
1853 //-- Voir BRepExtrema_ExtPF.cxx
1854 BRepAdaptor_Surface Surf(face);
1855 Standard_Real U1, U2, V1, V2;
1856 BRepTools::UVBounds(face, U1, U2, V1, V2);
1857 Extrema_ExtPS myExtrem;
1858 myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
1859 myExtrem.Perform(P);
1860 //----------------------------------------------------------
1861 //-- On cherche le point le plus proche , PUIS
1862 //-- On le classifie.
1863 Standard_Integer nbv = 0; // xpu
1864 Standard_Real MaxDist = RealLast();
1865 Standard_Integer indice = 0;
1866 if(myExtrem.IsDone()) {
1867 nbv = myExtrem.NbExt();
1868 for (Standard_Integer i = 1; i <= nbv; i++) {
1869 Standard_Real d = myExtrem.Value(i);
1879 Standard_Real U1,U2;
1880 myExtrem.Point(indice).Parameter(U1, U2);
1881 Puv.SetCoord(U1, U2);
1882 fc->Perform(face, Puv, Tol);
1885 fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
1886 //-- le resultat est TopAbs_OUT;