1 // Copyright (C) 2007-2008 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
20 // File : BLSURFPlugin_BLSURF.cxx
21 // Authors : Francis KLOSS (OCC) & Patrick LAUG (INRIA) & Lioka RAZAFINDRAZAKA (CEA)
22 // & Aurelien ALLEAUME (DISTENE)
23 // Size maps developement: Nicolas GEIMER (OCC) & Gilles DAVID (EURIWARE)
26 #include "BLSURFPlugin_BLSURF.hxx"
27 #include "BLSURFPlugin_Hypothesis.hxx"
30 #include <distene/api.h>
31 #include <distene/blsurf.h>
34 #include <structmember.h>
37 #include <Basics_Utils.hxx>
39 #include <SMESH_Gen.hxx>
40 #include <SMESH_Mesh.hxx>
41 #include <SMESH_ControlsDef.hxx>
43 #include <utilities.h>
51 // OPENCASCADE includes
52 #include <BRep_Tool.hxx>
54 #include <TopExp_Explorer.hxx>
56 #include <NCollection_Map.hxx>
58 #include <Geom_Surface.hxx>
59 #include <Handle_Geom_Surface.hxx>
60 #include <Geom2d_Curve.hxx>
61 #include <Handle_Geom2d_Curve.hxx>
62 #include <Geom_Curve.hxx>
63 #include <Handle_Geom_Curve.hxx>
64 #include <Handle_AIS_InteractiveObject.hxx>
65 #include <TopoDS_Vertex.hxx>
66 #include <TopoDS_Edge.hxx>
67 #include <TopoDS_Wire.hxx>
68 #include <TopoDS_Face.hxx>
70 #include <gp_Pnt2d.hxx>
71 #include <TopTools_IndexedMapOfShape.hxx>
72 #include <TopoDS_Shape.hxx>
73 #include <BRep_Builder.hxx>
74 #include <BRepTools.hxx>
76 #include <TopTools_DataMapOfShapeInteger.hxx>
77 #include <GProp_GProps.hxx>
78 #include <BRepGProp.hxx>
84 #include <Standard_ErrorHandler.hxx>
85 #include <GeomAPI_ProjectPointOnCurve.hxx>
86 #include <GeomAPI_ProjectPointOnSurf.hxx>
89 // #include <BRepClass_FaceClassifier.hxx>
90 #include <TopTools_MapOfShape.hxx>
92 /* ==================================
93 * =========== PYTHON ==============
94 * ==================================*/
103 PyStdOut_dealloc(PyStdOut *self)
109 PyStdOut_write(PyStdOut *self, PyObject *args)
113 if (!PyArg_ParseTuple(args, "t#:write",&c, &l))
117 *(self->out)=*(self->out)+c;
123 static PyMethodDef PyStdOut_methods[] = {
124 {"write", (PyCFunction)PyStdOut_write, METH_VARARGS,
125 PyDoc_STR("write(string) -> None")},
126 {NULL, NULL} /* sentinel */
129 static PyMemberDef PyStdOut_memberlist[] = {
130 {"softspace", T_INT, offsetof(PyStdOut, softspace), 0,
131 "flag indicating that a space needs to be printed; used by print"},
132 {NULL} /* Sentinel */
135 static PyTypeObject PyStdOut_Type = {
136 /* The ob_type field must be initialized in the module init function
137 * to be portable to Windows without using C++. */
138 PyObject_HEAD_INIT(NULL)
141 sizeof(PyStdOut), /*tp_basicsize*/
144 (destructor)PyStdOut_dealloc, /*tp_dealloc*/
151 0, /*tp_as_sequence*/
156 PyObject_GenericGetAttr, /*tp_getattro*/
157 /* softspace is writable: we must supply tp_setattro */
158 PyObject_GenericSetAttr, /* tp_setattro */
160 Py_TPFLAGS_DEFAULT, /*tp_flags*/
164 0, /*tp_richcompare*/
165 0, /*tp_weaklistoffset*/
168 PyStdOut_methods, /*tp_methods*/
169 PyStdOut_memberlist, /*tp_members*/
183 PyObject * newPyStdOut( std::string& out )
186 self = PyObject_New(PyStdOut, &PyStdOut_Type);
191 return (PyObject*)self;
195 ////////////////////////END PYTHON///////////////////////////
197 //////////////////MY MAPS////////////////////////////////////////
198 TopTools_IndexedMapOfShape FacesWithSizeMap;
199 std::map<int,string> FaceId2SizeMap;
200 TopTools_IndexedMapOfShape EdgesWithSizeMap;
201 std::map<int,string> EdgeId2SizeMap;
202 TopTools_IndexedMapOfShape VerticesWithSizeMap;
203 std::map<int,string> VertexId2SizeMap;
205 std::map<int,PyObject*> FaceId2PythonSmp;
206 std::map<int,PyObject*> EdgeId2PythonSmp;
207 std::map<int,PyObject*> VertexId2PythonSmp;
209 std::map<int,std::vector<double> > FaceId2AttractorCoords;
211 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
212 std::map< int, std::set< std::vector<double> > > FaceId2EnforcedVertexCoords;
214 bool HasSizeMapOnFace=false;
215 bool HasSizeMapOnEdge=false;
216 bool HasSizeMapOnVertex=false;
218 //=============================================================================
222 //=============================================================================
224 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId, int studyId,
226 : SMESH_2D_Algo(hypId, studyId, gen)
228 MESSAGE("BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF");
231 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
232 _compatibleHypothesis.push_back("BLSURF_Parameters");
233 _requireDescretBoundary = false;
234 _onlyUnaryInput = false;
237 smeshGen_i = SMESH_Gen_i::GetSMESHGen();
238 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
239 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
241 MESSAGE("studyid = " << _studyId);
244 myStudy = aStudyMgr->GetStudyByID(_studyId);
245 MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
247 /* Initialize the Python interpreter */
248 assert(Py_IsInitialized());
249 PyGILState_STATE gstate;
250 gstate = PyGILState_Ensure();
253 main_mod = PyImport_AddModule("__main__");
256 main_dict = PyModule_GetDict(main_mod);
258 PyRun_SimpleString("from math import *");
259 PyGILState_Release(gstate);
261 FacesWithSizeMap.Clear();
262 FaceId2SizeMap.clear();
263 EdgesWithSizeMap.Clear();
264 EdgeId2SizeMap.clear();
265 VerticesWithSizeMap.Clear();
266 VertexId2SizeMap.clear();
267 FaceId2PythonSmp.clear();
268 EdgeId2PythonSmp.clear();
269 VertexId2PythonSmp.clear();
270 FaceId2AttractorCoords.clear();
271 FacesWithEnforcedVertices.Clear();
272 FaceId2EnforcedVertexCoords.clear();
276 //=============================================================================
280 //=============================================================================
282 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
284 MESSAGE("BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF");
288 //=============================================================================
292 //=============================================================================
294 bool BLSURFPlugin_BLSURF::CheckHypothesis
296 const TopoDS_Shape& aShape,
297 SMESH_Hypothesis::Hypothesis_Status& aStatus)
301 list<const SMESHDS_Hypothesis*>::const_iterator itl;
302 const SMESHDS_Hypothesis* theHyp;
304 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
305 int nbHyp = hyps.size();
308 aStatus = SMESH_Hypothesis::HYP_OK;
309 return true; // can work with no hypothesis
313 theHyp = (*itl); // use only the first hypothesis
315 string hypName = theHyp->GetName();
317 if (hypName == "BLSURF_Parameters")
319 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
321 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
322 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
323 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
324 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
326 aStatus = SMESH_Hypothesis::HYP_OK;
329 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
331 return aStatus == SMESH_Hypothesis::HYP_OK;
334 //=============================================================================
336 * Pass parameters to BLSURF
338 //=============================================================================
340 inline std::string to_string(double d)
342 std::ostringstream o;
347 inline std::string to_string(int i)
349 std::ostringstream o;
354 double _smp_phy_size;
355 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
356 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
357 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
359 double my_u_min=1e6,my_v_min=1e6,my_u_max=-1e6,my_v_max=-1e6;
365 /////////////////////////////////////////////////////////
366 projectionPoint getProjectionPoint(const TopoDS_Face& face, const gp_XYZ& point)
368 projectionPoint myPoint;
369 Handle(Geom_Surface) surface = BRep_Tool::Surface(face);
370 GeomAPI_ProjectPointOnSurf projector( point, surface );
371 if ( !projector.IsDone() || projector.NbPoints()==0 )
372 throw "getProjectionPoint: Can't project";
374 Quantity_Parameter u,v;
375 projector.LowerDistanceParameters(u,v);
376 myPoint.uv = gp_XY(u,v);
377 gp_Pnt aPnt = projector.NearestPoint();
378 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
382 /////////////////////////////////////////////////////////
384 /////////////////////////////////////////////////////////
385 double getT(const TopoDS_Edge& edge, const gp_XYZ& point)
388 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, f,l);
389 GeomAPI_ProjectPointOnCurve projector( point, curve);
390 if ( projector.NbPoints() == 0 )
392 return projector.LowerDistanceParameter();
395 /////////////////////////////////////////////////////////
396 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
398 MESSAGE("BLSURFPlugin_BLSURF::entryToShape "<<entry );
399 GEOM::GEOM_Object_var aGeomObj;
400 TopoDS_Shape S = TopoDS_Shape();
401 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
402 SALOMEDS::GenericAttribute_var anAttr;
404 if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR")) {
405 SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
406 CORBA::String_var aVal = anIOR->Value();
407 CORBA::Object_var obj = myStudy->ConvertIORToObject(aVal);
408 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
410 if ( !aGeomObj->_is_nil() )
411 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
415 /////////////////////////////////////////////////////////
416 void createEnforcedVertexOnFace(TopoDS_Shape GeomShape, BLSURFPlugin_Hypothesis::TEnforcedVertexList enforcedVertexList)
419 std::vector<double> coords;
420 BLSURFPlugin_Hypothesis::TEnforcedVertex enforcedVertex;
421 BLSURFPlugin_Hypothesis::TEnforcedVertexList::const_iterator evlIt = enforcedVertexList.begin();
423 for( ; evlIt != enforcedVertexList.end() ; ++evlIt ) {
425 enforcedVertex = *evlIt;
426 xe = enforcedVertex[0];
427 ye = enforcedVertex[1];
428 ze = enforcedVertex[2];
429 MESSAGE("Enforced Vertex: " << xe << ", " << ye << ", " << ze);
430 // Get the (u,v) values of the enforced vertex on the face
431 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_XYZ(xe,ye,ze));
432 gp_XY uvPoint = myPoint.uv;
433 gp_XYZ xyzPoint = myPoint.xyz;
434 Standard_Real u0 = uvPoint.X();
435 Standard_Real v0 = uvPoint.Y();
436 Standard_Real x0 = xyzPoint.X();
437 Standard_Real y0 = xyzPoint.Y();
438 Standard_Real z0 = xyzPoint.Z();
439 MESSAGE("Projected Vertex: " << x0 << ", " << y0 << ", " << z0);
440 coords.push_back(u0);
441 coords.push_back(v0);
442 coords.push_back(x0);
443 coords.push_back(y0);
444 coords.push_back(z0);
447 if (! FacesWithEnforcedVertices.Contains(TopoDS::Face(GeomShape))) {
448 key = FacesWithEnforcedVertices.Add(TopoDS::Face(GeomShape));
451 key = FacesWithEnforcedVertices.FindIndex(TopoDS::Face(GeomShape));
454 // If a node is already created by an attractor, do not create enforced vertex
455 int attractorKey = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
456 bool sameAttractor = false;
457 if (attractorKey >= 0)
458 if (FaceId2AttractorCoords.count(attractorKey) > 0)
459 if (FaceId2AttractorCoords[attractorKey] == coords)
460 sameAttractor = true;
462 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
463 MESSAGE("Map of enf. vertex has key " << key)
464 MESSAGE("Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
466 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
468 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
469 MESSAGE("New Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
472 MESSAGE("Map of enf. vertex has not key " << key << ": creating it")
473 if (! sameAttractor) {
474 std::set< std::vector<double> > ens;
476 FaceId2EnforcedVertexCoords[key] = ens;
479 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
484 /////////////////////////////////////////////////////////
485 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction)
487 MESSAGE("Attractor function: "<< AttractorFunction);
488 double xa, ya, za; // Coordinates of attractor point
489 double a, b; // Attractor parameter
490 bool createNode=false; // To create a node on attractor projection
492 const char *sep = ";";
493 // atIt->second has the following pattern:
494 // ATTRACTOR(xa;ya;za;a;b)
496 // xa;ya;za : coordinates of attractor
497 // a : desired size on attractor
498 // b : distance of influence of attractor
500 // We search the parameters in the string
502 pos1 = AttractorFunction.find(sep);
503 if (pos1!=string::npos)
504 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
506 pos2 = AttractorFunction.find(sep, pos1+1);
507 if (pos2!=string::npos) {
508 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
512 pos2 = AttractorFunction.find(sep, pos1+1);
513 if (pos2!=string::npos) {
514 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
518 pos2 = AttractorFunction.find(sep, pos1+1);
519 if (pos2!=string::npos) {
520 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
524 pos2 = AttractorFunction.find(sep, pos1+1);
525 if (pos2!=string::npos) {
526 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
530 pos2 = AttractorFunction.find(")");
531 if (pos2!=string::npos) {
532 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
533 MESSAGE("createNode: " << createNodeStr);
534 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
537 // Get the (u,v) values of the attractor on the face
538 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_XYZ(xa,ya,za));
539 gp_XY uvPoint = myPoint.uv;
540 gp_XYZ xyzPoint = myPoint.xyz;
541 Standard_Real u0 = uvPoint.X();
542 Standard_Real v0 = uvPoint.Y();
543 Standard_Real x0 = xyzPoint.X();
544 Standard_Real y0 = xyzPoint.Y();
545 Standard_Real z0 = xyzPoint.Z();
546 std::vector<double> coords;
547 coords.push_back(u0);
548 coords.push_back(v0);
549 coords.push_back(x0);
550 coords.push_back(y0);
551 coords.push_back(z0);
552 // We construct the python function
553 ostringstream attractorFunctionStream;
554 attractorFunctionStream << "def f(u,v): return ";
555 attractorFunctionStream << _smp_phy_size << "-(" << _smp_phy_size <<"-" << a << ")";
556 attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
558 MESSAGE("Python function for attractor:" << std::endl << attractorFunctionStream.str());
561 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
562 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
565 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
567 FaceId2SizeMap[key] =attractorFunctionStream.str();
569 MESSAGE("Creating node on ("<<x0<<","<<y0<<","<<z0<<")");
570 FaceId2AttractorCoords[key] = coords;
574 /////////////////////////////////////////////////////////
576 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp, blsurf_session_t *bls)
578 int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
579 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
580 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
581 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
582 double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
583 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
584 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
585 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
586 bool _decimesh = BLSURFPlugin_Hypothesis::GetDefaultDecimesh();
587 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
590 MESSAGE("BLSURFPlugin_BLSURF::SetParameters");
591 _topology = (int) hyp->GetTopology();
592 _physicalMesh = (int) hyp->GetPhysicalMesh();
593 _phySize = hyp->GetPhySize();
594 _geometricMesh = (int) hyp->GetGeometricMesh();
595 _angleMeshS = hyp->GetAngleMeshS();
596 _angleMeshC = hyp->GetAngleMeshC();
597 _gradation = hyp->GetGradation();
598 _quadAllowed = hyp->GetQuadAllowed();
599 _decimesh = hyp->GetDecimesh();
600 _verb = hyp->GetVerbosity();
602 if ( hyp->GetPhyMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
603 blsurf_set_param(bls, "hphymin", to_string(hyp->GetPhyMin()).c_str());
604 if ( hyp->GetPhyMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
605 blsurf_set_param(bls, "hphymax", to_string(hyp->GetPhyMax()).c_str());
606 if ( hyp->GetGeoMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
607 blsurf_set_param(bls, "hgeomin", to_string(hyp->GetGeoMin()).c_str());
608 if ( hyp->GetGeoMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
609 blsurf_set_param(bls, "hgeomax", to_string(hyp->GetGeoMax()).c_str());
611 const BLSURFPlugin_Hypothesis::TOptionValues & opts = hyp->GetOptionValues();
612 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
613 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt )
614 if ( !opIt->second.empty() ) {
615 MESSAGE("blsurf_set_param(): " << opIt->first << " = " << opIt->second);
616 blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
620 MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
622 _smp_phy_size = _phySize;
623 blsurf_set_param(bls, "topo_points", _topology > 0 ? "1" : "0");
624 blsurf_set_param(bls, "topo_curves", _topology > 0 ? "1" : "0");
625 blsurf_set_param(bls, "topo_project", _topology > 0 ? "1" : "0");
626 blsurf_set_param(bls, "clean_boundary", _topology > 1 ? "1" : "0");
627 blsurf_set_param(bls, "close_boundary", _topology > 1 ? "1" : "0");
628 blsurf_set_param(bls, "hphy_flag", to_string(_physicalMesh).c_str());
629 // blsurf_set_param(bls, "hphy_flag", "2");
630 if ((to_string(_physicalMesh))=="2"){
631 TopoDS_Shape GeomShape;
632 TopAbs_ShapeEnum GeomType;
634 // Standard Size Maps
636 MESSAGE("Setting a Size Map");
637 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
638 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
639 for ( ; smIt != sizeMaps.end(); ++smIt ) {
640 if ( !smIt->second.empty() ) {
641 MESSAGE("blsurf_set_sizeMap(): " << smIt->first << " = " << smIt->second);
642 GeomShape = entryToShape(smIt->first);
643 GeomType = GeomShape.ShapeType();
644 MESSAGE("Geomtype is " << GeomType);
647 if (GeomType == TopAbs_COMPOUND){
648 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
650 if (it.Value().ShapeType() == TopAbs_FACE){
651 HasSizeMapOnFace = true;
652 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
653 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
656 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
657 // MESSAGE("Face with key " << key << " already in map");
659 FaceId2SizeMap[key] = smIt->second;
662 if (it.Value().ShapeType() == TopAbs_EDGE){
663 HasSizeMapOnEdge = true;
664 HasSizeMapOnFace = true;
665 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
666 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
669 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
670 // MESSAGE("Edge with key " << key << " already in map");
672 EdgeId2SizeMap[key] = smIt->second;
675 if (it.Value().ShapeType() == TopAbs_VERTEX){
676 HasSizeMapOnVertex = true;
677 HasSizeMapOnEdge = true;
678 HasSizeMapOnFace = true;
679 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
680 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
683 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
684 MESSAGE("Group of vertices with key " << key << " already in map");
686 MESSAGE("Group of vertices with key " << key << " has a size map: " << smIt->second);
687 VertexId2SizeMap[key] = smIt->second;
692 if (GeomType == TopAbs_FACE){
693 HasSizeMapOnFace = true;
694 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
695 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
698 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
699 // MESSAGE("Face with key " << key << " already in map");
701 FaceId2SizeMap[key] = smIt->second;
704 if (GeomType == TopAbs_EDGE){
705 HasSizeMapOnEdge = true;
706 HasSizeMapOnFace = true;
707 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
708 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
711 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
712 // MESSAGE("Edge with key " << key << " already in map");
714 EdgeId2SizeMap[key] = smIt->second;
717 if (GeomType == TopAbs_VERTEX){
718 HasSizeMapOnVertex = true;
719 HasSizeMapOnEdge = true;
720 HasSizeMapOnFace = true;
721 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
722 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
725 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
726 MESSAGE("Vertex with key " << key << " already in map");
728 MESSAGE("Vertex with key " << key << " has a size map: " << smIt->second);
729 VertexId2SizeMap[key] = smIt->second;
737 MESSAGE("Setting Attractors");
738 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
739 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
740 for ( ; atIt != attractors.end(); ++atIt ) {
741 if ( !atIt->second.empty() ) {
742 MESSAGE("blsurf_set_attractor(): " << atIt->first << " = " << atIt->second);
743 GeomShape = entryToShape(atIt->first);
744 GeomType = GeomShape.ShapeType();
746 if (GeomType == TopAbs_COMPOUND){
747 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
748 if (it.Value().ShapeType() == TopAbs_FACE){
749 HasSizeMapOnFace = true;
750 createAttractorOnFace(it.Value(), atIt->second);
755 if (GeomType == TopAbs_FACE){
756 HasSizeMapOnFace = true;
757 createAttractorOnFace(GeomShape, atIt->second);
760 if (GeomType == TopAbs_EDGE){
761 HasSizeMapOnEdge = true;
762 HasSizeMapOnFace = true;
763 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
765 if (GeomType == TopAbs_VERTEX){
766 HasSizeMapOnVertex = true;
767 HasSizeMapOnEdge = true;
768 HasSizeMapOnFace = true;
769 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
779 MESSAGE("Setting Enforced Vertices");
780 const BLSURFPlugin_Hypothesis::TEnforcedVertexMap enforcedVertexMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVertices(hyp);
781 BLSURFPlugin_Hypothesis::TEnforcedVertexMap::const_iterator enfIt = enforcedVertexMap.begin();
782 for ( ; enfIt != enforcedVertexMap.end(); ++enfIt ) {
783 if ( !enfIt->second.empty() ) {
784 GeomShape = entryToShape(enfIt->first);
785 GeomType = GeomShape.ShapeType();
787 if (GeomType == TopAbs_COMPOUND){
788 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
789 if (it.Value().ShapeType() == TopAbs_FACE){
790 HasSizeMapOnFace = true;
791 createEnforcedVertexOnFace(it.Value(), enfIt->second);
796 if (GeomType == TopAbs_FACE){
797 HasSizeMapOnFace = true;
798 createEnforcedVertexOnFace(GeomShape, enfIt->second);
803 // if (HasSizeMapOnFace){
804 // In all size map cases (hphy_flag = 2), at least map on face must be defined
805 MESSAGE("Setting Size Map on FACES ");
806 blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
809 if (HasSizeMapOnEdge){
810 MESSAGE("Setting Size Map on EDGES ");
811 blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
813 if (HasSizeMapOnVertex){
814 MESSAGE("Setting Size Map on VERTICES ");
815 blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
818 blsurf_set_param(bls, "hphydef", to_string(_phySize).c_str());
819 blsurf_set_param(bls, "hgeo_flag", to_string(_geometricMesh).c_str());
820 blsurf_set_param(bls, "relax_size", _decimesh ? "0": to_string(_geometricMesh).c_str());
821 blsurf_set_param(bls, "angle_meshs", to_string(_angleMeshS).c_str());
822 blsurf_set_param(bls, "angle_meshc", to_string(_angleMeshC).c_str());
823 blsurf_set_param(bls, "gradation", to_string(_gradation).c_str());
824 blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
825 blsurf_set_param(bls, "element", _quadAllowed ? "q1.0" : "p1");
826 blsurf_set_param(bls, "verb", to_string(_verb).c_str());
829 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
830 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
831 real *duu, real *duv, real *dvv, void *user_data);
832 status_t message_cb(message_t *msg, void *user_data);
833 status_t interrupt_cb(integer *interrupt_status, void *user_data);
835 //=============================================================================
839 //=============================================================================
841 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
843 MESSAGE("BLSURFPlugin_BLSURF::Compute");
845 // if (aShape.ShapeType() == TopAbs_COMPOUND) {
846 // MESSAGE(" the shape is a COMPOUND");
849 // MESSAGE(" the shape is UNKNOWN");
852 // Fix problem with locales
853 Kernel_Utils::Localizer loc;
855 /* create a distene context (generic object) */
856 status_t status = STATUS_ERROR;
858 context_t *ctx = context_new();
860 /* Set the message callback in the working context */
861 context_set_message_callback(ctx, message_cb, &_comment);
862 context_set_interrupt_callback(ctx, interrupt_cb, NULL);
864 /* create the CAD object we will work on. It is associated to the context ctx. */
865 cad_t *c = cad_new(ctx);
867 blsurf_session_t *bls = blsurf_session_new(ctx);
869 FacesWithSizeMap.Clear();
870 FaceId2SizeMap.clear();
871 EdgesWithSizeMap.Clear();
872 EdgeId2SizeMap.clear();
873 VerticesWithSizeMap.Clear();
874 VertexId2SizeMap.clear();
876 MESSAGE("BEGIN SetParameters");
877 SetParameters(_hypothesis, bls);
878 MESSAGE("END SetParameters");
880 /* Now fill the CAD object with data from your CAD
881 * environement. This is the most complex part of a successfull
885 // needed to prevent the opencascade memory managmement from freeing things
886 vector<Handle(Geom2d_Curve)> curves;
887 vector<Handle(Geom_Surface)> surfaces;
892 TopTools_IndexedMapOfShape fmap;
893 TopTools_IndexedMapOfShape emap;
894 TopTools_IndexedMapOfShape pmap;
897 FaceId2PythonSmp.clear();
899 EdgeId2PythonSmp.clear();
901 VertexId2PythonSmp.clear();
903 assert(Py_IsInitialized());
904 PyGILState_STATE gstate;
905 gstate = PyGILState_Ensure();
907 string theSizeMapStr;
909 /****************************************************************************************
911 *****************************************************************************************/
913 string bad_end = "return";
915 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
916 TopoDS_Face f=TopoDS::Face(face_iter.Current());
918 if (fmap.FindIndex(f) > 0)
923 surfaces.push_back(BRep_Tool::Surface(f));
925 /* create an object representing the face for blsurf */
926 /* where face_id is an integer identifying the face.
927 * surf_function is the function that defines the surface
928 * (For this face, it will be called by blsurf with your_face_object_ptr
931 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
933 /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
934 cad_face_set_tag(fce, iface);
936 /* Set face orientation (optional if you want a well oriented output mesh)*/
937 if(f.Orientation() != TopAbs_FORWARD){
938 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
940 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
943 if (HasSizeMapOnFace){
944 std::cout << "A size map is defined on a face" << std::endl;
946 faceKey = FacesWithSizeMap.FindIndex(f);
948 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
949 theSizeMapStr = FaceId2SizeMap[faceKey];
950 // check if function ends with "return"
951 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
953 // Expr To Python function, verification is performed at validation in GUI
954 PyObject * obj = NULL;
955 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
957 PyObject * func = NULL;
958 func = PyObject_GetAttrString(main_mod, "f");
959 FaceId2PythonSmp[iface]=func;
960 FaceId2SizeMap.erase(faceKey);
963 // Specific size map = Attractor
964 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
966 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
967 if (attractor_iter->first == faceKey) {
968 MESSAGE("Face indice: " << iface);
969 MESSAGE("Adding attractor");
971 double xyzCoords[3] = {attractor_iter->second[2],
972 attractor_iter->second[3],
973 attractor_iter->second[4]};
975 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
976 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
977 BRepClass_FaceClassifier scl(f,P,1e-7);
978 // scl.Perform() is bugged. The function was rewritten
980 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
981 TopAbs_State result = scl.State();
982 MESSAGE("Position of point on face: "<<result);
983 if ( result == TopAbs_OUT )
984 MESSAGE("Point is out of face: node is not created");
985 if ( result == TopAbs_UNKNOWN )
986 MESSAGE("Point position on face is unknown: node is not created");
987 if ( result == TopAbs_ON )
988 MESSAGE("Point is on border of face: node is not created");
989 if ( result == TopAbs_IN )
991 // Point is inside face and not on border
992 MESSAGE("Point is in face: node is created");
993 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
995 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << iatt);
996 cad_point_t* point_p = cad_point_new(fce, iatt, uvCoords);
997 cad_point_set_tag(point_p, iatt);
999 FaceId2AttractorCoords.erase(faceKey);
1003 // Enforced Vertices
1004 faceKey = FacesWithEnforcedVertices.FindIndex(f);
1005 std::map<int,std::set<std::vector<double> > >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1006 if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1007 std::cout << "Some enforced vertices are defined" << std::endl;
1009 std::set<std::vector<double> > evl;
1010 // std::vector<double> ev;
1011 MESSAGE("Face indice: " << iface);
1012 MESSAGE("Adding enforced vertices");
1013 evl = evmIt->second;
1014 MESSAGE("Number of vertices to add: "<< evl.size())
1015 std::set< std::vector<double> >::const_iterator evlIt = evl.begin();
1016 for (; evlIt != evl.end(); ++evlIt) {
1018 // for (int i=0; i<evl.size() ; i++) {
1021 // double xyzCoords[3] = {ev[2], ev[3], ev[4]};
1022 double xyzCoords[3] = {evlIt->at(0), evlIt->at(3), evlIt->at(4)};
1023 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1024 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1025 BRepClass_FaceClassifier scl(f,P,1e-7);
1026 // scl.Perform() is bugged. The function was rewritten
1028 BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1029 TopAbs_State result = scl.State();
1030 MESSAGE("Position of point on face: "<<result);
1031 if ( result == TopAbs_OUT )
1032 MESSAGE("Point is out of face: node is not created");
1033 if ( result == TopAbs_UNKNOWN )
1034 MESSAGE("Point position on face is unknown: node is not created");
1035 if ( result == TopAbs_ON )
1036 MESSAGE("Point is on border of face: node is not created");
1037 if ( result == TopAbs_IN )
1039 // Point is inside face and not on border
1040 MESSAGE("Point is in face: node is created");
1041 // double uvCoords[2] = {ev[0],ev[1]};
1042 double uvCoords[2] = {evlIt->at(0),evlIt->at(1)};
1044 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1045 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1046 cad_point_set_tag(point_p, ienf);
1049 FaceId2EnforcedVertexCoords.erase(faceKey);
1052 std::cout << "No enforced vertex defined" << std::endl;
1056 /****************************************************************************************
1058 now create the edges associated to this face
1059 *****************************************************************************************/
1061 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1062 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1063 int ic = emap.FindIndex(e);
1068 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1070 if (HasSizeMapOnEdge){
1071 edgeKey = EdgesWithSizeMap.FindIndex(e);
1072 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1073 MESSAGE("Sizemap defined on edge with index " << edgeKey);
1074 theSizeMapStr = EdgeId2SizeMap[edgeKey];
1075 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1077 // Expr To Python function, verification is performed at validation in GUI
1078 PyObject * obj = NULL;
1079 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1081 PyObject * func = NULL;
1082 func = PyObject_GetAttrString(main_mod, "f");
1083 EdgeId2PythonSmp[ic]=func;
1084 EdgeId2SizeMap.erase(edgeKey);
1088 /* attach the edge to the current blsurf face */
1089 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1091 /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1092 cad_edge_set_tag(edg, ic);
1094 /* by default, an edge does not necessalry appear in the resulting mesh,
1095 unless the following property is set :
1097 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1099 /* by default an edge is a boundary edge */
1100 if (e.Orientation() == TopAbs_INTERNAL)
1101 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1105 gp_Pnt2d e0 = curves.back()->Value(tmin);
1106 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1107 Standard_Real d1=0,d2=0;
1110 /****************************************************************************************
1112 *****************************************************************************************/
1114 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1115 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1119 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1122 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1124 *ip = pmap.FindIndex(v);
1128 //vertexKey = VerticesWithSizeMap.FindIndex(v);
1129 if (HasSizeMapOnVertex){
1130 vertexKey = VerticesWithSizeMap.FindIndex(v);
1131 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1132 theSizeMapStr = VertexId2SizeMap[vertexKey];
1133 //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1134 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1136 // Expr To Python function, verification is performed at validation in GUI
1137 PyObject * obj = NULL;
1138 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1140 PyObject * func = NULL;
1141 func = PyObject_GetAttrString(main_mod, "f");
1142 VertexId2PythonSmp[*ip]=func;
1143 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
1148 // should not happen
1149 MESSAGE("An edge does not have 2 extremities.");
1152 // This defines the curves extremity connectivity
1153 cad_edge_set_extremities(edg, ip1, ip2);
1154 /* set the tag (color) to the same value as the extremity id : */
1155 cad_edge_set_extremities_tag(edg, ip1, ip2);
1158 cad_edge_set_extremities(edg, ip2, ip1);
1159 cad_edge_set_extremities_tag(edg, ip2, ip1);
1166 PyGILState_Release(gstate);
1168 blsurf_data_set_cad(bls, c);
1170 std::cout << std::endl;
1171 std::cout << "Beginning of Surface Mesh generation" << std::endl;
1172 std::cout << std::endl;
1174 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1176 feclearexcept( FE_ALL_EXCEPT );
1177 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1183 status = blsurf_compute_mesh(bls);
1186 catch ( std::exception& exc ) {
1187 _comment += exc.what();
1189 catch (Standard_Failure& ex) {
1190 _comment += ex.DynamicType()->Name();
1191 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1193 _comment += ex.GetMessageString();
1197 if ( _comment.empty() )
1198 _comment = "Exception in blsurf_compute_mesh()";
1200 if ( status != STATUS_OK) {
1201 // There was an error while meshing
1202 blsurf_session_delete(bls);
1204 context_delete(ctx);
1206 return error(_comment);
1209 std::cout << std::endl;
1210 std::cout << "End of Surface Mesh generation" << std::endl;
1211 std::cout << std::endl;
1214 blsurf_data_get_mesh(bls, &msh);
1216 blsurf_session_delete(bls);
1218 context_delete(ctx);
1220 return error(_comment);
1224 /* retrieve mesh data (see distene/mesh.h) */
1225 integer nv, ne, nt, nq, vtx[4], tag;
1228 mesh_get_vertex_count(msh, &nv);
1229 mesh_get_edge_count(msh, &ne);
1230 mesh_get_triangle_count(msh, &nt);
1231 mesh_get_quadrangle_count(msh, &nq);
1234 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1235 SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1236 bool* tags = new bool[nv+1];
1238 /* enumerated vertices */
1239 for(int iv=1;iv<=nv;iv++) {
1240 mesh_get_vertex_coordinates(msh, iv, xyz);
1241 mesh_get_vertex_tag(msh, iv, &tag);
1242 // Issue 0020656. Use vertex coordinates
1244 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex(pmap(tag)));
1245 xyz[0] = p.X(); xyz[1] = p.Y(); xyz[2] = p.Z();
1247 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1248 // internal point are tagged to zero
1250 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1257 /* enumerate edges */
1258 for(int it=1;it<=ne;it++) {
1259 mesh_get_edge_vertices(msh, it, vtx);
1260 SMDS_MeshEdge* edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1261 mesh_get_edge_tag(msh, it, &tag);
1264 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1265 tags[vtx[0]] = false;
1268 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1269 tags[vtx[1]] = false;
1271 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1275 /* enumerate triangles */
1276 for(int it=1;it<=nt;it++) {
1277 mesh_get_triangle_vertices(msh, it, vtx);
1278 SMDS_MeshFace* tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1279 mesh_get_triangle_tag(msh, it, &tag);
1280 meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1282 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1283 tags[vtx[0]] = false;
1286 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1287 tags[vtx[1]] = false;
1290 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1291 tags[vtx[2]] = false;
1295 /* enumerate quadrangles */
1296 for(int it=1;it<=nq;it++) {
1297 mesh_get_quadrangle_vertices(msh, it, vtx);
1298 SMDS_MeshFace* quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1299 mesh_get_quadrangle_tag(msh, it, &tag);
1300 meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1302 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1303 tags[vtx[0]] = false;
1306 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1307 tags[vtx[1]] = false;
1310 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1311 tags[vtx[2]] = false;
1314 meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1315 tags[vtx[3]] = false;
1321 /* release the mesh object */
1322 blsurf_data_regain_mesh(bls, msh);
1324 /* clean up everything */
1325 blsurf_session_delete(bls);
1328 context_delete(ctx);
1330 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1332 if ( oldFEFlags > 0 )
1333 feenableexcept( oldFEFlags );
1334 feclearexcept( FE_ALL_EXCEPT );
1338 std::cout << "FacesWithSizeMap" << std::endl;
1339 FacesWithSizeMap.Statistics(std::cout);
1340 std::cout << "EdgesWithSizeMap" << std::endl;
1341 EdgesWithSizeMap.Statistics(std::cout);
1342 std::cout << "VerticesWithSizeMap" << std::endl;
1343 VerticesWithSizeMap.Statistics(std::cout);
1344 std::cout << "FacesWithEnforcedVertices" << std::endl;
1345 FacesWithEnforcedVertices.Statistics(std::cout);
1351 //=============================================================================
1355 //=============================================================================
1357 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1358 const TopoDS_Edge edge = TopoDS::Edge(ed);
1360 gp_Pnt pnt(node->X(), node->Y(), node->Z());
1362 Standard_Real p0 = 0.0;
1363 Standard_Real p1 = 1.0;
1364 TopLoc_Location loc;
1365 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
1367 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
1368 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
1371 if ( proj.NbPoints() > 0 )
1373 pa = (double)proj.LowerDistanceParameter();
1374 // Issue 0020656. Move node if it is too far from edge
1375 gp_Pnt curve_pnt = curve->Value( pa );
1376 double dist2 = pnt.SquareDistance( curve_pnt );
1377 double tol = BRep_Tool::Tolerance( edge );
1378 if ( 1e-12 < dist2 && dist2 <= 2*tol*tol ) // large enough and within tolerance
1380 curve_pnt.Transform( loc );
1381 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
1384 // GProp_GProps LProps;
1385 // BRepGProp::LinearProperties(ed, LProps);
1386 // double lg = (double)LProps.Mass();
1388 meshDS->SetNodeOnEdge(node, edge, pa);
1391 //=============================================================================
1395 //=============================================================================
1397 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
1402 //=============================================================================
1406 //=============================================================================
1408 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
1413 //=============================================================================
1417 //=============================================================================
1419 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
1421 return hyp.SaveTo( save );
1424 //=============================================================================
1428 //=============================================================================
1430 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
1432 return hyp.LoadFrom( load );
1435 /* Curve definition function See cad_curv_t in file distene/cad.h for
1437 * NOTE : if when your CAD systems evaluates second
1438 * order derivatives it also computes first order derivatives and
1439 * function evaluation, you can optimize this example by making only
1440 * one CAD call and filling the necessary uv, dt, dtt arrays.
1442 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
1444 /* t is given. It contains the t (time) 1D parametric coordintaes
1445 of the point PreCAD/BLSurf is querying on the curve */
1447 /* user_data identifies the edge PreCAD/BLSurf is querying
1448 * (see cad_edge_new later in this example) */
1449 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
1452 /* BLSurf is querying the function evaluation */
1455 uv[0]=P.X(); uv[1]=P.Y();
1459 /* query for the first order derivatives */
1462 dt[0]=V1.X(); dt[1]=V1.Y();
1466 /* query for the second order derivatives */
1469 dtt[0]=V2.X(); dtt[1]=V2.Y();
1475 /* Surface definition function.
1476 * See cad_surf_t in file distene/cad.h for more information.
1477 * NOTE : if when your CAD systems evaluates second order derivatives it also
1478 * computes first order derivatives and function evaluation, you can optimize
1479 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
1482 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1483 real *duu, real *duv, real *dvv, void *user_data)
1485 /* uv[2] is given. It contains the u,v coordinates of the point
1486 * PreCAD/BLSurf is querying on the surface */
1488 /* user_data identifies the face PreCAD/BLSurf is querying (see
1489 * cad_face_new later in this example)*/
1490 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
1494 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
1495 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
1502 geometry->D1(uv[0],uv[1],P,D1U,D1V);
1503 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
1504 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
1507 if(duu && duv && dvv){
1511 gp_Vec D2U,D2V,D2UV;
1513 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
1514 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
1515 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
1516 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
1523 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
1526 if (my_u_min > uv[0]) {
1529 if (my_v_min > uv[1]) {
1532 if (my_u_max < uv[0]) {
1535 if (my_v_max < uv[1]) {
1540 if (FaceId2PythonSmp.count(face_id) != 0){
1541 PyObject * pyresult = NULL;
1542 PyObject* new_stderr = NULL;
1543 assert(Py_IsInitialized());
1544 PyGILState_STATE gstate;
1545 gstate = PyGILState_Ensure();
1546 pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],"(f,f)",uv[0],uv[1]);
1548 if ( pyresult == NULL){
1550 string err_description="";
1551 new_stderr = newPyStdOut(err_description);
1552 PySys_SetObject("stderr", new_stderr);
1554 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1555 Py_DECREF(new_stderr);
1556 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
1557 result = *((double*)user_data);
1560 result = PyFloat_AsDouble(pyresult);
1561 Py_DECREF(pyresult);
1564 //MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
1565 PyGILState_Release(gstate);
1568 *size = *((double*)user_data);
1573 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
1575 if (EdgeId2PythonSmp.count(edge_id) != 0){
1576 PyObject * pyresult = NULL;
1577 PyObject* new_stderr = NULL;
1578 assert(Py_IsInitialized());
1579 PyGILState_STATE gstate;
1580 gstate = PyGILState_Ensure();
1581 pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],"(f)",t);
1583 if ( pyresult == NULL){
1585 string err_description="";
1586 new_stderr = newPyStdOut(err_description);
1587 PySys_SetObject("stderr", new_stderr);
1589 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1590 Py_DECREF(new_stderr);
1591 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
1592 result = *((double*)user_data);
1595 result = PyFloat_AsDouble(pyresult);
1596 Py_DECREF(pyresult);
1599 PyGILState_Release(gstate);
1602 *size = *((double*)user_data);
1607 status_t size_on_vertex(integer point_id, real *size, void *user_data)
1609 if (VertexId2PythonSmp.count(point_id) != 0){
1610 PyObject * pyresult = NULL;
1611 PyObject* new_stderr = NULL;
1612 assert(Py_IsInitialized());
1613 PyGILState_STATE gstate;
1614 gstate = PyGILState_Ensure();
1615 pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],"");
1617 if ( pyresult == NULL){
1619 string err_description="";
1620 new_stderr = newPyStdOut(err_description);
1621 PySys_SetObject("stderr", new_stderr);
1623 PySys_SetObject("stderr", PySys_GetObject("__stderr__"));
1624 Py_DECREF(new_stderr);
1625 MESSAGE("Can't evaluate f()" << " error is " << err_description);
1626 result = *((double*)user_data);
1629 result = PyFloat_AsDouble(pyresult);
1630 Py_DECREF(pyresult);
1633 PyGILState_Release(gstate);
1636 *size = *((double*)user_data);
1642 * The following function will be called for PreCAD/BLSurf message
1643 * printing. See context_set_message_callback (later in this
1644 * template) for how to set user_data.
1646 status_t message_cb(message_t *msg, void *user_data)
1648 integer errnumber = 0;
1650 message_get_number(msg, &errnumber);
1651 message_get_description(msg, &desc);
1652 if ( errnumber < 0 ) {
1653 string * error = (string*)user_data;
1654 // if ( !error->empty() )
1656 // remove ^A from the tail
1657 int len = strlen( desc );
1658 while (len > 0 && desc[len-1] != '\n')
1660 error->append( desc, len );
1663 std::cout << desc << std::endl;
1668 /* This is the interrupt callback. PreCAD/BLSurf will call this
1669 * function regularily. See the file distene/interrupt.h
1671 status_t interrupt_cb(integer *interrupt_status, void *user_data)
1673 integer you_want_to_continue = 1;
1675 if(you_want_to_continue)
1676 *interrupt_status = INTERRUPT_CONTINUE;
1677 else /* you want to stop BLSurf */
1678 *interrupt_status = INTERRUPT_STOP;
1683 //=============================================================================
1687 //=============================================================================
1688 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
1689 const TopoDS_Shape& aShape,
1690 MapShapeNbElems& aResMap)
1692 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
1693 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
1694 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
1695 //double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
1696 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
1698 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
1699 _phySize = _hypothesis->GetPhySize();
1700 //_geometricMesh = (int) hyp->GetGeometricMesh();
1701 //_angleMeshS = hyp->GetAngleMeshS();
1702 _angleMeshC = _hypothesis->GetAngleMeshC();
1705 bool IsQuadratic = false;
1710 TopTools_DataMapOfShapeInteger EdgesMap;
1711 double fullLen = 0.0;
1712 double fullNbSeg = 0;
1713 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
1714 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
1715 if( EdgesMap.IsBound(E) )
1717 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
1718 double aLen = SMESH_Algo::EdgeLength(E);
1721 if(_physicalMesh==1) {
1722 nb1d = (int)( aLen/_phySize + 1 );
1727 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
1728 double fullAng = 0.0;
1729 double dp = (l-f)/200;
1734 for(int j=2; j<=200; j++) {
1737 fullAng += fabs(V1.Angle(V2));
1741 nb1d = (int)( fullAng/_angleMeshC + 1 );
1744 std::vector<int> aVec(SMDSEntity_Last);
1745 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1746 if( IsQuadratic > 0 ) {
1747 aVec[SMDSEntity_Node] = 2*nb1d - 1;
1748 aVec[SMDSEntity_Quad_Edge] = nb1d;
1751 aVec[SMDSEntity_Node] = nb1d - 1;
1752 aVec[SMDSEntity_Edge] = nb1d;
1754 aResMap.insert(std::make_pair(sm,aVec));
1755 EdgesMap.Bind(E,nb1d);
1757 double ELen = fullLen/fullNbSeg;
1761 // try to evaluate as in MEFISTO
1762 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
1763 TopoDS_Face F = TopoDS::Face( exp.Current() );
1764 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
1766 BRepGProp::SurfaceProperties(F,G);
1767 double anArea = G.Mass();
1769 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
1770 nb1d += EdgesMap.Find(exp1.Current());
1772 int nbFaces = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
1773 int nbNodes = (int) ( ( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
1774 std::vector<int> aVec(SMDSEntity_Last);
1775 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1777 int nb1d_in = (nbFaces*3 - nb1d) / 2;
1778 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
1779 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
1782 aVec[SMDSEntity_Node] = nbNodes;
1783 aVec[SMDSEntity_Triangle] = nbFaces;
1785 aResMap.insert(std::make_pair(sm,aVec));
1792 BRepGProp::VolumeProperties(aShape,G);
1793 double aVolume = G.Mass();
1794 double tetrVol = 0.1179*ELen*ELen*ELen;
1795 int nbVols = int(aVolume/tetrVol);
1796 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
1797 std::vector<int> aVec(SMDSEntity_Last);
1798 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
1800 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
1801 aVec[SMDSEntity_Quad_Tetra] = nbVols;
1804 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
1805 aVec[SMDSEntity_Tetra] = nbVols;
1807 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
1808 aResMap.insert(std::make_pair(sm,aVec));
1813 //=============================================================================
1815 * Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
1816 * Following line was added:
1817 * myExtrem.Perform(P);
1819 //=============================================================================
1820 void BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
1821 const TopoDS_Face& face,
1823 const Standard_Real Tol)
1825 //-- Voir BRepExtrema_ExtPF.cxx
1826 BRepAdaptor_Surface Surf(face);
1827 Standard_Real U1, U2, V1, V2;
1828 BRepTools::UVBounds(face, U1, U2, V1, V2);
1829 Extrema_ExtPS myExtrem;
1830 myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
1831 myExtrem.Perform(P);
1832 //----------------------------------------------------------
1833 //-- On cherche le point le plus proche , PUIS
1834 //-- On le classifie.
1835 Standard_Integer nbv = 0; // xpu
1836 Standard_Real MaxDist = RealLast();
1837 Standard_Integer indice = 0;
1838 if(myExtrem.IsDone()) {
1839 nbv = myExtrem.NbExt();
1840 for (Standard_Integer i = 1; i <= nbv; i++) {
1841 Standard_Real d = myExtrem.Value(i);
1851 Standard_Real U1,U2;
1852 myExtrem.Point(indice).Parameter(U1, U2);
1853 Puv.SetCoord(U1, U2);
1854 fc->Perform(face, Puv, Tol);
1857 fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
1858 //-- le resultat est TopAbs_OUT;