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 coords.push_back(u0);
442 coords.push_back(v0);
443 coords.push_back(x0);
444 coords.push_back(y0);
445 coords.push_back(z0);
448 if (! FacesWithEnforcedVertices.Contains(TopoDS::Face(GeomShape))) {
449 key = FacesWithEnforcedVertices.Add(TopoDS::Face(GeomShape));
452 key = FacesWithEnforcedVertices.FindIndex(TopoDS::Face(GeomShape));
455 // If a node is already created by an attractor, do not create enforced vertex
456 int attractorKey = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
457 bool sameAttractor = false;
458 if (attractorKey >= 0)
459 if (FaceId2AttractorCoords.count(attractorKey) > 0)
460 if (FaceId2AttractorCoords[attractorKey] == coords)
461 sameAttractor = true;
463 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
464 MESSAGE("Map of enf. vertex has key " << key)
465 MESSAGE("Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
467 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
469 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
470 MESSAGE("New Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
473 MESSAGE("Map of enf. vertex has not key " << key << ": creating it")
474 if (! sameAttractor) {
475 std::set< std::vector<double> > ens;
477 FaceId2EnforcedVertexCoords[key] = ens;
480 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
485 /////////////////////////////////////////////////////////
486 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction)
488 MESSAGE("Attractor function: "<< AttractorFunction);
489 double xa, ya, za; // Coordinates of attractor point
490 double a, b; // Attractor parameter
491 bool createNode=false; // To create a node on attractor projection
493 const char *sep = ";";
494 // atIt->second has the following pattern:
495 // ATTRACTOR(xa;ya;za;a;b)
497 // xa;ya;za : coordinates of attractor
498 // a : desired size on attractor
499 // b : distance of influence of attractor
501 // We search the parameters in the string
503 pos1 = AttractorFunction.find(sep);
504 if (pos1!=string::npos)
505 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
507 pos2 = AttractorFunction.find(sep, pos1+1);
508 if (pos2!=string::npos) {
509 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
513 pos2 = AttractorFunction.find(sep, pos1+1);
514 if (pos2!=string::npos) {
515 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
519 pos2 = AttractorFunction.find(sep, pos1+1);
520 if (pos2!=string::npos) {
521 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
525 pos2 = AttractorFunction.find(sep, pos1+1);
526 if (pos2!=string::npos) {
527 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
531 pos2 = AttractorFunction.find(")");
532 if (pos2!=string::npos) {
533 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
534 MESSAGE("createNode: " << createNodeStr);
535 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
538 // Get the (u,v) values of the attractor on the face
539 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_XYZ(xa,ya,za));
540 gp_XY uvPoint = myPoint.uv;
541 gp_XYZ xyzPoint = myPoint.xyz;
542 Standard_Real u0 = uvPoint.X();
543 Standard_Real v0 = uvPoint.Y();
544 Standard_Real x0 = xyzPoint.X();
545 Standard_Real y0 = xyzPoint.Y();
546 Standard_Real z0 = xyzPoint.Z();
547 std::vector<double> coords;
548 coords.push_back(u0);
549 coords.push_back(v0);
550 coords.push_back(x0);
551 coords.push_back(y0);
552 coords.push_back(z0);
553 // We construct the python function
554 ostringstream attractorFunctionStream;
555 attractorFunctionStream << "def f(u,v): return ";
556 attractorFunctionStream << _smp_phy_size << "-(" << _smp_phy_size <<"-" << a << ")";
557 attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
559 MESSAGE("Python function for attractor:" << std::endl << attractorFunctionStream.str());
562 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
563 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
566 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
568 FaceId2SizeMap[key] =attractorFunctionStream.str();
570 MESSAGE("Creating node on ("<<x0<<","<<y0<<","<<z0<<")");
571 FaceId2AttractorCoords[key] = coords;
575 /////////////////////////////////////////////////////////
577 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp, blsurf_session_t *bls)
579 int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
580 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
581 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
582 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
583 double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
584 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
585 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
586 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
587 bool _decimesh = BLSURFPlugin_Hypothesis::GetDefaultDecimesh();
588 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
591 MESSAGE("BLSURFPlugin_BLSURF::SetParameters");
592 _topology = (int) hyp->GetTopology();
593 _physicalMesh = (int) hyp->GetPhysicalMesh();
594 _phySize = hyp->GetPhySize();
595 _geometricMesh = (int) hyp->GetGeometricMesh();
596 _angleMeshS = hyp->GetAngleMeshS();
597 _angleMeshC = hyp->GetAngleMeshC();
598 _gradation = hyp->GetGradation();
599 _quadAllowed = hyp->GetQuadAllowed();
600 _decimesh = hyp->GetDecimesh();
601 _verb = hyp->GetVerbosity();
603 if ( hyp->GetPhyMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
604 blsurf_set_param(bls, "hphymin", to_string(hyp->GetPhyMin()).c_str());
605 if ( hyp->GetPhyMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
606 blsurf_set_param(bls, "hphymax", to_string(hyp->GetPhyMax()).c_str());
607 if ( hyp->GetGeoMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
608 blsurf_set_param(bls, "hgeomin", to_string(hyp->GetGeoMin()).c_str());
609 if ( hyp->GetGeoMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
610 blsurf_set_param(bls, "hgeomax", to_string(hyp->GetGeoMax()).c_str());
612 const BLSURFPlugin_Hypothesis::TOptionValues & opts = hyp->GetOptionValues();
613 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
614 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt )
615 if ( !opIt->second.empty() ) {
616 MESSAGE("blsurf_set_param(): " << opIt->first << " = " << opIt->second);
617 blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
621 MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
623 _smp_phy_size = _phySize;
624 blsurf_set_param(bls, "topo_points", _topology > 0 ? "1" : "0");
625 blsurf_set_param(bls, "topo_curves", _topology > 0 ? "1" : "0");
626 blsurf_set_param(bls, "topo_project", _topology > 0 ? "1" : "0");
627 blsurf_set_param(bls, "clean_boundary", _topology > 1 ? "1" : "0");
628 blsurf_set_param(bls, "close_boundary", _topology > 1 ? "1" : "0");
629 blsurf_set_param(bls, "hphy_flag", to_string(_physicalMesh).c_str());
630 // blsurf_set_param(bls, "hphy_flag", "2");
631 if ((to_string(_physicalMesh))=="2"){
632 TopoDS_Shape GeomShape;
633 TopAbs_ShapeEnum GeomType;
635 // Standard Size Maps
637 MESSAGE("Setting a Size Map");
638 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
639 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
640 for ( ; smIt != sizeMaps.end(); ++smIt ) {
641 if ( !smIt->second.empty() ) {
642 MESSAGE("blsurf_set_sizeMap(): " << smIt->first << " = " << smIt->second);
643 GeomShape = entryToShape(smIt->first);
644 GeomType = GeomShape.ShapeType();
645 MESSAGE("Geomtype is " << GeomType);
648 if (GeomType == TopAbs_COMPOUND){
649 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
651 if (it.Value().ShapeType() == TopAbs_FACE){
652 HasSizeMapOnFace = true;
653 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
654 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
657 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
658 // MESSAGE("Face with key " << key << " already in map");
660 FaceId2SizeMap[key] = smIt->second;
663 if (it.Value().ShapeType() == TopAbs_EDGE){
664 HasSizeMapOnEdge = true;
665 HasSizeMapOnFace = true;
666 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
667 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
670 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
671 // MESSAGE("Edge with key " << key << " already in map");
673 EdgeId2SizeMap[key] = smIt->second;
676 if (it.Value().ShapeType() == TopAbs_VERTEX){
677 HasSizeMapOnVertex = true;
678 HasSizeMapOnEdge = true;
679 HasSizeMapOnFace = true;
680 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
681 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
684 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
685 MESSAGE("Group of vertices with key " << key << " already in map");
687 MESSAGE("Group of vertices with key " << key << " has a size map: " << smIt->second);
688 VertexId2SizeMap[key] = smIt->second;
693 if (GeomType == TopAbs_FACE){
694 HasSizeMapOnFace = true;
695 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
696 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
699 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
700 // MESSAGE("Face with key " << key << " already in map");
702 FaceId2SizeMap[key] = smIt->second;
705 if (GeomType == TopAbs_EDGE){
706 HasSizeMapOnEdge = true;
707 HasSizeMapOnFace = true;
708 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
709 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
712 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
713 // MESSAGE("Edge with key " << key << " already in map");
715 EdgeId2SizeMap[key] = smIt->second;
718 if (GeomType == TopAbs_VERTEX){
719 HasSizeMapOnVertex = true;
720 HasSizeMapOnEdge = true;
721 HasSizeMapOnFace = true;
722 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
723 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
726 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
727 MESSAGE("Vertex with key " << key << " already in map");
729 MESSAGE("Vertex with key " << key << " has a size map: " << smIt->second);
730 VertexId2SizeMap[key] = smIt->second;
738 MESSAGE("Setting Attractors");
739 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
740 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
741 for ( ; atIt != attractors.end(); ++atIt ) {
742 if ( !atIt->second.empty() ) {
743 MESSAGE("blsurf_set_attractor(): " << atIt->first << " = " << atIt->second);
744 GeomShape = entryToShape(atIt->first);
745 GeomType = GeomShape.ShapeType();
747 if (GeomType == TopAbs_COMPOUND){
748 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
749 if (it.Value().ShapeType() == TopAbs_FACE){
750 HasSizeMapOnFace = true;
751 createAttractorOnFace(it.Value(), atIt->second);
756 if (GeomType == TopAbs_FACE){
757 HasSizeMapOnFace = true;
758 createAttractorOnFace(GeomShape, atIt->second);
761 if (GeomType == TopAbs_EDGE){
762 HasSizeMapOnEdge = true;
763 HasSizeMapOnFace = true;
764 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
766 if (GeomType == TopAbs_VERTEX){
767 HasSizeMapOnVertex = true;
768 HasSizeMapOnEdge = true;
769 HasSizeMapOnFace = true;
770 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
780 MESSAGE("Setting Enforced Vertices");
781 const BLSURFPlugin_Hypothesis::TEnforcedVertexMap enforcedVertexMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVertices(hyp);
782 BLSURFPlugin_Hypothesis::TEnforcedVertexMap::const_iterator enfIt = enforcedVertexMap.begin();
783 for ( ; enfIt != enforcedVertexMap.end(); ++enfIt ) {
784 if ( !enfIt->second.empty() ) {
785 GeomShape = entryToShape(enfIt->first);
786 GeomType = GeomShape.ShapeType();
788 if (GeomType == TopAbs_COMPOUND){
789 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
790 if (it.Value().ShapeType() == TopAbs_FACE){
791 HasSizeMapOnFace = true;
792 createEnforcedVertexOnFace(it.Value(), enfIt->second);
797 if (GeomType == TopAbs_FACE){
798 HasSizeMapOnFace = true;
799 createEnforcedVertexOnFace(GeomShape, enfIt->second);
804 // if (HasSizeMapOnFace){
805 // In all size map cases (hphy_flag = 2), at least map on face must be defined
806 MESSAGE("Setting Size Map on FACES ");
807 blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
810 if (HasSizeMapOnEdge){
811 MESSAGE("Setting Size Map on EDGES ");
812 blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
814 if (HasSizeMapOnVertex){
815 MESSAGE("Setting Size Map on VERTICES ");
816 blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
819 blsurf_set_param(bls, "hphydef", to_string(_phySize).c_str());
820 blsurf_set_param(bls, "hgeo_flag", to_string(_geometricMesh).c_str());
821 blsurf_set_param(bls, "relax_size", _decimesh ? "0": to_string(_geometricMesh).c_str());
822 blsurf_set_param(bls, "angle_meshs", to_string(_angleMeshS).c_str());
823 blsurf_set_param(bls, "angle_meshc", to_string(_angleMeshC).c_str());
824 blsurf_set_param(bls, "gradation", to_string(_gradation).c_str());
825 blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
826 blsurf_set_param(bls, "element", _quadAllowed ? "q1.0" : "p1");
827 blsurf_set_param(bls, "verb", to_string(_verb).c_str());
830 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
831 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
832 real *duu, real *duv, real *dvv, void *user_data);
833 status_t message_cb(message_t *msg, void *user_data);
834 status_t interrupt_cb(integer *interrupt_status, void *user_data);
836 //=============================================================================
840 //=============================================================================
842 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
844 MESSAGE("BLSURFPlugin_BLSURF::Compute");
846 // if (aShape.ShapeType() == TopAbs_COMPOUND) {
847 // MESSAGE(" the shape is a COMPOUND");
850 // MESSAGE(" the shape is UNKNOWN");
853 // Fix problem with locales
854 Kernel_Utils::Localizer loc;
856 /* create a distene context (generic object) */
857 status_t status = STATUS_ERROR;
859 context_t *ctx = context_new();
861 /* Set the message callback in the working context */
862 context_set_message_callback(ctx, message_cb, &_comment);
863 context_set_interrupt_callback(ctx, interrupt_cb, NULL);
865 /* create the CAD object we will work on. It is associated to the context ctx. */
866 cad_t *c = cad_new(ctx);
868 blsurf_session_t *bls = blsurf_session_new(ctx);
870 FacesWithSizeMap.Clear();
871 FaceId2SizeMap.clear();
872 EdgesWithSizeMap.Clear();
873 EdgeId2SizeMap.clear();
874 VerticesWithSizeMap.Clear();
875 VertexId2SizeMap.clear();
877 MESSAGE("BEGIN SetParameters");
878 SetParameters(_hypothesis, bls);
879 MESSAGE("END SetParameters");
881 /* Now fill the CAD object with data from your CAD
882 * environement. This is the most complex part of a successfull
886 // needed to prevent the opencascade memory managmement from freeing things
887 vector<Handle(Geom2d_Curve)> curves;
888 vector<Handle(Geom_Surface)> surfaces;
893 TopTools_IndexedMapOfShape fmap;
894 TopTools_IndexedMapOfShape emap;
895 TopTools_IndexedMapOfShape pmap;
898 FaceId2PythonSmp.clear();
900 EdgeId2PythonSmp.clear();
902 VertexId2PythonSmp.clear();
904 assert(Py_IsInitialized());
905 PyGILState_STATE gstate;
906 gstate = PyGILState_Ensure();
908 string theSizeMapStr;
910 /****************************************************************************************
912 *****************************************************************************************/
914 string bad_end = "return";
916 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
917 TopoDS_Face f=TopoDS::Face(face_iter.Current());
919 if (fmap.FindIndex(f) > 0)
924 surfaces.push_back(BRep_Tool::Surface(f));
926 /* create an object representing the face for blsurf */
927 /* where face_id is an integer identifying the face.
928 * surf_function is the function that defines the surface
929 * (For this face, it will be called by blsurf with your_face_object_ptr
932 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
934 /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
935 cad_face_set_tag(fce, iface);
937 /* Set face orientation (optional if you want a well oriented output mesh)*/
938 if(f.Orientation() != TopAbs_FORWARD){
939 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
941 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
944 if (HasSizeMapOnFace){
945 std::cout << "A size map is defined on a face" << std::endl;
947 faceKey = FacesWithSizeMap.FindIndex(f);
949 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
950 theSizeMapStr = FaceId2SizeMap[faceKey];
951 // check if function ends with "return"
952 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
954 // Expr To Python function, verification is performed at validation in GUI
955 PyObject * obj = NULL;
956 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
958 PyObject * func = NULL;
959 func = PyObject_GetAttrString(main_mod, "f");
960 FaceId2PythonSmp[iface]=func;
961 FaceId2SizeMap.erase(faceKey);
964 // Specific size map = Attractor
965 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
967 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
968 if (attractor_iter->first == faceKey) {
969 MESSAGE("Face indice: " << iface);
970 MESSAGE("Adding attractor");
972 double xyzCoords[3] = {attractor_iter->second[2],
973 attractor_iter->second[3],
974 attractor_iter->second[4]};
976 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
977 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
978 BRepClass_FaceClassifier scl(f,P,1e-7);
979 // scl.Perform() is bugged. The function was rewritten
981 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
982 TopAbs_State result = scl.State();
983 MESSAGE("Position of point on face: "<<result);
984 if ( result == TopAbs_OUT )
985 MESSAGE("Point is out of face: node is not created");
986 if ( result == TopAbs_UNKNOWN )
987 MESSAGE("Point position on face is unknown: node is not created");
988 if ( result == TopAbs_ON )
989 MESSAGE("Point is on border of face: node is not created");
990 if ( result == TopAbs_IN )
992 // Point is inside face and not on border
993 MESSAGE("Point is in face: node is created");
994 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
996 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << iatt);
997 cad_point_t* point_p = cad_point_new(fce, iatt, uvCoords);
998 cad_point_set_tag(point_p, iatt);
1000 FaceId2AttractorCoords.erase(faceKey);
1004 // Enforced Vertices
1005 faceKey = FacesWithEnforcedVertices.FindIndex(f);
1006 std::map<int,std::set<std::vector<double> > >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1007 if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1008 std::cout << "Some enforced vertices are defined" << std::endl;
1010 std::set<std::vector<double> > evl;
1011 // std::vector<double> ev;
1012 MESSAGE("Face indice: " << iface);
1013 MESSAGE("Adding enforced vertices");
1014 evl = evmIt->second;
1015 MESSAGE("Number of vertices to add: "<< evl.size())
1016 std::set< std::vector<double> >::const_iterator evlIt = evl.begin();
1017 for (; evlIt != evl.end(); ++evlIt) {
1019 // for (int i=0; i<evl.size() ; i++) {
1022 // double xyzCoords[3] = {ev[2], ev[3], ev[4]};
1023 double xyzCoords[3] = {evlIt->at(0), evlIt->at(3), evlIt->at(4)};
1024 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1025 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1026 BRepClass_FaceClassifier scl(f,P,1e-7);
1027 // scl.Perform() is bugged. The function was rewritten
1029 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1030 TopAbs_State result = scl.State();
1031 MESSAGE("Position of point on face: "<<result);
1032 if ( result == TopAbs_OUT )
1033 MESSAGE("Point is out of face: node is not created");
1034 if ( result == TopAbs_UNKNOWN )
1035 MESSAGE("Point position on face is unknown: node is not created");
1036 if ( result == TopAbs_ON )
1037 MESSAGE("Point is on border of face: node is not created");
1038 if ( result == TopAbs_IN )
1040 // Point is inside face and not on border
1041 MESSAGE("Point is in face: node is created");
1042 // double uvCoords[2] = {ev[0],ev[1]};
1043 double uvCoords[2] = {evlIt->at(0),evlIt->at(1)};
1045 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1046 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1047 cad_point_set_tag(point_p, ienf);
1050 FaceId2EnforcedVertexCoords.erase(faceKey);
1053 std::cout << "No enforced vertex defined" << std::endl;
1057 /****************************************************************************************
1059 now create the edges associated to this face
1060 *****************************************************************************************/
1062 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1063 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1064 int ic = emap.FindIndex(e);
1069 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1071 if (HasSizeMapOnEdge){
1072 edgeKey = EdgesWithSizeMap.FindIndex(e);
1073 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1074 MESSAGE("Sizemap defined on edge with index " << edgeKey);
1075 theSizeMapStr = EdgeId2SizeMap[edgeKey];
1076 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1078 // Expr To Python function, verification is performed at validation in GUI
1079 PyObject * obj = NULL;
1080 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1082 PyObject * func = NULL;
1083 func = PyObject_GetAttrString(main_mod, "f");
1084 EdgeId2PythonSmp[ic]=func;
1085 EdgeId2SizeMap.erase(edgeKey);
1089 /* attach the edge to the current blsurf face */
1090 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1092 /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1093 cad_edge_set_tag(edg, ic);
1095 /* by default, an edge does not necessalry appear in the resulting mesh,
1096 unless the following property is set :
1098 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1100 /* by default an edge is a boundary edge */
1101 if (e.Orientation() == TopAbs_INTERNAL)
1102 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1106 gp_Pnt2d e0 = curves.back()->Value(tmin);
1107 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1108 Standard_Real d1=0,d2=0;
1111 /****************************************************************************************
1113 *****************************************************************************************/
1115 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1116 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1120 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1123 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1125 *ip = pmap.FindIndex(v);
1129 //vertexKey = VerticesWithSizeMap.FindIndex(v);
1130 if (HasSizeMapOnVertex){
1131 vertexKey = VerticesWithSizeMap.FindIndex(v);
1132 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1133 theSizeMapStr = VertexId2SizeMap[vertexKey];
1134 //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1135 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1137 // Expr To Python function, verification is performed at validation in GUI
1138 PyObject * obj = NULL;
1139 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1141 PyObject * func = NULL;
1142 func = PyObject_GetAttrString(main_mod, "f");
1143 VertexId2PythonSmp[*ip]=func;
1144 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
1149 // should not happen
1150 MESSAGE("An edge does not have 2 extremities.");
1153 // This defines the curves extremity connectivity
1154 cad_edge_set_extremities(edg, ip1, ip2);
1155 /* set the tag (color) to the same value as the extremity id : */
1156 cad_edge_set_extremities_tag(edg, ip1, ip2);
1159 cad_edge_set_extremities(edg, ip2, ip1);
1160 cad_edge_set_extremities_tag(edg, ip2, ip1);
1167 PyGILState_Release(gstate);
1169 blsurf_data_set_cad(bls, c);
1171 std::cout << std::endl;
1172 std::cout << "Beginning of Surface Mesh generation" << std::endl;
1173 std::cout << std::endl;
1175 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1177 feclearexcept( FE_ALL_EXCEPT );
1178 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1184 status = blsurf_compute_mesh(bls);
1187 catch ( std::exception& exc ) {
1188 _comment += exc.what();
1190 catch (Standard_Failure& ex) {
1191 _comment += ex.DynamicType()->Name();
1192 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1194 _comment += ex.GetMessageString();
1198 if ( _comment.empty() )
1199 _comment = "Exception in blsurf_compute_mesh()";
1201 if ( status != STATUS_OK) {
1202 // There was an error while meshing
1203 blsurf_session_delete(bls);
1205 context_delete(ctx);
1207 return error(_comment);
1210 std::cout << std::endl;
1211 std::cout << "End of Surface Mesh generation" << std::endl;
1212 std::cout << std::endl;
1215 blsurf_data_get_mesh(bls, &msh);
1217 blsurf_session_delete(bls);
1219 context_delete(ctx);
1221 return error(_comment);
1225 /* retrieve mesh data (see distene/mesh.h) */
1226 integer nv, ne, nt, nq, vtx[4], tag;
1229 mesh_get_vertex_count(msh, &nv);
1230 mesh_get_edge_count(msh, &ne);
1231 mesh_get_triangle_count(msh, &nt);
1232 mesh_get_quadrangle_count(msh, &nq);
1235 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1236 SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1237 bool* tags = new bool[nv+1];
1239 /* enumerated vertices */
1240 for(int iv=1;iv<=nv;iv++) {
1241 mesh_get_vertex_coordinates(msh, iv, xyz);
1242 mesh_get_vertex_tag(msh, iv, &tag);
1243 // Issue 0020656. Use vertex coordinates
1245 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex(pmap(tag)));
1246 xyz[0] = p.X(); xyz[1] = p.Y(); xyz[2] = p.Z();
1248 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1249 // internal point are tagged to zero
1251 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1258 /* enumerate edges */
1259 for(int it=1;it<=ne;it++) {
1260 mesh_get_edge_vertices(msh, it, vtx);
1261 SMDS_MeshEdge* edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1262 mesh_get_edge_tag(msh, it, &tag);
1265 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1266 tags[vtx[0]] = false;
1269 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1270 tags[vtx[1]] = false;
1272 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1276 /* enumerate triangles */
1277 for(int it=1;it<=nt;it++) {
1278 mesh_get_triangle_vertices(msh, it, vtx);
1279 SMDS_MeshFace* tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1280 mesh_get_triangle_tag(msh, it, &tag);
1281 meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1283 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1284 tags[vtx[0]] = false;
1287 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1288 tags[vtx[1]] = false;
1291 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1292 tags[vtx[2]] = false;
1296 /* enumerate quadrangles */
1297 for(int it=1;it<=nq;it++) {
1298 mesh_get_quadrangle_vertices(msh, it, vtx);
1299 SMDS_MeshFace* quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1300 mesh_get_quadrangle_tag(msh, it, &tag);
1301 meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1303 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1304 tags[vtx[0]] = false;
1307 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1308 tags[vtx[1]] = false;
1311 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1312 tags[vtx[2]] = false;
1315 meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1316 tags[vtx[3]] = false;
1322 /* release the mesh object */
1323 blsurf_data_regain_mesh(bls, msh);
1325 /* clean up everything */
1326 blsurf_session_delete(bls);
1329 context_delete(ctx);
1331 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1333 if ( oldFEFlags > 0 )
1334 feenableexcept( oldFEFlags );
1335 feclearexcept( FE_ALL_EXCEPT );
1339 std::cout << "FacesWithSizeMap" << std::endl;
1340 FacesWithSizeMap.Statistics(std::cout);
1341 std::cout << "EdgesWithSizeMap" << std::endl;
1342 EdgesWithSizeMap.Statistics(std::cout);
1343 std::cout << "VerticesWithSizeMap" << std::endl;
1344 VerticesWithSizeMap.Statistics(std::cout);
1345 std::cout << "FacesWithEnforcedVertices" << std::endl;
1346 FacesWithEnforcedVertices.Statistics(std::cout);
1352 //=============================================================================
1356 //=============================================================================
1358 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1359 const TopoDS_Edge edge = TopoDS::Edge(ed);
1361 gp_Pnt pnt(node->X(), node->Y(), node->Z());
1363 Standard_Real p0 = 0.0;
1364 Standard_Real p1 = 1.0;
1365 TopLoc_Location loc;
1366 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1368 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1369 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1372 if ( proj.NbPoints() > 0 )
1374 pa = (double)proj.LowerDistanceParameter();
1375 // Issue 0020656. Move node if it is too far from edge
1376 gp_Pnt curve_pnt = curve->Value( pa );
1377 double dist2 = pnt.SquareDistance( curve_pnt );
1378 double tol = BRep_Tool::Tolerance( edge );
1379 if ( 1e-12 < dist2 && dist2 <= 2*tol*tol ) // large enough and within tolerance
1381 curve_pnt.Transform( loc );
1382 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1385 // GProp_GProps LProps;
1386 // BRepGProp::LinearProperties(ed, LProps);
1387 // double lg = (double)LProps.Mass();
1389 meshDS->SetNodeOnEdge(node, edge, pa);
1392 //=============================================================================
1396 //=============================================================================
1398 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
1403 //=============================================================================
1407 //=============================================================================
1409 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
1414 //=============================================================================
1418 //=============================================================================
1420 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
1422 return hyp.SaveTo( save );
1425 //=============================================================================
1429 //=============================================================================
1431 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
1433 return hyp.LoadFrom( load );
1436 /* Curve definition function See cad_curv_t in file distene/cad.h for
1438 * NOTE : if when your CAD systems evaluates second
1439 * order derivatives it also computes first order derivatives and
1440 * function evaluation, you can optimize this example by making only
1441 * one CAD call and filling the necessary uv, dt, dtt arrays.
1443 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
1445 /* t is given. It contains the t (time) 1D parametric coordintaes
1446 of the point PreCAD/BLSurf is querying on the curve */
1448 /* user_data identifies the edge PreCAD/BLSurf is querying
1449 * (see cad_edge_new later in this example) */
1450 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
1453 /* BLSurf is querying the function evaluation */
1456 uv[0]=P.X(); uv[1]=P.Y();
1460 /* query for the first order derivatives */
1463 dt[0]=V1.X(); dt[1]=V1.Y();
1467 /* query for the second order derivatives */
1470 dtt[0]=V2.X(); dtt[1]=V2.Y();
1476 /* Surface definition function.
1477 * See cad_surf_t in file distene/cad.h for more information.
1478 * NOTE : if when your CAD systems evaluates second order derivatives it also
1479 * computes first order derivatives and function evaluation, you can optimize
1480 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
1483 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1484 real *duu, real *duv, real *dvv, void *user_data)
1486 /* uv[2] is given. It contains the u,v coordinates of the point
1487 * PreCAD/BLSurf is querying on the surface */
1489 /* user_data identifies the face PreCAD/BLSurf is querying (see
1490 * cad_face_new later in this example)*/
1491 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
1495 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
1496 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
1503 geometry->D1(uv[0],uv[1],P,D1U,D1V);
1504 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
1505 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
1508 if(duu && duv && dvv){
1512 gp_Vec D2U,D2V,D2UV;
1514 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
1515 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
1516 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
1517 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
1524 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
1527 if (my_u_min > uv[0]) {
1530 if (my_v_min > uv[1]) {
1533 if (my_u_max < uv[0]) {
1536 if (my_v_max < uv[1]) {
1541 if (FaceId2PythonSmp.count(face_id) != 0){
1542 PyObject * pyresult = NULL;
1543 PyObject* new_stderr = NULL;
1544 assert(Py_IsInitialized());
1545 PyGILState_STATE gstate;
1546 gstate = PyGILState_Ensure();
1547 pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],"(f,f)",uv[0],uv[1]);
1549 if ( pyresult == NULL){
1551 string err_description="";
1552 new_stderr = newPyStdOut(err_description);
1553 PySys_SetObject("stderr", new_stderr);
1555 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1556 Py_DECREF(new_stderr);
1557 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
1558 result = *((double*)user_data);
1561 result = PyFloat_AsDouble(pyresult);
1562 Py_DECREF(pyresult);
1565 //MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1566 PyGILState_Release(gstate);
1569 *size = *((double*)user_data);
1574 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
1576 if (EdgeId2PythonSmp.count(edge_id) != 0){
1577 PyObject * pyresult = NULL;
1578 PyObject* new_stderr = NULL;
1579 assert(Py_IsInitialized());
1580 PyGILState_STATE gstate;
1581 gstate = PyGILState_Ensure();
1582 pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],"(f)",t);
1584 if ( pyresult == NULL){
1586 string err_description="";
1587 new_stderr = newPyStdOut(err_description);
1588 PySys_SetObject("stderr", new_stderr);
1590 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1591 Py_DECREF(new_stderr);
1592 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
1593 result = *((double*)user_data);
1596 result = PyFloat_AsDouble(pyresult);
1597 Py_DECREF(pyresult);
1600 PyGILState_Release(gstate);
1603 *size = *((double*)user_data);
1608 status_t size_on_vertex(integer point_id, real *size, void *user_data)
1610 if (VertexId2PythonSmp.count(point_id) != 0){
1611 PyObject * pyresult = NULL;
1612 PyObject* new_stderr = NULL;
1613 assert(Py_IsInitialized());
1614 PyGILState_STATE gstate;
1615 gstate = PyGILState_Ensure();
1616 pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],"");
1618 if ( pyresult == NULL){
1620 string err_description="";
1621 new_stderr = newPyStdOut(err_description);
1622 PySys_SetObject("stderr", new_stderr);
1624 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1625 Py_DECREF(new_stderr);
1626 MESSAGE("Can't evaluate f()" << " error is " << err_description);
1627 result = *((double*)user_data);
1630 result = PyFloat_AsDouble(pyresult);
1631 Py_DECREF(pyresult);
1634 PyGILState_Release(gstate);
1637 *size = *((double*)user_data);
1643 * The following function will be called for PreCAD/BLSurf message
1644 * printing. See context_set_message_callback (later in this
1645 * template) for how to set user_data.
1647 status_t message_cb(message_t *msg, void *user_data)
1649 integer errnumber = 0;
1651 message_get_number(msg, &errnumber);
1652 message_get_description(msg, &desc);
1653 if ( errnumber < 0 ) {
1654 string * error = (string*)user_data;
1655 // if ( !error->empty() )
1657 // remove ^A from the tail
1658 int len = strlen( desc );
1659 while (len > 0 && desc[len-1] != '\n')
1661 error->append( desc, len );
1664 std::cout << desc << std::endl;
1669 /* This is the interrupt callback. PreCAD/BLSurf will call this
1670 * function regularily. See the file distene/interrupt.h
1672 status_t interrupt_cb(integer *interrupt_status, void *user_data)
1674 integer you_want_to_continue = 1;
1676 if(you_want_to_continue)
1677 *interrupt_status = INTERRUPT_CONTINUE;
1678 else /* you want to stop BLSurf */
1679 *interrupt_status = INTERRUPT_STOP;
1684 //=============================================================================
1688 //=============================================================================
1689 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
1690 const TopoDS_Shape& aShape,
1691 MapShapeNbElems& aResMap)
1693 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
1694 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
1695 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
1696 //double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
1697 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
1698 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
1700 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
1701 _phySize = _hypothesis->GetPhySize();
1702 //_geometricMesh = (int) hyp->GetGeometricMesh();
1703 //_angleMeshS = hyp->GetAngleMeshS();
1704 _angleMeshC = _hypothesis->GetAngleMeshC();
1705 _quadAllowed = _hypothesis->GetQuadAllowed();
1708 bool IsQuadratic = false;
1713 TopTools_DataMapOfShapeInteger EdgesMap;
1714 double fullLen = 0.0;
1715 double fullNbSeg = 0;
1716 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1717 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
1718 if( EdgesMap.IsBound(E) )
1720 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
1721 double aLen = SMESH_Algo::EdgeLength(E);
1724 if(_physicalMesh==1) {
1725 nb1d = (int)( aLen/_phySize + 1 );
1730 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
1731 double fullAng = 0.0;
1732 double dp = (l-f)/200;
1737 for(int j=2; j<=200; j++) {
1740 fullAng += fabs(V1.Angle(V2));
1744 nb1d = (int)( fullAng/_angleMeshC + 1 );
1747 std::vector<int> aVec(SMDSEntity_Last);
1748 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1749 if( IsQuadratic > 0 ) {
1750 aVec[SMDSEntity_Node] = 2*nb1d - 1;
1751 aVec[SMDSEntity_Quad_Edge] = nb1d;
1754 aVec[SMDSEntity_Node] = nb1d - 1;
1755 aVec[SMDSEntity_Edge] = nb1d;
1757 aResMap.insert(std::make_pair(sm,aVec));
1758 EdgesMap.Bind(E,nb1d);
1760 double ELen = fullLen/fullNbSeg;
1764 // try to evaluate as in MEFISTO
1765 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1766 TopoDS_Face F = TopoDS::Face( exp.Current() );
1767 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1769 BRepGProp::SurfaceProperties(F,G);
1770 double anArea = G.Mass();
1772 std::vector<int> nb1dVec;
1773 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
1774 int nbSeg = EdgesMap.Find(exp1.Current());
1776 nb1dVec.push_back( nbSeg );
1779 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
1780 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
1783 if ( nb1dVec.size() == 4 ) // quadrangle geom face
1785 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
1787 nbNodes = (n1 + 1) * (n2 + 1);
1792 nbTria = nbQuad = nbTria / 3 + 1;
1795 std::vector<int> aVec(SMDSEntity_Last,0);
1797 int nb1d_in = (nbTria*3 - nb1d) / 2;
1798 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
1799 aVec[SMDSEntity_Quad_Triangle] = nbTria;
1800 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
1803 aVec[SMDSEntity_Node] = nbNodes;
1804 aVec[SMDSEntity_Triangle] = nbTria;
1805 aVec[SMDSEntity_Quadrangle] = nbQuad;
1807 aResMap.insert(std::make_pair(sm,aVec));
1814 BRepGProp::VolumeProperties(aShape,G);
1815 double aVolume = G.Mass();
1816 double tetrVol = 0.1179*ELen*ELen*ELen;
1817 int nbVols = int(aVolume/tetrVol);
1818 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
1819 std::vector<int> aVec(SMDSEntity_Last);
1820 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1822 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
1823 aVec[SMDSEntity_Quad_Tetra] = nbVols;
1826 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
1827 aVec[SMDSEntity_Tetra] = nbVols;
1829 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1830 aResMap.insert(std::make_pair(sm,aVec));
1835 //=============================================================================
1837 * Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
1838 * Following line was added:
1839 * myExtrem.Perform(P);
1841 //=============================================================================
1842 void BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
1843 const TopoDS_Face& face,
1845 const Standard_Real Tol)
1847 //-- Voir BRepExtrema_ExtPF.cxx
1848 BRepAdaptor_Surface Surf(face);
1849 Standard_Real U1, U2, V1, V2;
1850 BRepTools::UVBounds(face, U1, U2, V1, V2);
1851 Extrema_ExtPS myExtrem;
1852 myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
1853 myExtrem.Perform(P);
1854 //----------------------------------------------------------
1855 //-- On cherche le point le plus proche , PUIS
1856 //-- On le classifie.
1857 Standard_Integer nbv = 0; // xpu
1858 Standard_Real MaxDist = RealLast();
1859 Standard_Integer indice = 0;
1860 if(myExtrem.IsDone()) {
1861 nbv = myExtrem.NbExt();
1862 for (Standard_Integer i = 1; i <= nbv; i++) {
1863 Standard_Real d = myExtrem.Value(i);
1873 Standard_Real U1,U2;
1874 myExtrem.Point(indice).Parameter(U1, U2);
1875 Puv.SetCoord(U1, U2);
1876 fc->Perform(face, Puv, Tol);
1879 fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
1880 //-- le resultat est TopAbs_OUT;