1 // Copyright (C) 2007-2012 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"
28 #include "BLSURFPlugin_Attractor.hxx"
31 #include <distene/api.h>
32 #include <distene/blsurf.h>
33 #include <distene/precad.h>
36 #include <structmember.h>
39 #include <Basics_Utils.hxx>
40 #include <Basics_OCCTVersion.hxx>
42 #include <SMESHDS_Group.hxx>
43 #include <SMESH_Gen.hxx>
44 #include <SMESH_Group.hxx>
45 #include <SMESH_Mesh.hxx>
46 #include <SMESH_MeshEditor.hxx>
47 #include <SMESH_MesherHelper.hxx>
48 #include <StdMeshers_FaceSide.hxx>
50 #include <utilities.h>
58 // OPENCASCADE includes
59 #include <BRep_Tool.hxx>
61 #include <TopExp_Explorer.hxx>
63 #include <NCollection_Map.hxx>
65 #include <Geom_Surface.hxx>
66 #include <Handle_Geom_Surface.hxx>
67 #include <Geom2d_Curve.hxx>
68 #include <Handle_Geom2d_Curve.hxx>
69 #include <Geom_Curve.hxx>
70 #include <Handle_Geom_Curve.hxx>
71 #include <Handle_AIS_InteractiveObject.hxx>
72 #include <TopoDS_Vertex.hxx>
73 #include <TopoDS_Edge.hxx>
74 #include <TopoDS_Wire.hxx>
75 #include <TopoDS_Face.hxx>
77 #include <gp_Pnt2d.hxx>
79 #include <TopTools_IndexedMapOfShape.hxx>
80 #include <TopoDS_Shape.hxx>
81 #include <BRep_Builder.hxx>
82 #include <BRepBuilderAPI_MakeVertex.hxx>
83 #include <BRepTools.hxx>
85 #include <TopTools_DataMapOfShapeInteger.hxx>
86 #include <GProp_GProps.hxx>
87 #include <BRepGProp.hxx>
93 #include <Standard_ErrorHandler.hxx>
94 #include <GeomAPI_ProjectPointOnCurve.hxx>
95 #include <GeomAPI_ProjectPointOnSurf.hxx>
98 // #include <BRepClass_FaceClassifier.hxx>
99 #include <TopTools_MapOfShape.hxx>
101 /* ==================================
102 * =========== PYTHON ==============
103 * ==================================*/
112 PyStdOut_dealloc(PyStdOut *self)
118 PyStdOut_write(PyStdOut *self, PyObject *args)
122 if (!PyArg_ParseTuple(args, "t#:write",&c, &l))
126 *(self->out)=*(self->out)+c;
132 static PyMethodDef PyStdOut_methods[] = {
133 {"write", (PyCFunction)PyStdOut_write, METH_VARARGS,
134 PyDoc_STR("write(string) -> None")},
135 {NULL, NULL} /* sentinel */
138 static PyMemberDef PyStdOut_memberlist[] = {
139 {(char*)"softspace", T_INT, offsetof(PyStdOut, softspace), 0,
140 (char*)"flag indicating that a space needs to be printed; used by print"},
141 {NULL} /* Sentinel */
144 static PyTypeObject PyStdOut_Type = {
145 /* The ob_type field must be initialized in the module init function
146 * to be portable to Windows without using C++. */
147 PyObject_HEAD_INIT(NULL)
150 sizeof(PyStdOut), /*tp_basicsize*/
153 (destructor)PyStdOut_dealloc, /*tp_dealloc*/
160 0, /*tp_as_sequence*/
165 PyObject_GenericGetAttr, /*tp_getattro*/
166 /* softspace is writable: we must supply tp_setattro */
167 PyObject_GenericSetAttr, /* tp_setattro */
169 Py_TPFLAGS_DEFAULT, /*tp_flags*/
173 0, /*tp_richcompare*/
174 0, /*tp_weaklistoffset*/
177 PyStdOut_methods, /*tp_methods*/
178 PyStdOut_memberlist, /*tp_members*/
192 PyObject * newPyStdOut( std::string& out )
195 self = PyObject_New(PyStdOut, &PyStdOut_Type);
200 return (PyObject*)self;
204 ////////////////////////END PYTHON///////////////////////////
206 //////////////////MY MAPS////////////////////////////////////////
207 TopTools_IndexedMapOfShape FacesWithSizeMap;
208 std::map<int,string> FaceId2SizeMap;
209 TopTools_IndexedMapOfShape EdgesWithSizeMap;
210 std::map<int,string> EdgeId2SizeMap;
211 TopTools_IndexedMapOfShape VerticesWithSizeMap;
212 std::map<int,string> VertexId2SizeMap;
214 std::map<int,PyObject*> FaceId2PythonSmp;
215 std::map<int,PyObject*> EdgeId2PythonSmp;
216 std::map<int,PyObject*> VertexId2PythonSmp;
218 std::map<int,std::vector<double> > FaceId2AttractorCoords;
219 std::map<int,BLSURFPlugin_Attractor*> FaceId2ClassAttractor;
220 std::map<int,BLSURFPlugin_Attractor*> FaceIndex2ClassAttractor;
222 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
223 std::map< int, BLSURFPlugin_Hypothesis::TEnfVertexCoordsList > FaceId2EnforcedVertexCoords;
224 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexCoords > EnfVertexCoords2ProjVertex;
225 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList > EnfVertexCoords2EnfVertexList;
227 bool HasSizeMapOnFace=false;
228 bool HasSizeMapOnEdge=false;
229 bool HasSizeMapOnVertex=false;
230 //bool HasAttractorOnFace=false;
232 //=============================================================================
236 //=============================================================================
238 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId, int studyId,
240 : SMESH_2D_Algo(hypId, studyId, gen)
242 MESSAGE("BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF");
245 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
246 _compatibleHypothesis.push_back(BLSURFPlugin_Hypothesis::GetHypType());
247 _requireDiscreteBoundary = false;
248 _onlyUnaryInput = false;
250 _supportSubmeshes = true;
252 smeshGen_i = SMESH_Gen_i::GetSMESHGen();
253 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
254 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
256 MESSAGE("studyid = " << _studyId);
259 myStudy = aStudyMgr->GetStudyByID(_studyId);
261 MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
263 /* Initialize the Python interpreter */
264 assert(Py_IsInitialized());
265 PyGILState_STATE gstate;
266 gstate = PyGILState_Ensure();
269 main_mod = PyImport_AddModule("__main__");
272 main_dict = PyModule_GetDict(main_mod);
274 PyRun_SimpleString("from math import *");
275 PyGILState_Release(gstate);
277 FacesWithSizeMap.Clear();
278 FaceId2SizeMap.clear();
279 EdgesWithSizeMap.Clear();
280 EdgeId2SizeMap.clear();
281 VerticesWithSizeMap.Clear();
282 VertexId2SizeMap.clear();
283 FaceId2PythonSmp.clear();
284 EdgeId2PythonSmp.clear();
285 VertexId2PythonSmp.clear();
286 FaceId2AttractorCoords.clear();
287 FaceId2ClassAttractor.clear();
288 FaceIndex2ClassAttractor.clear();
289 FacesWithEnforcedVertices.Clear();
290 FaceId2EnforcedVertexCoords.clear();
291 EnfVertexCoords2ProjVertex.clear();
292 EnfVertexCoords2EnfVertexList.clear();
294 #ifdef WITH_SMESH_CANCEL_COMPUTE
295 _compute_canceled = false;
299 //=============================================================================
303 //=============================================================================
305 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
307 MESSAGE("BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF");
311 //=============================================================================
315 //=============================================================================
317 bool BLSURFPlugin_BLSURF::CheckHypothesis
319 const TopoDS_Shape& aShape,
320 SMESH_Hypothesis::Hypothesis_Status& aStatus)
324 list<const SMESHDS_Hypothesis*>::const_iterator itl;
325 const SMESHDS_Hypothesis* theHyp;
327 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
328 int nbHyp = hyps.size();
331 aStatus = SMESH_Hypothesis::HYP_OK;
332 return true; // can work with no hypothesis
336 theHyp = (*itl); // use only the first hypothesis
338 string hypName = theHyp->GetName();
340 if (hypName == "BLSURF_Parameters")
342 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
344 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
345 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
346 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
347 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
349 aStatus = SMESH_Hypothesis::HYP_OK;
352 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
354 return aStatus == SMESH_Hypothesis::HYP_OK;
357 //=============================================================================
359 * Pass parameters to BLSURF
361 //=============================================================================
363 inline std::string to_string(double d)
365 std::ostringstream o;
370 inline std::string to_string(int i)
372 std::ostringstream o;
377 double _smp_phy_size;
378 // #if BLSURF_VERSION_LONG >= "3.1.1"
379 // // sizemap_t *geo_sizemap_e, *geo_sizemap_f;
380 // sizemap_t *iso_sizemap_p, *iso_sizemap_e, *iso_sizemap_f;
381 // // sizemap_t *clean_geo_sizemap_e, *clean_geo_sizemap_f;
382 // sizemap_t *clean_iso_sizemap_p, *clean_iso_sizemap_e, *clean_iso_sizemap_f;
384 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
385 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
386 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
388 double my_u_min=1e6,my_v_min=1e6,my_u_max=-1e6,my_v_max=-1e6;
394 /////////////////////////////////////////////////////////
395 projectionPoint getProjectionPoint(const TopoDS_Face& face, const gp_Pnt& point)
397 projectionPoint myPoint;
398 Handle(Geom_Surface) surface = BRep_Tool::Surface(face);
399 GeomAPI_ProjectPointOnSurf projector( point, surface );
400 if ( !projector.IsDone() || projector.NbPoints()==0 )
401 throw "getProjectionPoint: Can't project";
403 Quantity_Parameter u,v;
404 projector.LowerDistanceParameters(u,v);
405 myPoint.uv = gp_XY(u,v);
406 gp_Pnt aPnt = projector.NearestPoint();
407 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
411 /////////////////////////////////////////////////////////
413 /////////////////////////////////////////////////////////
414 double getT(const TopoDS_Edge& edge, const gp_Pnt& point)
417 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, f,l);
418 GeomAPI_ProjectPointOnCurve projector( point, curve);
419 if ( projector.NbPoints() == 0 )
421 return projector.LowerDistanceParameter();
424 /////////////////////////////////////////////////////////
425 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
427 MESSAGE("BLSURFPlugin_BLSURF::entryToShape "<<entry );
428 GEOM::GEOM_Object_var aGeomObj;
429 TopoDS_Shape S = TopoDS_Shape();
430 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
431 SALOMEDS::GenericAttribute_var anAttr;
433 if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR")) {
434 SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
435 CORBA::String_var aVal = anIOR->Value();
436 CORBA::Object_var obj = myStudy->ConvertIORToObject(aVal);
437 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
439 if ( !aGeomObj->_is_nil() )
440 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
444 void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
446 BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
451 // Get the (u,v) values of the enforced vertex on the face
452 projectionPoint myPoint = getProjectionPoint(faceShape,aPnt);
454 MESSAGE("Enforced Vertex: " << aPnt.X() << ", " << aPnt.Y() << ", " << aPnt.Z());
455 MESSAGE("Projected Vertex: " << myPoint.xyz.X() << ", " << myPoint.xyz.Y() << ", " << myPoint.xyz.Z());
456 MESSAGE("Parametric coordinates: " << myPoint.uv.X() << ", " << myPoint.uv.Y() );
458 enf_coords.push_back(aPnt.X());
459 enf_coords.push_back(aPnt.Y());
460 enf_coords.push_back(aPnt.Z());
462 coords.push_back(myPoint.uv.X());
463 coords.push_back(myPoint.uv.Y());
464 coords.push_back(myPoint.xyz.X());
465 coords.push_back(myPoint.xyz.Y());
466 coords.push_back(myPoint.xyz.Z());
468 s_coords.push_back(myPoint.xyz.X());
469 s_coords.push_back(myPoint.xyz.Y());
470 s_coords.push_back(myPoint.xyz.Z());
472 // Save pair projected vertex / enf vertex
473 MESSAGE("Storing pair projected vertex / enf vertex:");
474 MESSAGE("("<< myPoint.xyz.X() << ", " << myPoint.xyz.Y() << ", " << myPoint.xyz.Z() <<") / (" << aPnt.X() << ", " << aPnt.Y() << ", " << aPnt.Z()<<")");
475 EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
476 MESSAGE("Group name is: \"" << enfVertex->grpName << "\"");
477 pair<BLSURFPlugin_Hypothesis::TEnfVertexList::iterator,bool> ret;
478 BLSURFPlugin_Hypothesis::TEnfVertexList::iterator it;
479 ret = EnfVertexCoords2EnfVertexList[s_coords].insert(enfVertex);
480 if (ret.second == false) {
482 (*it)->grpName = enfVertex->grpName;
486 if (! FacesWithEnforcedVertices.Contains(faceShape)) {
487 key = FacesWithEnforcedVertices.Add(faceShape);
490 key = FacesWithEnforcedVertices.FindIndex(faceShape);
493 // If a node is already created by an attractor, do not create enforced vertex
494 int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
495 bool sameAttractor = false;
496 if (attractorKey >= 0)
497 if (FaceId2AttractorCoords.count(attractorKey) > 0)
498 if (FaceId2AttractorCoords[attractorKey] == coords)
499 sameAttractor = true;
501 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
502 MESSAGE("Map of enf. vertex has key " << key)
503 MESSAGE("Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
505 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
507 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
508 MESSAGE("New Enf. vertex list size is: " << FaceId2EnforcedVertexCoords[key].size())
511 MESSAGE("Map of enf. vertex has not key " << key << ": creating it")
512 if (! sameAttractor) {
513 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList ens;
515 FaceId2EnforcedVertexCoords[key] = ens;
518 MESSAGE("An attractor node is already defined: I don't add the enforced vertex");
522 /////////////////////////////////////////////////////////
523 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
525 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
528 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
530 for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
531 enfVertex = *enfVertexListIt;
532 // Case of manual coords
533 if (enfVertex->coords.size() != 0) {
534 aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
535 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
538 // Case of geom vertex coords
539 if (enfVertex->geomEntry != "") {
540 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
541 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
542 if (GeomType == TopAbs_VERTEX){
543 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(GeomShape));
544 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
547 if (GeomType == TopAbs_COMPOUND){
548 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
549 if (it.Value().ShapeType() == TopAbs_VERTEX){
550 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
551 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
559 /////////////////////////////////////////////////////////
560 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction)
562 MESSAGE("Attractor function: "<< AttractorFunction);
563 double xa, ya, za; // Coordinates of attractor point
564 double a, b; // Attractor parameter
566 bool createNode=false; // To create a node on attractor projection
568 const char *sep = ";";
569 // atIt->second has the following pattern:
570 // ATTRACTOR(xa;ya;za;a;b;True|False;d)
572 // xa;ya;za : coordinates of attractor
573 // a : desired size on attractor
574 // b : distance of influence of attractor
575 // d : distance until which the size remains constant
577 // We search the parameters in the string
579 pos1 = AttractorFunction.find(sep);
580 if (pos1!=string::npos)
581 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
583 pos2 = AttractorFunction.find(sep, pos1+1);
584 if (pos2!=string::npos) {
585 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
589 pos2 = AttractorFunction.find(sep, pos1+1);
590 if (pos2!=string::npos) {
591 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
595 pos2 = AttractorFunction.find(sep, pos1+1);
596 if (pos2!=string::npos) {
597 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
601 pos2 = AttractorFunction.find(sep, pos1+1);
602 if (pos2!=string::npos) {
603 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
607 pos2 = AttractorFunction.find(sep, pos1+1);
608 if (pos2!=string::npos) {
609 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
610 MESSAGE("createNode: " << createNodeStr);
611 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
615 pos2 = AttractorFunction.find(")");
616 if (pos2!=string::npos) {
617 d = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
620 // Get the (u,v) values of the attractor on the face
621 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
622 gp_XY uvPoint = myPoint.uv;
623 gp_XYZ xyzPoint = myPoint.xyz;
624 Standard_Real u0 = uvPoint.X();
625 Standard_Real v0 = uvPoint.Y();
626 Standard_Real x0 = xyzPoint.X();
627 Standard_Real y0 = xyzPoint.Y();
628 Standard_Real z0 = xyzPoint.Z();
629 std::vector<double> coords;
630 coords.push_back(u0);
631 coords.push_back(v0);
632 coords.push_back(x0);
633 coords.push_back(y0);
634 coords.push_back(z0);
635 // We construct the python function
636 ostringstream attractorFunctionStream;
637 attractorFunctionStream << "def f(u,v): return ";
638 attractorFunctionStream << _smp_phy_size << "-(" << _smp_phy_size <<"-" << a << ")";
639 //attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
640 // rnc: make possible to keep the size constant until
641 // a defined distance. Distance is expressed as the positiv part
642 // of r-d where r is the distance to (u0,v0)
643 attractorFunctionStream << "*exp(-(0.5*(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"+abs(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"))/(" << b << "))**2)";
645 MESSAGE("Python function for attractor:" << std::endl << attractorFunctionStream.str());
648 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
649 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
652 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
654 FaceId2SizeMap[key] =attractorFunctionStream.str();
656 MESSAGE("Creating node on ("<<x0<<","<<y0<<","<<z0<<")");
657 FaceId2AttractorCoords[key] = coords;
659 // // Test for new attractors
660 // gp_Pnt myP(xyzPoint);
661 // TopoDS_Vertex myV = BRepBuilderAPI_MakeVertex(myP);
662 // BLSURFPlugin_Attractor myAttractor(TopoDS::Face(GeomShape),myV,200);
663 // myAttractor.SetParameters(a, _smp_phy_size, b, d);
664 // myAttractor.SetType(1);
665 // FaceId2ClassAttractor[key] = myAttractor;
666 // if(FaceId2ClassAttractor[key].GetFace().IsNull()){
667 // MESSAGE("face nulle ");
670 // MESSAGE("face OK");
672 // if (FaceId2ClassAttractor[key].GetAttractorShape().IsNull()){
673 // MESSAGE("pas de point");
676 // MESSAGE("point OK");
679 /////////////////////////////////////////////////////////
681 void BLSURFPlugin_BLSURF::SetParameters(
682 // #if BLSURF_VERSION_LONG >= "3.1.1"
685 const BLSURFPlugin_Hypothesis* hyp,
686 blsurf_session_t * bls,
687 precad_session_t * pcs,
693 // Clear map so that it is not stored in the algorithm with old enforced vertices in it
694 EnfVertexCoords2EnfVertexList.clear();
696 int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
697 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
698 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
699 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
700 double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
701 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
702 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
703 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
704 bool _decimesh = BLSURFPlugin_Hypothesis::GetDefaultDecimesh();
705 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
708 int _precadMergeEdges = BLSURFPlugin_Hypothesis::GetDefaultPreCADMergeEdges();
709 int _precadRemoveNanoEdges = BLSURFPlugin_Hypothesis::GetDefaultPreCADRemoveNanoEdges();
710 int _precadDiscardInput = BLSURFPlugin_Hypothesis::GetDefaultPreCADDiscardInput();
711 double _precadEpsNano = BLSURFPlugin_Hypothesis::GetDefaultPreCADEpsNano();
715 MESSAGE("BLSURFPlugin_BLSURF::SetParameters");
716 _topology = (int) hyp->GetTopology();
717 _physicalMesh = (int) hyp->GetPhysicalMesh();
718 _phySize = hyp->GetPhySize();
719 _geometricMesh = (int) hyp->GetGeometricMesh();
720 _angleMeshS = hyp->GetAngleMeshS();
721 _angleMeshC = hyp->GetAngleMeshC();
722 _gradation = hyp->GetGradation();
723 _quadAllowed = hyp->GetQuadAllowed();
724 _decimesh = hyp->GetDecimesh();
725 _verb = hyp->GetVerbosity();
726 if ( hyp->GetPhyMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
727 blsurf_set_param(bls, "hphymin", to_string(hyp->GetPhyMin()).c_str());
728 if ( hyp->GetPhyMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
729 blsurf_set_param(bls, "hphymax", to_string(hyp->GetPhyMax()).c_str());
730 if ( hyp->GetGeoMin() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
731 blsurf_set_param(bls, "hgeomin", to_string(hyp->GetGeoMin()).c_str());
732 if ( hyp->GetGeoMax() != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
733 blsurf_set_param(bls, "hgeomax", to_string(hyp->GetGeoMax()).c_str());
735 const BLSURFPlugin_Hypothesis::TOptionValues & opts = hyp->GetOptionValues();
736 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
737 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt )
738 if ( !opIt->second.empty() ) {
739 MESSAGE("blsurf_set_param(): " << opIt->first << " = " << opIt->second);
740 blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
743 const BLSURFPlugin_Hypothesis::TOptionValues & preCADopts = hyp->GetPreCADOptionValues();
744 for ( opIt = preCADopts.begin(); opIt != preCADopts.end(); ++opIt )
745 if ( !opIt->second.empty() ) {
746 if (_topology == BLSURFPlugin_Hypothesis::PreCAD) {
747 MESSAGE("precad_set_param(): " << opIt->first << " = " << opIt->second);
748 blsurf_set_param(bls, opIt->first.c_str(), opIt->second.c_str());
753 _precadMergeEdges = hyp->GetPreCADMergeEdges();
754 _precadRemoveNanoEdges = hyp->GetPreCADRemoveNanoEdges();
755 _precadDiscardInput = hyp->GetPreCADDiscardInput();
756 _precadEpsNano = hyp->GetPreCADEpsNano();
759 //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
760 // GetDefaultPhySize() sometimes leads to computation failure
761 _phySize = mesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
762 MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
766 if (_topology == BLSURFPlugin_Hypothesis::PreCAD) {
768 precad_set_param(pcs, "verbose", to_string(_verb).c_str());
769 precad_set_param(pcs, "merge_edges", _precadMergeEdges ? "1" : "0");
770 precad_set_param(pcs, "remove_nano_edges", _precadRemoveNanoEdges ? "1" : "0");
771 precad_set_param(pcs, "discard_input_topology", _precadDiscardInput ? "1" : "0");
772 if ( _precadEpsNano != ::BLSURFPlugin_Hypothesis::undefinedDouble() )
773 precad_set_param(pcs, "eps_nano", to_string(_precadEpsNano).c_str());
776 _smp_phy_size = _phySize;
777 blsurf_set_param(bls, "topo_points", _topology == BLSURFPlugin_Hypothesis::Process || _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
778 blsurf_set_param(bls, "topo_curves", _topology == BLSURFPlugin_Hypothesis::Process || _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
779 blsurf_set_param(bls, "topo_project", _topology == BLSURFPlugin_Hypothesis::Process || _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
780 blsurf_set_param(bls, "clean_boundary", _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
781 blsurf_set_param(bls, "close_boundary", _topology == BLSURFPlugin_Hypothesis::Process2 ? "1" : "0");
782 blsurf_set_param(bls, "hphy_flag", to_string(_physicalMesh).c_str());
783 blsurf_set_param(bls, "hphydef", to_string(_phySize).c_str());
784 blsurf_set_param(bls, "hgeo_flag", to_string(_geometricMesh).c_str());
785 blsurf_set_param(bls, "relax_size", _decimesh ? "0": to_string(_geometricMesh).c_str());
786 blsurf_set_param(bls, "angle_meshs", to_string(_angleMeshS).c_str());
787 blsurf_set_param(bls, "angle_meshc", to_string(_angleMeshC).c_str());
788 blsurf_set_param(bls, "gradation", to_string(_gradation).c_str());
789 blsurf_set_param(bls, "patch_independent", _decimesh ? "1" : "0");
790 blsurf_set_param(bls, "element", _quadAllowed ? "q1.0" : "p1");
791 blsurf_set_param(bls, "verb", to_string(_verb).c_str());
793 if (_physicalMesh == BLSURFPlugin_Hypothesis::SizeMap){
794 TopoDS_Shape GeomShape;
795 TopoDS_Shape AttShape;
796 TopAbs_ShapeEnum GeomType;
798 // Standard Size Maps
800 MESSAGE("Setting a Size Map");
801 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
802 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
803 for ( ; smIt != sizeMaps.end(); ++smIt ) {
804 if ( !smIt->second.empty() ) {
805 MESSAGE("blsurf_set_sizeMap(): " << smIt->first << " = " << smIt->second);
806 GeomShape = entryToShape(smIt->first);
807 GeomType = GeomShape.ShapeType();
808 MESSAGE("Geomtype is " << GeomType);
811 if (GeomType == TopAbs_COMPOUND){
812 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
814 if (it.Value().ShapeType() == TopAbs_FACE){
815 HasSizeMapOnFace = true;
816 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
817 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
820 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
821 // MESSAGE("Face with key " << key << " already in map");
823 FaceId2SizeMap[key] = smIt->second;
826 if (it.Value().ShapeType() == TopAbs_EDGE){
827 HasSizeMapOnEdge = true;
828 HasSizeMapOnFace = true;
829 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
830 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
833 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
834 // MESSAGE("Edge with key " << key << " already in map");
836 EdgeId2SizeMap[key] = smIt->second;
839 if (it.Value().ShapeType() == TopAbs_VERTEX){
840 HasSizeMapOnVertex = true;
841 HasSizeMapOnEdge = true;
842 HasSizeMapOnFace = true;
843 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
844 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
847 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
848 MESSAGE("Group of vertices with key " << key << " already in map");
850 MESSAGE("Group of vertices with key " << key << " has a size map: " << smIt->second);
851 VertexId2SizeMap[key] = smIt->second;
856 if (GeomType == TopAbs_FACE){
857 HasSizeMapOnFace = true;
858 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
859 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
862 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
863 // MESSAGE("Face with key " << key << " already in map");
865 FaceId2SizeMap[key] = smIt->second;
868 if (GeomType == TopAbs_EDGE){
869 HasSizeMapOnEdge = true;
870 HasSizeMapOnFace = true;
871 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
872 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
875 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
876 // MESSAGE("Edge with key " << key << " already in map");
878 EdgeId2SizeMap[key] = smIt->second;
881 if (GeomType == TopAbs_VERTEX){
882 HasSizeMapOnVertex = true;
883 HasSizeMapOnEdge = true;
884 HasSizeMapOnFace = true;
885 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
886 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
889 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
890 MESSAGE("Vertex with key " << key << " already in map");
892 MESSAGE("Vertex with key " << key << " has a size map: " << smIt->second);
893 VertexId2SizeMap[key] = smIt->second;
901 // TODO appeler le constructeur des attracteurs directement ici
902 MESSAGE("Setting Attractors");
903 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
904 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
905 for ( ; atIt != attractors.end(); ++atIt ) {
906 if ( !atIt->second.empty() ) {
907 MESSAGE("blsurf_set_attractor(): " << atIt->first << " = " << atIt->second);
908 GeomShape = entryToShape(atIt->first);
909 GeomType = GeomShape.ShapeType();
911 if (GeomType == TopAbs_COMPOUND){
912 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
913 if (it.Value().ShapeType() == TopAbs_FACE){
914 HasSizeMapOnFace = true;
915 createAttractorOnFace(it.Value(), atIt->second);
920 if (GeomType == TopAbs_FACE){
921 HasSizeMapOnFace = true;
922 createAttractorOnFace(GeomShape, atIt->second);
925 if (GeomType == TopAbs_EDGE){
926 HasSizeMapOnEdge = true;
927 HasSizeMapOnFace = true;
928 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
930 if (GeomType == TopAbs_VERTEX){
931 HasSizeMapOnVertex = true;
932 HasSizeMapOnEdge = true;
933 HasSizeMapOnFace = true;
934 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
941 // temporary commented out for testing
943 // - Fill in the BLSURFPlugin_Hypothesis::TAttractorMap map in the hypothesis
944 // - Uncomment and complete this part to construct the attractors from the attractor shape and the passed parameters on each face of the map
945 // - To do this use the public methodss: SetParameters(several double parameters) and SetType(int type)
947 // - Construct the attractors with an empty dist. map in the hypothesis
948 // - build the map here for each face with an attractor set and only if the attractor shape as changed since the last call to _buildmap()
949 // -> define a bool _mapbuilt in the class that is set to false by default and set to true when calling _buildmap() OK
951 const BLSURFPlugin_Hypothesis::TAttractorMap class_attractors = BLSURFPlugin_Hypothesis::GetClassAttractorEntries(hyp);
953 BLSURFPlugin_Hypothesis::TAttractorMap::const_iterator AtIt = class_attractors.begin();
954 for ( ; AtIt != class_attractors.end(); ++AtIt ) {
955 if ( !AtIt->second->Empty() ) {
956 // MESSAGE("blsurf_set_attractor(): " << AtIt->first << " = " << AtIt->second);
957 GeomShape = entryToShape(AtIt->first);
958 AttShape = AtIt->second->GetAttractorShape();
959 GeomType = GeomShape.ShapeType();
961 // if (GeomType == TopAbs_COMPOUND){
962 // for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
963 // if (it.Value().ShapeType() == TopAbs_FACE){
964 // HasAttractorOnFace = true;
965 // myAttractor = BLSURFPluginAttractor(GeomShape, AttShape);
970 if (GeomType == TopAbs_FACE
971 && (AttShape.ShapeType() == TopAbs_VERTEX
972 || AttShape.ShapeType() == TopAbs_EDGE
973 || AttShape.ShapeType() == TopAbs_WIRE
974 || AttShape.ShapeType() == TopAbs_COMPOUND) ){
975 HasSizeMapOnFace = true;
977 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape)) ) {
978 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
981 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
982 // MESSAGE("Face with key " << key << " already in map");
985 FaceId2ClassAttractor[key] = AtIt->second;
988 MESSAGE("Wrong shape type !!")
998 MESSAGE("Setting Enforced Vertices");
999 const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
1000 BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
1001 for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
1002 if ( !enfIt->second.empty() ) {
1003 GeomShape = entryToShape(enfIt->first);
1004 GeomType = GeomShape.ShapeType();
1006 if (GeomType == TopAbs_COMPOUND){
1007 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1008 if (it.Value().ShapeType() == TopAbs_FACE){
1009 HasSizeMapOnFace = true;
1010 createEnforcedVertexOnFace(it.Value(), enfIt->second);
1015 if (GeomType == TopAbs_FACE){
1016 HasSizeMapOnFace = true;
1017 createEnforcedVertexOnFace(GeomShape, enfIt->second);
1022 // Internal vertices
1023 bool useInternalVertexAllFaces = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFaces(hyp);
1024 if (useInternalVertexAllFaces) {
1025 std::string grpName = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFacesGroup(hyp);
1026 MESSAGE("Setting Internal Enforced Vertices");
1027 GeomShape = mesh.GetShapeToMesh();
1029 TopExp_Explorer exp (GeomShape, TopAbs_FACE);
1030 for (; exp.More(); exp.Next()){
1031 MESSAGE("Iterating shapes. Shape type is " << exp.Current().ShapeType());
1032 TopExp_Explorer exp_face (exp.Current(), TopAbs_VERTEX, TopAbs_EDGE);
1033 for (; exp_face.More(); exp_face.Next())
1035 // Get coords of vertex
1036 // Check if current coords is already in enfVertexList
1037 // If coords not in enfVertexList, add new enfVertex
1038 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(exp_face.Current()));
1039 MESSAGE("Found vertex on face at " << aPnt.X() <<", "<<aPnt.Y()<<", "<<aPnt.Z());
1040 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex = new BLSURFPlugin_Hypothesis::TEnfVertex();
1041 enfVertex->coords.push_back(aPnt.X());
1042 enfVertex->coords.push_back(aPnt.Y());
1043 enfVertex->coords.push_back(aPnt.Z());
1044 enfVertex->name = "";
1045 enfVertex->faceEntries.clear();
1046 enfVertex->geomEntry = "";
1047 enfVertex->grpName = grpName;
1048 enfVertex->vertex = TopoDS::Vertex( exp_face.Current() );
1049 _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()), aPnt, enfVertex);
1050 HasSizeMapOnFace = true;
1055 MESSAGE("Setting Size Map on FACES ");
1056 // #if BLSURF_VERSION_LONG < "3.1.1"
1057 blsurf_data_set_sizemap_iso_cad_face(bls, size_on_surface, &_smp_phy_size);
1060 // iso_sizemap_f = sizemap_new(c, distene_sizemap_type_iso_cad_face, (void *)size_on_surface, NULL);
1062 // clean_iso_sizemap_f = sizemap_new(c, distene_sizemap_type_iso_cad_face, (void *)size_on_surface, NULL);
1065 if (HasSizeMapOnEdge){
1066 MESSAGE("Setting Size Map on EDGES ");
1067 // #if BLSURF_VERSION_LONG < "3.1.1"
1068 blsurf_data_set_sizemap_iso_cad_edge(bls, size_on_edge, &_smp_phy_size);
1071 // iso_sizemap_e = sizemap_new(c, distene_sizemap_type_iso_cad_edge, (void *)size_on_edge, NULL);
1073 // clean_iso_sizemap_e = sizemap_new(c, distene_sizemap_type_iso_cad_edge, (void *)size_on_edge, NULL);
1076 if (HasSizeMapOnVertex){
1077 MESSAGE("Setting Size Map on VERTICES ");
1078 // #if BLSURF_VERSION_LONG < "3.1.1"
1079 blsurf_data_set_sizemap_iso_cad_point(bls, size_on_vertex, &_smp_phy_size);
1082 // iso_sizemap_p = sizemap_new(c, distene_sizemap_type_iso_cad_point, (void *)size_on_vertex, NULL);
1084 // clean_iso_sizemap_p = sizemap_new(c, distene_sizemap_type_iso_cad_point, (void *)size_on_vertex, NULL);
1093 * \brief Class correctly terminating usage of BLSURF library at destruction
1095 class BLSURF_Cleaner
1098 blsurf_session_t* _bls;
1102 BLSURF_Cleaner(context_t * ctx,
1103 blsurf_session_t* bls,
1114 blsurf_session_delete(_bls);
1116 // #if BLSURF_VERSION_LONG >= "3.1.1"
1117 // // if(geo_sizemap_e)
1118 // // distene_sizemap_delete(geo_sizemap_e);
1119 // // if(geo_sizemap_f)
1120 // // distene_sizemap_delete(geo_sizemap_f);
1121 // if(iso_sizemap_p)
1122 // distene_sizemap_delete(iso_sizemap_p);
1123 // if(iso_sizemap_e)
1124 // distene_sizemap_delete(iso_sizemap_e);
1125 // if(iso_sizemap_f)
1126 // distene_sizemap_delete(iso_sizemap_f);
1128 // // if(clean_geo_sizemap_e)
1129 // // distene_sizemap_delete(clean_geo_sizemap_e);
1130 // // if(clean_geo_sizemap_f)
1131 // // distene_sizemap_delete(clean_geo_sizemap_f);
1132 // if(clean_iso_sizemap_p)
1133 // distene_sizemap_delete(clean_iso_sizemap_p);
1134 // if(clean_iso_sizemap_e)
1135 // distene_sizemap_delete(clean_iso_sizemap_e);
1136 // if(clean_iso_sizemap_f)
1137 // distene_sizemap_delete(clean_iso_sizemap_f);
1142 context_delete(_ctx);
1147 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
1148 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1149 real *duu, real *duv, real *dvv, void *user_data);
1150 status_t message_cb(message_t *msg, void *user_data);
1151 status_t interrupt_cb(integer *interrupt_status, void *user_data);
1153 //=============================================================================
1157 //=============================================================================
1159 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape) {
1161 MESSAGE("BLSURFPlugin_BLSURF::Compute");
1163 // Fix problem with locales
1164 Kernel_Utils::Localizer aLocalizer;
1166 /* create a distene context (generic object) */
1167 status_t status = STATUS_ERROR;
1169 context_t *ctx = context_new();
1171 /* Set the message callback in the working context */
1172 context_set_message_callback(ctx, message_cb, &_comment);
1173 #ifdef WITH_SMESH_CANCEL_COMPUTE
1174 _compute_canceled = false;
1175 context_set_interrupt_callback(ctx, interrupt_cb, this);
1177 context_set_interrupt_callback(ctx, interrupt_cb, NULL);
1180 /* create the CAD object we will work on. It is associated to the context ctx. */
1181 cad_t *c = cad_new(ctx);
1182 dcad_t *dcad = dcad_new(c);
1184 FacesWithSizeMap.Clear();
1185 FaceId2SizeMap.clear();
1186 FaceId2ClassAttractor.clear();
1187 FaceIndex2ClassAttractor.clear();
1188 EdgesWithSizeMap.Clear();
1189 EdgeId2SizeMap.clear();
1190 VerticesWithSizeMap.Clear();
1191 VertexId2SizeMap.clear();
1193 SMESH_MesherHelper helper( aMesh );
1194 // do not call helper.IsQuadraticSubMesh() because submeshes
1195 // may be cleaned and helper.myTLinkNodeMap gets invalid in such a case
1196 const bool haveQudraticSubMesh = SMESH_MesherHelper( aMesh ).IsQuadraticSubMesh( aShape );
1197 helper.SetIsQuadratic( haveQudraticSubMesh );
1198 bool needMerge = false;
1199 set< SMESH_subMesh* > edgeSubmeshes;
1200 set< SMESH_subMesh* >& mergeSubmeshes = edgeSubmeshes;
1202 /* Now fill the CAD object with data from your CAD
1203 * environement. This is the most complex part of a successfull
1208 // If user requests it, send the CAD through Distene preprocessor : PreCAD
1209 cad_t *cleanc = NULL;
1210 precad_session_t *pcs = precad_session_new(ctx);
1211 precad_data_set_cad(pcs, c);
1213 blsurf_session_t *bls = blsurf_session_new(ctx);
1215 // an object that correctly deletes all blsurf objects at destruction
1216 BLSURF_Cleaner cleaner( ctx,bls,c,dcad );
1218 MESSAGE("BEGIN SetParameters");
1219 bool use_precad = false;
1221 // #if BLSURF_VERSION_LONG >= "3.1.1"
1224 _hypothesis, bls, pcs, aMesh, &use_precad);
1225 MESSAGE("END SetParameters");
1227 // needed to prevent the opencascade memory managmement from freeing things
1228 vector<Handle(Geom2d_Curve)> curves;
1229 vector<Handle(Geom_Surface)> surfaces;
1234 TopTools_IndexedMapOfShape fmap;
1235 TopTools_IndexedMapOfShape emap;
1236 TopTools_IndexedMapOfShape pmap;
1239 FaceId2PythonSmp.clear();
1241 EdgeId2PythonSmp.clear();
1243 VertexId2PythonSmp.clear();
1245 assert(Py_IsInitialized());
1246 PyGILState_STATE gstate;
1247 gstate = PyGILState_Ensure();
1249 string theSizeMapStr;
1251 /****************************************************************************************
1253 *****************************************************************************************/
1255 string bad_end = "return";
1257 TopTools_IndexedMapOfShape _map;
1258 TopExp::MapShapes(aShape,TopAbs_VERTEX,_map);
1259 int ienf = _map.Extent();
1261 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next()) {
1262 TopoDS_Face f=TopoDS::Face(face_iter.Current());
1264 // make INTERNAL face oriented FORWARD (issue 0020993)
1265 if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
1266 f.Orientation(TopAbs_FORWARD);
1268 if (fmap.FindIndex(f) > 0)
1273 surfaces.push_back(BRep_Tool::Surface(f));
1275 /* create an object representing the face for blsurf */
1276 /* where face_id is an integer identifying the face.
1277 * surf_function is the function that defines the surface
1278 * (For this face, it will be called by blsurf with your_face_object_ptr
1279 * as last parameter.
1281 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
1283 /* by default a face has no tag (color). The following call sets it to the same value as the face_id : */
1284 cad_face_set_tag(fce, iface);
1286 /* Set face orientation (optional if you want a well oriented output mesh)*/
1287 if(f.Orientation() != TopAbs_FORWARD){
1288 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
1290 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
1293 if (HasSizeMapOnFace && !use_precad){
1294 // MESSAGE("A size map is defined on a face")
1295 // std::cout << "A size map is defined on a face" << std::endl;
1297 faceKey = FacesWithSizeMap.FindIndex(f);
1300 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end()){
1301 MESSAGE("A size map is defined on face :"<<faceKey)
1302 theSizeMapStr = FaceId2SizeMap[faceKey];
1303 // check if function ends with "return"
1304 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1306 // Expr To Python function, verification is performed at validation in GUI
1307 PyObject * obj = NULL;
1308 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1310 PyObject * func = NULL;
1311 func = PyObject_GetAttrString(main_mod, "f");
1312 FaceId2PythonSmp[iface]=func;
1313 FaceId2SizeMap.erase(faceKey);
1316 // Specific size map = Attractor
1317 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
1319 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
1320 if (attractor_iter->first == faceKey) {
1321 MESSAGE("Face indice: " << iface);
1322 MESSAGE("Adding attractor");
1324 double xyzCoords[3] = {attractor_iter->second[2],
1325 attractor_iter->second[3],
1326 attractor_iter->second[4]};
1328 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1329 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1330 BRepClass_FaceClassifier scl(f,P,1e-7);
1331 // scl.Perform() is bugged. The function was rewritten
1332 // BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1333 // OCC 6.5.2: scl.Perform() is not bugged anymore
1334 scl.Perform(f, P, 1e-7);
1335 TopAbs_State result = scl.State();
1336 MESSAGE("Position of point on face: "<<result);
1337 if ( result == TopAbs_OUT )
1338 MESSAGE("Point is out of face: node is not created");
1339 if ( result == TopAbs_UNKNOWN )
1340 MESSAGE("Point position on face is unknown: node is not created");
1341 if ( result == TopAbs_ON )
1342 MESSAGE("Point is on border of face: node is not created");
1343 if ( result == TopAbs_IN )
1345 // Point is inside face and not on border
1346 MESSAGE("Point is in face: node is created");
1347 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
1349 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1350 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1351 cad_point_set_tag(point_p, ienf);
1353 FaceId2AttractorCoords.erase(faceKey);
1358 std::map<int,BLSURFPlugin_Attractor* >::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
1359 if (clAttractor_iter != FaceId2ClassAttractor.end()){
1360 MESSAGE("Face indice: " << iface);
1361 MESSAGE("Adding attractor");
1362 FaceIndex2ClassAttractor[iface]=clAttractor_iter->second;
1363 FaceId2ClassAttractor.erase(clAttractor_iter);
1367 // Enforced Vertices
1368 faceKey = FacesWithEnforcedVertices.FindIndex(f);
1369 std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
1370 if (evmIt != FaceId2EnforcedVertexCoords.end()) {
1371 MESSAGE("Some enforced vertices are defined");
1372 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl;
1373 MESSAGE("Face indice: " << iface);
1374 MESSAGE("Adding enforced vertices");
1375 evl = evmIt->second;
1376 MESSAGE("Number of vertices to add: "<< evl.size());
1377 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
1378 for (; evlIt != evl.end(); ++evlIt) {
1379 BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
1380 xyzCoords.push_back(evlIt->at(2));
1381 xyzCoords.push_back(evlIt->at(3));
1382 xyzCoords.push_back(evlIt->at(4));
1383 MESSAGE("Check position of vertex =(" << xyzCoords[0] << "," << xyzCoords[1] << "," << xyzCoords[2] << ")");
1384 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
1385 BRepClass_FaceClassifier scl(f,P,1e-7);
1386 // OCC 6.3sp6 : scl.Perform() is bugged. The function was rewritten
1387 // BRepClass_FaceClassifierPerform(&scl,f,P,1e-7);
1388 // OCC 6.5.2: scl.Perform() is not bugged anymore
1389 scl.Perform(f, P, 1e-7);
1390 TopAbs_State result = scl.State();
1391 MESSAGE("Position of point on face: "<<result);
1392 if ( result == TopAbs_OUT ) {
1393 MESSAGE("Point is out of face: node is not created");
1394 if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1395 EnfVertexCoords2ProjVertex.erase(xyzCoords);
1396 EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1399 if ( result == TopAbs_UNKNOWN ) {
1400 MESSAGE("Point position on face is unknown: node is not created");
1401 if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1402 EnfVertexCoords2ProjVertex.erase(xyzCoords);
1403 EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1406 if ( result == TopAbs_ON ) {
1407 MESSAGE("Point is on border of face: node is not created");
1408 if (EnfVertexCoords2ProjVertex.find(xyzCoords) != EnfVertexCoords2ProjVertex.end()) {
1409 EnfVertexCoords2ProjVertex.erase(xyzCoords);
1410 EnfVertexCoords2EnfVertexList.erase(xyzCoords);
1413 if ( result == TopAbs_IN )
1415 // Point is inside face and not on border
1416 MESSAGE("Point is in face: node is created");
1417 double uvCoords[2] = {evlIt->at(0),evlIt->at(1)};
1419 MESSAGE("Add cad point on (u,v)=(" << uvCoords[0] << "," << uvCoords[1] << ") with id = " << ienf);
1420 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
1422 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(xyzCoords);
1423 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end() &&
1424 !enfCoordsIt->second.empty() )
1426 TopoDS_Vertex v = (*enfCoordsIt->second.begin())->vertex;
1427 if ( v.IsNull() ) v = (*enfCoordsIt->second.rbegin())->vertex;
1428 if ( !v.IsNull() ) {
1429 tag = pmap.Add( v );
1430 mergeSubmeshes.insert( aMesh.GetSubMesh( v ));
1431 //if ( tag != pmap.Extent() )
1435 if ( tag == 0 ) tag = ienf;
1436 cad_point_set_tag(point_p, tag);
1439 FaceId2EnforcedVertexCoords.erase(faceKey);
1442 // std::cout << "No enforced vertex defined" << std::endl;
1446 /****************************************************************************************
1448 now create the edges associated to this face
1449 *****************************************************************************************/
1451 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next()) {
1452 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
1453 int ic = emap.FindIndex(e);
1458 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
1460 if (HasSizeMapOnEdge){
1461 edgeKey = EdgesWithSizeMap.FindIndex(e);
1462 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end()) {
1463 MESSAGE("Sizemap defined on edge with index " << edgeKey);
1464 theSizeMapStr = EdgeId2SizeMap[edgeKey];
1465 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1467 // Expr To Python function, verification is performed at validation in GUI
1468 PyObject * obj = NULL;
1469 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1471 PyObject * func = NULL;
1472 func = PyObject_GetAttrString(main_mod, "f");
1473 EdgeId2PythonSmp[ic]=func;
1474 EdgeId2SizeMap.erase(edgeKey);
1478 /* attach the edge to the current blsurf face */
1479 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
1481 /* by default an edge has no tag (color). The following call sets it to the same value as the edge_id : */
1482 cad_edge_set_tag(edg, ic);
1484 /* by default, an edge does not necessalry appear in the resulting mesh,
1485 unless the following property is set :
1487 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
1489 /* by default an edge is a boundary edge */
1490 if (e.Orientation() == TopAbs_INTERNAL)
1491 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
1493 // pass existing nodes of sub-meshes to BLSURF
1494 SMESH_subMesh* sm = aMesh.GetSubMesh( e );
1495 if ( !sm->IsEmpty() )
1497 edgeSubmeshes.insert( sm );
1499 StdMeshers_FaceSide edgeOfFace( f, e, &aMesh, e.Orientation() == TopAbs_FORWARD,
1500 /*ignoreMedium=*/haveQudraticSubMesh);
1501 if ( edgeOfFace.MissVertexNode() )
1502 return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
1504 const int nbNodes = edgeOfFace.NbPoints();
1506 dcad_edge_discretization_t *dedge;
1507 dcad_get_edge_discretization(dcad, edg, &dedge);
1508 dcad_edge_discretization_set_vertex_count( dedge, nbNodes );
1510 const std::vector<UVPtStruct>& nodeData = edgeOfFace.GetUVPtStruct();
1511 for ( int iN = 0; iN < nbNodes; ++iN )
1513 const UVPtStruct& nData = nodeData[ iN ];
1514 double t = nData.param;
1515 real uv[2] = { nData.u, nData.v };
1516 SMESH_TNodeXYZ nXYZ( nData.node );
1517 dcad_edge_discretization_set_vertex_coordinates( dedge, iN+1, t, uv, nXYZ._xyz );
1519 dcad_edge_discretization_set_property(dedge, DISTENE_DCAD_PROPERTY_REQUIRED);
1522 /****************************************************************************************
1524 *****************************************************************************************/
1528 gp_Pnt2d e0 = curves.back()->Value(tmin);
1529 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
1530 Standard_Real d1=0,d2=0;
1533 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
1534 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
1538 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1541 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
1543 *ip = pmap.FindIndex(v);
1547 if (HasSizeMapOnVertex){
1548 vertexKey = VerticesWithSizeMap.FindIndex(v);
1549 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
1550 theSizeMapStr = VertexId2SizeMap[vertexKey];
1551 //MESSAGE("VertexId2SizeMap[faceKey]: " << VertexId2SizeMap[vertexKey]);
1552 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
1554 // Expr To Python function, verification is performed at validation in GUI
1555 PyObject * obj = NULL;
1556 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
1558 PyObject * func = NULL;
1559 func = PyObject_GetAttrString(main_mod, "f");
1560 VertexId2PythonSmp[*ip]=func;
1561 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
1566 // should not happen
1567 MESSAGE("An edge does not have 2 extremities.");
1570 // This defines the curves extremity connectivity
1571 cad_edge_set_extremities(edg, ip1, ip2);
1572 /* set the tag (color) to the same value as the extremity id : */
1573 cad_edge_set_extremities_tag(edg, ip1, ip2);
1576 cad_edge_set_extremities(edg, ip2, ip1);
1577 cad_edge_set_extremities_tag(edg, ip2, ip1);
1583 // Clear mesh from already meshed edges if possible else
1584 // remember that merge is needed
1585 set< SMESH_subMesh* >::iterator smIt = edgeSubmeshes.begin();
1586 for ( ; smIt != edgeSubmeshes.end(); ++smIt )
1588 SMESH_subMesh* sm = *smIt;
1589 SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
1590 /*complexFirst=*/false);
1591 while ( subsmIt->more() )
1593 sm = subsmIt->next();
1594 if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() )
1596 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1599 const SMDS_MeshNode* n = nIt->next();
1600 if ( n->NbInverseElements( SMDSAbs_Face ) > 0 )
1603 // add existing medium nodes to helper
1604 if ( aMesh.NbEdges( ORDER_QUADRATIC ) > 0 )
1606 SMDS_ElemIteratorPtr edgeIt = smDS->GetElements();
1607 while ( edgeIt->more() )
1608 helper.AddTLinks( static_cast<const SMDS_MeshEdge*>(edgeIt->next()));
1613 sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
1620 PyGILState_Release(gstate);
1623 /* Now launch the PreCAD process */
1624 status = precad_process(pcs);
1625 if(status != STATUS_OK){
1626 cout << "PreCAD processing failed with error code " << status << "\n";
1627 // Now we can delete the PreCAD session
1628 precad_session_delete(pcs);
1631 // retrieve the pre-processed CAD object
1632 cleanc = precad_new_cad(pcs);
1634 cout << "Unable to retrieve PreCAD result \n";
1636 cout << "PreCAD processing successfull \n";
1638 // #if BLSURF_VERSION_LONG >= "3.1.1"
1639 // /* We can now get the updated sizemaps (if any) */
1640 // // if(geo_sizemap_e)
1641 // // clean_geo_sizemap_e = precad_new_sizemap(pcs, geo_sizemap_e);
1643 // // if(geo_sizemap_f)
1644 // // clean_geo_sizemap_f = precad_new_sizemap(pcs, geo_sizemap_f);
1646 // if(iso_sizemap_p)
1647 // clean_iso_sizemap_p = precad_new_sizemap(pcs, iso_sizemap_p);
1649 // if(iso_sizemap_e)
1650 // clean_iso_sizemap_e = precad_new_sizemap(pcs, iso_sizemap_e);
1652 // if(iso_sizemap_f)
1653 // clean_iso_sizemap_f = precad_new_sizemap(pcs, iso_sizemap_f);
1656 // Now we can delete the PreCAD session
1657 precad_session_delete(pcs);
1661 blsurf_data_set_dcad(bls, dcad);
1663 // Give the pre-processed CAD object to the current BLSurf session
1664 blsurf_data_set_cad(bls, cleanc);
1667 // Use the original one
1668 blsurf_data_set_cad(bls, c);
1671 // #if BLSURF_VERSION_LONG >= "3.1.1"
1672 // blsurf_data_set_sizemap(bls, clean_iso_sizemap_p);
1673 // blsurf_data_set_sizemap(bls, clean_iso_sizemap_e);
1674 // blsurf_data_set_sizemap(bls, clean_iso_sizemap_f);
1677 std::cout << std::endl;
1678 std::cout << "Beginning of Surface Mesh generation" << std::endl;
1679 std::cout << std::endl;
1681 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1683 feclearexcept( FE_ALL_EXCEPT );
1684 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1690 status = blsurf_compute_mesh(bls);
1693 catch ( std::exception& exc ) {
1694 _comment += exc.what();
1696 catch (Standard_Failure& ex) {
1697 _comment += ex.DynamicType()->Name();
1698 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
1700 _comment += ex.GetMessageString();
1704 if ( _comment.empty() )
1705 _comment = "Exception in blsurf_compute_mesh()";
1707 if ( status != STATUS_OK) {
1708 // There was an error while meshing
1709 //return error(_comment);
1712 std::cout << std::endl;
1713 std::cout << "End of Surface Mesh generation" << std::endl;
1714 std::cout << std::endl;
1717 blsurf_data_get_mesh(bls, &msh);
1719 /* release the mesh object */
1720 blsurf_data_regain_mesh(bls, msh);
1721 return error(_comment);
1724 std::string GMFFileName = BLSURFPlugin_Hypothesis::GetDefaultGMFFile();
1726 GMFFileName = _hypothesis->GetGMFFile();
1727 if (GMFFileName != "") {
1728 // bool GMFFileMode = _hypothesis->GetGMFFileMode();
1729 bool asciiFound = (GMFFileName.find(".mesh",GMFFileName.length()-5) != std::string::npos);
1730 bool binaryFound = (GMFFileName.find(".meshb",GMFFileName.length()-6) != std::string::npos);
1731 if (!asciiFound && !binaryFound)
1732 GMFFileName.append(".mesh");
1737 GMFFileName.append("b");
1739 GMFFileName.append(".meshb");
1744 GMFFileName.append(".mesh");
1747 mesh_write_mesh(msh, GMFFileName.c_str());
1750 /* retrieve mesh data (see distene/mesh.h) */
1751 integer nv, ne, nt, nq, vtx[4], tag;
1754 mesh_get_vertex_count(msh, &nv);
1755 mesh_get_edge_count(msh, &ne);
1756 mesh_get_triangle_count(msh, &nt);
1757 mesh_get_quadrangle_count(msh, &nq);
1760 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1761 SMDS_MeshNode** nodes = new SMDS_MeshNode*[nv+1];
1762 bool* tags = new bool[nv+1];
1764 /* enumerated vertices */
1765 for(int iv=1;iv<=nv;iv++) {
1766 mesh_get_vertex_coordinates(msh, iv, xyz);
1767 mesh_get_vertex_tag(msh, iv, &tag);
1768 // Issue 0020656. Use vertex coordinates
1769 if ( tag > 0 && tag <= pmap.Extent() ) {
1770 TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
1771 double tol = BRep_Tool::Tolerance( v );
1772 gp_Pnt p = BRep_Tool::Pnt( v );
1773 if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
1774 xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
1776 tag = 0; // enforced or attracted vertex
1778 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
1780 // Create group of enforced vertices if requested
1781 BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
1783 projVertex.push_back((double)xyz[0]);
1784 projVertex.push_back((double)xyz[1]);
1785 projVertex.push_back((double)xyz[2]);
1786 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
1787 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end()) {
1788 MESSAGE("Found enforced vertex @ " << xyz[0] << ", " << xyz[1] << ", " << xyz[2]);
1789 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
1790 BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
1791 for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
1792 currentEnfVertex = (*enfListIt);
1793 if (currentEnfVertex->grpName != "") {
1794 bool groupDone = false;
1795 SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
1796 MESSAGE("currentEnfVertex->grpName: " << currentEnfVertex->grpName);
1797 MESSAGE("Parsing the groups of the mesh");
1798 while (grIt->more()) {
1799 SMESH_Group * group = grIt->next();
1800 if ( !group ) continue;
1801 MESSAGE("Group: " << group->GetName());
1802 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
1803 if ( !groupDS ) continue;
1804 MESSAGE("group->SMDSGroup().GetType(): " << (groupDS->GetType()));
1805 MESSAGE("group->SMDSGroup().GetType()==SMDSAbs_Node: " << (groupDS->GetType()==SMDSAbs_Node));
1806 MESSAGE("currentEnfVertex->grpName.compare(group->GetStoreName())==0: " << (currentEnfVertex->grpName.compare(group->GetName())==0));
1807 if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
1808 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
1809 aGroupDS->SMDSGroup().Add(nodes[iv]);
1810 MESSAGE("Node ID: " << nodes[iv]->GetID());
1811 // How can I inform the hypothesis ?
1812 // _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
1814 MESSAGE("Successfully added enforced vertex to existing group " << currentEnfVertex->grpName);
1821 SMESH_Group* aGroup = aMesh.AddGroup(SMDSAbs_Node, currentEnfVertex->grpName.c_str(), groupId);
1822 aGroup->SetName( currentEnfVertex->grpName.c_str() );
1823 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
1824 aGroupDS->SMDSGroup().Add(nodes[iv]);
1825 MESSAGE("Successfully created enforced vertex group " << currentEnfVertex->grpName);
1829 throw SALOME_Exception(LOCALIZED("An enforced vertex node was not added to a group"));
1832 MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
1837 // internal point are tagged to zero
1838 if(tag > 0 && tag <= pmap.Extent() ){
1839 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
1846 /* enumerate edges */
1847 for(int it=1;it<=ne;it++) {
1849 mesh_get_edge_vertices(msh, it, vtx);
1850 mesh_get_edge_tag(msh, it, &tag);
1852 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
1853 tags[vtx[0]] = false;
1856 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
1857 tags[vtx[1]] = false;
1859 SMDS_MeshEdge* edg = helper.AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
1860 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
1863 /* enumerate triangles */
1864 for(int it=1;it<=nt;it++) {
1865 mesh_get_triangle_vertices(msh, it, vtx);
1866 SMDS_MeshFace* tri = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
1867 mesh_get_triangle_tag(msh, it, &tag);
1868 meshDS->SetMeshElementOnShape(tri, TopoDS::Face(fmap(tag)));
1870 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1871 tags[vtx[0]] = false;
1874 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1875 tags[vtx[1]] = false;
1878 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1879 tags[vtx[2]] = false;
1883 /* enumerate quadrangles */
1884 for(int it=1;it<=nq;it++) {
1885 mesh_get_quadrangle_vertices(msh, it, vtx);
1886 SMDS_MeshFace* quad = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
1887 mesh_get_quadrangle_tag(msh, it, &tag);
1888 meshDS->SetMeshElementOnShape(quad, TopoDS::Face(fmap(tag)));
1890 meshDS->SetNodeOnFace(nodes[vtx[0]], TopoDS::Face(fmap(tag)));
1891 tags[vtx[0]] = false;
1894 meshDS->SetNodeOnFace(nodes[vtx[1]], TopoDS::Face(fmap(tag)));
1895 tags[vtx[1]] = false;
1898 meshDS->SetNodeOnFace(nodes[vtx[2]], TopoDS::Face(fmap(tag)));
1899 tags[vtx[2]] = false;
1902 meshDS->SetNodeOnFace(nodes[vtx[3]], TopoDS::Face(fmap(tag)));
1903 tags[vtx[3]] = false;
1907 // SetIsAlwaysComputed( true ) to sub-meshes of EDGEs w/o mesh
1908 TopLoc_Location loc; double f,l;
1909 for (int i = 1; i <= emap.Extent(); i++)
1910 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
1911 sm->SetIsAlwaysComputed( true );
1912 for (int i = 1; i <= pmap.Extent(); i++)
1913 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( pmap( i )))
1914 if ( !sm->IsMeshComputed() )
1915 sm->SetIsAlwaysComputed( true );
1919 set< SMESH_subMesh* >::iterator smIt = mergeSubmeshes.begin();
1920 for ( ; smIt != mergeSubmeshes.end(); ++smIt )
1922 SMESH_subMesh* sm = *smIt;
1923 SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
1924 /*complexFirst=*/false);
1925 TIDSortedNodeSet nodesOnEdge;
1926 double mergeTolerance = 1e-7, tol;
1927 while ( subsmIt->more() )
1929 // get nodes from an edge
1930 sm = subsmIt->next();
1931 if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() )
1933 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1934 while ( nIt->more() )
1935 nodesOnEdge.insert( nIt->next() );
1937 // get max tolerance
1938 TopoDS_Shape subShape = sm->GetSubShape();
1939 if ( subShape.ShapeType() == TopAbs_EDGE )
1940 tol = BRep_Tool::Tolerance( TopoDS::Edge( subShape ));
1942 tol = BRep_Tool::Tolerance( TopoDS::Vertex( subShape ));
1943 mergeTolerance = Max( tol, mergeTolerance );
1945 // find nodes to merge
1946 SMESH_MeshEditor::TListOfListOfNodes nodeGroupsToMerge;
1947 SMESH_MeshEditor editor( &aMesh );
1948 editor.FindCoincidentNodes( nodesOnEdge, mergeTolerance*2, nodeGroupsToMerge );
1950 editor.MergeNodes( nodeGroupsToMerge );
1956 /* release the mesh object */
1957 blsurf_data_regain_mesh(bls, msh);
1959 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1961 if ( oldFEFlags > 0 )
1962 feenableexcept( oldFEFlags );
1963 feclearexcept( FE_ALL_EXCEPT );
1967 std::cout << "FacesWithSizeMap" << std::endl;
1968 FacesWithSizeMap.Statistics(std::cout);
1969 std::cout << "EdgesWithSizeMap" << std::endl;
1970 EdgesWithSizeMap.Statistics(std::cout);
1971 std::cout << "VerticesWithSizeMap" << std::endl;
1972 VerticesWithSizeMap.Statistics(std::cout);
1973 std::cout << "FacesWithEnforcedVertices" << std::endl;
1974 FacesWithEnforcedVertices.Statistics(std::cout);
1977 MESSAGE("END OF BLSURFPlugin_BLSURF::Compute()");
1981 #ifdef WITH_SMESH_CANCEL_COMPUTE
1982 void BLSURFPlugin_BLSURF::CancelCompute()
1984 _compute_canceled = true;
1988 //=============================================================================
1992 //=============================================================================
1994 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS, SMDS_MeshNode* node, const TopoDS_Shape& ed) {
1995 const TopoDS_Edge edge = TopoDS::Edge(ed);
1997 gp_Pnt pnt(node->X(), node->Y(), node->Z());
1999 Standard_Real p0 = 0.0;
2000 Standard_Real p1 = 1.0;
2001 TopLoc_Location loc;
2002 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
2004 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
2005 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
2008 if ( proj.NbPoints() > 0 )
2010 pa = (double)proj.LowerDistanceParameter();
2011 // Issue 0020656. Move node if it is too far from edge
2012 gp_Pnt curve_pnt = curve->Value( pa );
2013 double dist2 = pnt.SquareDistance( curve_pnt );
2014 double tol = BRep_Tool::Tolerance( edge );
2015 if ( 1e-14 < dist2 && dist2 <= 1000*tol ) // large enough and within tolerance
2017 curve_pnt.Transform( loc );
2018 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
2021 // GProp_GProps LProps;
2022 // BRepGProp::LinearProperties(ed, LProps);
2023 // double lg = (double)LProps.Mass();
2025 meshDS->SetNodeOnEdge(node, edge, pa);
2028 //=============================================================================
2032 //=============================================================================
2034 ostream & BLSURFPlugin_BLSURF::SaveTo(ostream & save)
2039 //=============================================================================
2043 //=============================================================================
2045 istream & BLSURFPlugin_BLSURF::LoadFrom(istream & load)
2050 //=============================================================================
2054 //=============================================================================
2056 ostream & operator << (ostream & save, BLSURFPlugin_BLSURF & hyp)
2058 return hyp.SaveTo( save );
2061 //=============================================================================
2065 //=============================================================================
2067 istream & operator >> (istream & load, BLSURFPlugin_BLSURF & hyp)
2069 return hyp.LoadFrom( load );
2072 /* Curve definition function See cad_curv_t in file distene/cad.h for
2074 * NOTE : if when your CAD systems evaluates second
2075 * order derivatives it also computes first order derivatives and
2076 * function evaluation, you can optimize this example by making only
2077 * one CAD call and filling the necessary uv, dt, dtt arrays.
2079 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
2081 /* t is given. It contains the t (time) 1D parametric coordintaes
2082 of the point PreCAD/BLSurf is querying on the curve */
2084 /* user_data identifies the edge PreCAD/BLSurf is querying
2085 * (see cad_edge_new later in this example) */
2086 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
2089 /* BLSurf is querying the function evaluation */
2092 uv[0]=P.X(); uv[1]=P.Y();
2096 /* query for the first order derivatives */
2099 dt[0]=V1.X(); dt[1]=V1.Y();
2103 /* query for the second order derivatives */
2106 dtt[0]=V2.X(); dtt[1]=V2.Y();
2112 /* Surface definition function.
2113 * See cad_surf_t in file distene/cad.h for more information.
2114 * NOTE : if when your CAD systems evaluates second order derivatives it also
2115 * computes first order derivatives and function evaluation, you can optimize
2116 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
2119 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
2120 real *duu, real *duv, real *dvv, void *user_data)
2122 /* uv[2] is given. It contains the u,v coordinates of the point
2123 * PreCAD/BLSurf is querying on the surface */
2125 /* user_data identifies the face PreCAD/BLSurf is querying (see
2126 * cad_face_new later in this example)*/
2127 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
2131 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
2132 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
2139 geometry->D1(uv[0],uv[1],P,D1U,D1V);
2140 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
2141 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
2144 if(duu && duv && dvv){
2148 gp_Vec D2U,D2V,D2UV;
2150 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
2151 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
2152 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
2153 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
2160 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
2163 if (my_u_min > uv[0]) {
2166 if (my_v_min > uv[1]) {
2169 if (my_u_max < uv[0]) {
2172 if (my_v_max < uv[1]) {
2176 //MESSAGE("size_on_surface")
2177 if (FaceId2PythonSmp.count(face_id) != 0){
2178 //MESSAGE("A size map is used to calculate size on face : "<<face_id)
2179 PyObject * pyresult = NULL;
2180 PyObject* new_stderr = NULL;
2181 assert(Py_IsInitialized());
2182 PyGILState_STATE gstate;
2183 gstate = PyGILState_Ensure();
2184 pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
2186 if ( pyresult == NULL){
2188 string err_description="";
2189 new_stderr = newPyStdOut(err_description);
2190 PySys_SetObject((char*)"stderr", new_stderr);
2192 PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
2193 Py_DECREF(new_stderr);
2194 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
2195 result = *((double*)user_data);
2198 result = PyFloat_AsDouble(pyresult);
2199 Py_DECREF(pyresult);
2201 // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
2203 PyGILState_Release(gstate);
2205 else if (FaceIndex2ClassAttractor.count(face_id) !=0 && !FaceIndex2ClassAttractor[face_id]->Empty()){
2206 // MESSAGE("attractor used on face :"<<face_id)
2207 // MESSAGE("List of attractor is not empty")
2208 // MESSAGE("Attractor empty : "<< FaceIndex2ClassAttractor[face_id]->Empty())
2209 double result = FaceIndex2ClassAttractor[face_id]->GetSize(uv[0],uv[1]);
2211 // MESSAGE("f(" << uv[0] << "," << uv[1] << ")" << " = " << result);
2214 // MESSAGE("List of attractor is empty !!!")
2215 *size = *((double*)user_data);
2220 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
2222 if (EdgeId2PythonSmp.count(edge_id) != 0){
2223 PyObject * pyresult = NULL;
2224 PyObject* new_stderr = NULL;
2225 assert(Py_IsInitialized());
2226 PyGILState_STATE gstate;
2227 gstate = PyGILState_Ensure();
2228 pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
2230 if ( pyresult == NULL){
2232 string err_description="";
2233 new_stderr = newPyStdOut(err_description);
2234 PySys_SetObject((char*)"stderr", new_stderr);
2236 PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
2237 Py_DECREF(new_stderr);
2238 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
2239 result = *((double*)user_data);
2242 result = PyFloat_AsDouble(pyresult);
2243 Py_DECREF(pyresult);
2246 PyGILState_Release(gstate);
2249 *size = *((double*)user_data);
2254 status_t size_on_vertex(integer point_id, real *size, void *user_data)
2256 if (VertexId2PythonSmp.count(point_id) != 0){
2257 PyObject * pyresult = NULL;
2258 PyObject* new_stderr = NULL;
2259 assert(Py_IsInitialized());
2260 PyGILState_STATE gstate;
2261 gstate = PyGILState_Ensure();
2262 pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
2264 if ( pyresult == NULL){
2266 string err_description="";
2267 new_stderr = newPyStdOut(err_description);
2268 PySys_SetObject((char*)"stderr", new_stderr);
2270 PySys_SetObject((char*)"stderr", PySys_GetObject((char*)"__stderr__"));
2271 Py_DECREF(new_stderr);
2272 MESSAGE("Can't evaluate f()" << " error is " << err_description);
2273 result = *((double*)user_data);
2276 result = PyFloat_AsDouble(pyresult);
2277 Py_DECREF(pyresult);
2280 PyGILState_Release(gstate);
2283 *size = *((double*)user_data);
2289 * The following function will be called for PreCAD/BLSurf message
2290 * printing. See context_set_message_callback (later in this
2291 * template) for how to set user_data.
2293 status_t message_cb(message_t *msg, void *user_data)
2295 integer errnumber = 0;
2297 message_get_number(msg, &errnumber);
2298 message_get_description(msg, &desc);
2300 if ( errnumber < 0 || err.find("license") != string::npos ) {
2301 string * error = (string*)user_data;
2302 // if ( !error->empty() )
2304 // remove ^A from the tail
2305 int len = strlen( desc );
2306 while (len > 0 && desc[len-1] != '\n')
2308 error->append( desc, len );
2311 std::cout << desc << std::endl;
2316 /* This is the interrupt callback. PreCAD/BLSurf will call this
2317 * function regularily. See the file distene/interrupt.h
2319 status_t interrupt_cb(integer *interrupt_status, void *user_data)
2321 integer you_want_to_continue = 1;
2322 #ifdef WITH_SMESH_CANCEL_COMPUTE
2323 BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
2324 you_want_to_continue = !tmp->computeCanceled();
2327 if(you_want_to_continue)
2329 *interrupt_status = INTERRUPT_CONTINUE;
2332 else /* you want to stop BLSurf */
2334 *interrupt_status = INTERRUPT_STOP;
2335 return STATUS_ERROR;
2339 //=============================================================================
2343 //=============================================================================
2344 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
2345 const TopoDS_Shape& aShape,
2346 MapShapeNbElems& aResMap)
2348 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
2349 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize();
2350 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
2351 //double _angleMeshS = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshS();
2352 double _angleMeshC = BLSURFPlugin_Hypothesis::GetDefaultAngleMeshC();
2353 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
2355 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
2356 _phySize = _hypothesis->GetPhySize();
2357 //_geometricMesh = (int) hyp->GetGeometricMesh();
2358 //_angleMeshS = hyp->GetAngleMeshS();
2359 _angleMeshC = _hypothesis->GetAngleMeshC();
2360 _quadAllowed = _hypothesis->GetQuadAllowed();
2362 //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
2363 // GetDefaultPhySize() sometimes leads to computation failure
2364 _phySize = aMesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
2365 MESSAGE("BLSURFPlugin_BLSURF::SetParameters using defaults");
2368 bool IsQuadratic = false;
2373 TopTools_DataMapOfShapeInteger EdgesMap;
2374 double fullLen = 0.0;
2375 double fullNbSeg = 0;
2376 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2377 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
2378 if( EdgesMap.IsBound(E) )
2380 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
2381 double aLen = SMESH_Algo::EdgeLength(E);
2384 if(_physicalMesh==1) {
2385 nb1d = (int)( aLen/_phySize + 1 );
2390 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
2391 double fullAng = 0.0;
2392 double dp = (l-f)/200;
2397 for(int j=2; j<=200; j++) {
2400 fullAng += fabs(V1.Angle(V2));
2404 nb1d = (int)( fullAng/_angleMeshC + 1 );
2407 std::vector<int> aVec(SMDSEntity_Last);
2408 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2409 if( IsQuadratic > 0 ) {
2410 aVec[SMDSEntity_Node] = 2*nb1d - 1;
2411 aVec[SMDSEntity_Quad_Edge] = nb1d;
2414 aVec[SMDSEntity_Node] = nb1d - 1;
2415 aVec[SMDSEntity_Edge] = nb1d;
2417 aResMap.insert(std::make_pair(sm,aVec));
2418 EdgesMap.Bind(E,nb1d);
2420 double ELen = fullLen/fullNbSeg;
2424 // try to evaluate as in MEFISTO
2425 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2426 TopoDS_Face F = TopoDS::Face( exp.Current() );
2427 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2429 BRepGProp::SurfaceProperties(F,G);
2430 double anArea = G.Mass();
2432 std::vector<int> nb1dVec;
2433 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
2434 int nbSeg = EdgesMap.Find(exp1.Current());
2436 nb1dVec.push_back( nbSeg );
2439 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
2440 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
2443 if ( nb1dVec.size() == 4 ) // quadrangle geom face
2445 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
2447 nbNodes = (n1 + 1) * (n2 + 1);
2452 nbTria = nbQuad = nbTria / 3 + 1;
2455 std::vector<int> aVec(SMDSEntity_Last,0);
2457 int nb1d_in = (nbTria*3 - nb1d) / 2;
2458 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
2459 aVec[SMDSEntity_Quad_Triangle] = nbTria;
2460 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
2463 aVec[SMDSEntity_Node] = nbNodes;
2464 aVec[SMDSEntity_Triangle] = nbTria;
2465 aVec[SMDSEntity_Quadrangle] = nbQuad;
2467 aResMap.insert(std::make_pair(sm,aVec));
2474 BRepGProp::VolumeProperties(aShape,G);
2475 double aVolume = G.Mass();
2476 double tetrVol = 0.1179*ELen*ELen*ELen;
2477 int nbVols = int(aVolume/tetrVol);
2478 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
2479 std::vector<int> aVec(SMDSEntity_Last);
2480 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2482 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
2483 aVec[SMDSEntity_Quad_Tetra] = nbVols;
2486 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
2487 aVec[SMDSEntity_Tetra] = nbVols;
2489 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2490 aResMap.insert(std::make_pair(sm,aVec));
2495 //=============================================================================
2497 * Rewritting of the BRepClass_FaceClassifier::Perform function which is bugged (CAS 6.3sp6)
2498 * Following line was added:
2499 * myExtrem.Perform(P);
2501 //=============================================================================
2502 void BLSURFPlugin_BLSURF::BRepClass_FaceClassifierPerform(BRepClass_FaceClassifier* fc,
2503 const TopoDS_Face& face,
2505 const Standard_Real Tol)
2507 //-- Voir BRepExtrema_ExtPF.cxx
2508 BRepAdaptor_Surface Surf(face);
2509 Standard_Real U1, U2, V1, V2;
2510 BRepTools::UVBounds(face, U1, U2, V1, V2);
2511 Extrema_ExtPS myExtrem;
2512 myExtrem.Initialize(Surf, U1, U2, V1, V2, Tol, Tol);
2513 myExtrem.Perform(P);
2514 //----------------------------------------------------------
2515 //-- On cherche le point le plus proche , PUIS
2516 //-- On le classifie.
2517 Standard_Integer nbv = 0; // xpu
2518 Standard_Real MaxDist = RealLast();
2519 Standard_Integer indice = 0;
2520 if (myExtrem.IsDone()) {
2521 nbv = myExtrem.NbExt();
2522 for (Standard_Integer i = 1; i <= nbv; i++) {
2523 #if OCC_VERSION_LARGE > 0x06040000 // Porting to OCCT6.5.1
2524 Standard_Real d = myExtrem.SquareDistance(i);
2526 Standard_Real d = myExtrem.Value(i);
2537 Standard_Real U1,U2;
2538 myExtrem.Point(indice).Parameter(U1, U2);
2539 Puv.SetCoord(U1, U2);
2540 fc->Perform(face, Puv, Tol);
2543 fc->Perform(face, gp_Pnt2d(U1-1.0,V1 - 1.0), Tol); //-- NYI etc BUG PAS BEAU En attendant l acces a rejected
2544 //-- le resultat est TopAbs_OUT;