1 // Copyright (C) 2007-2016 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, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 // File : BLSURFPlugin_BLSURF.cxx
22 // Authors : Francis KLOSS (OCC) & Patrick LAUG (INRIA) & Lioka RAZAFINDRAZAKA (CEA)
23 // & Aurelien ALLEAUME (DISTENE)
24 // Size maps developement: Nicolas GEIMER (OCC) & Gilles DAVID (EURIWARE)
27 #include "BLSURFPlugin_BLSURF.hxx"
28 #include "BLSURFPlugin_Hypothesis.hxx"
29 #include "BLSURFPlugin_Attractor.hxx"
32 #include <meshgems/meshgems.h>
33 #include <meshgems/cadsurf.h>
36 #include <structmember.h>
39 #include <Basics_Utils.hxx>
40 #include <Basics_OCCTVersion.hxx>
42 #include <SMDS_EdgePosition.hxx>
43 #include <SMESHDS_Group.hxx>
44 #include <SMESH_Gen.hxx>
45 #include <SMESH_Group.hxx>
46 #include <SMESH_Mesh.hxx>
47 #include <SMESH_MeshEditor.hxx>
48 #include <SMESH_MesherHelper.hxx>
49 #include <StdMeshers_FaceSide.hxx>
50 #include <StdMeshers_ViscousLayers2D.hxx>
51 #include <SMESH_File.hxx>
53 #include <utilities.h>
61 // OPENCASCADE includes
62 #include <BRepBuilderAPI_MakeFace.hxx>
63 #include <BRepBuilderAPI_MakePolygon.hxx>
64 #include <BRepBuilderAPI_MakeWire.hxx>
65 #include <BRepGProp.hxx>
66 #include <BRepTools.hxx>
67 #include <BRep_Builder.hxx>
68 #include <BRep_Tool.hxx>
69 #include <GProp_GProps.hxx>
70 #include <Geom2d_Curve.hxx>
71 #include <Geom2d_Line.hxx>
72 #include <GeomAPI_ProjectPointOnCurve.hxx>
73 #include <GeomAPI_ProjectPointOnSurf.hxx>
74 #include <Geom_Curve.hxx>
75 #include <Geom_Surface.hxx>
76 #include <NCollection_DataMap.hxx>
77 #include <NCollection_Map.hxx>
78 #include <Standard_ErrorHandler.hxx>
80 #include <TopExp_Explorer.hxx>
81 #include <TopTools_DataMapOfShapeInteger.hxx>
82 #include <TopTools_IndexedMapOfShape.hxx>
83 #include <TopTools_MapOfShape.hxx>
85 #include <TopoDS_Compound.hxx>
86 #include <TopoDS_Edge.hxx>
87 #include <TopoDS_Face.hxx>
88 #include <TopoDS_Shape.hxx>
89 #include <TopoDS_Vertex.hxx>
90 #include <TopoDS_Wire.hxx>
92 #include <gp_Pnt2d.hxx>
102 /* ==================================
103 * =========== PYTHON ==============
104 * ==================================*/
115 PyStdOut_dealloc(PyStdOut *self)
121 PyStdOut_write(PyStdOut *self, PyObject *args)
125 if (!PyArg_ParseTuple(args, "t#:write",&c, &l))
128 *(self->out)=*(self->out)+c;
134 static PyMethodDef PyStdOut_methods[] = {
135 {"write", (PyCFunction)PyStdOut_write, METH_VARARGS,
136 PyDoc_STR("write(string) -> None")},
137 {NULL, NULL} /* sentinel */
140 static PyMemberDef PyStdOut_memberlist[] = {
141 {(char*)"softspace", T_INT, offsetof(PyStdOut, softspace), 0,
142 (char*)"flag indicating that a space needs to be printed; used by print"},
143 {NULL} /* Sentinel */
146 static PyTypeObject PyStdOut_Type = {
147 /* The ob_type field must be initialized in the module init function
148 * to be portable to Windows without using C++. */
149 PyObject_HEAD_INIT(NULL)
152 sizeof(PyStdOut), /*tp_basicsize*/
155 (destructor)PyStdOut_dealloc, /*tp_dealloc*/
162 0, /*tp_as_sequence*/
167 PyObject_GenericGetAttr, /*tp_getattro*/
168 /* softspace is writable: we must supply tp_setattro */
169 PyObject_GenericSetAttr, /* tp_setattro */
171 Py_TPFLAGS_DEFAULT, /*tp_flags*/
175 0, /*tp_richcompare*/
176 0, /*tp_weaklistoffset*/
179 PyStdOut_methods, /*tp_methods*/
180 PyStdOut_memberlist, /*tp_members*/
194 PyObject * newPyStdOut( std::string& out )
196 PyStdOut* self = PyObject_New(PyStdOut, &PyStdOut_Type);
201 return (PyObject*)self;
206 ////////////////////////END PYTHON///////////////////////////
208 //////////////////MY MAPS////////////////////////////////////////
210 TopTools_IndexedMapOfShape FacesWithSizeMap;
211 std::map<int,string> FaceId2SizeMap;
212 TopTools_IndexedMapOfShape EdgesWithSizeMap;
213 std::map<int,string> EdgeId2SizeMap;
214 TopTools_IndexedMapOfShape VerticesWithSizeMap;
215 std::map<int,string> VertexId2SizeMap;
217 std::map<int,PyObject*> FaceId2PythonSmp;
218 std::map<int,PyObject*> EdgeId2PythonSmp;
219 std::map<int,PyObject*> VertexId2PythonSmp;
221 typedef std::map<int, std::vector< BLSURFPlugin_Attractor* > > TId2ClsAttractorVec;
222 TId2ClsAttractorVec FaceId2ClassAttractor;
223 TId2ClsAttractorVec FaceIndex2ClassAttractor;
224 std::map<int,std::vector<double> > FaceId2AttractorCoords;
227 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
228 std::map< int, BLSURFPlugin_Hypothesis::TEnfVertexCoordsList > FaceId2EnforcedVertexCoords;
229 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexCoords > EnfVertexCoords2ProjVertex;
230 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList > EnfVertexCoords2EnfVertexList;
231 SMESH_MesherHelper* theHelper;
233 bool HasSizeMapOnFace=false;
234 bool HasSizeMapOnEdge=false;
235 bool HasSizeMapOnVertex=false;
236 //bool HasAttractorOnFace=false;
238 //=============================================================================
242 //=============================================================================
244 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId,
248 : SMESH_2D_Algo(hypId, studyId, gen)
250 _name = theHasGEOM ? "MG-CADSurf" : "MG-CADSurf_NOGEOM";//"BLSURF";
251 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
252 _compatibleHypothesis.push_back(BLSURFPlugin_Hypothesis::GetHypType(theHasGEOM));
254 _compatibleHypothesis.push_back(StdMeshers_ViscousLayers2D::GetHypType());
255 _requireDiscreteBoundary = false;
256 _onlyUnaryInput = false;
258 _supportSubmeshes = true;
259 _requireShape = theHasGEOM;
261 smeshGen_i = SMESH_Gen_i::GetSMESHGen();
262 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
263 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
266 myStudy = aStudyMgr->GetStudyByID(_studyId);
268 /* Initialize the Python interpreter */
269 assert(Py_IsInitialized());
270 PyGILState_STATE gstate;
271 gstate = PyGILState_Ensure();
274 main_mod = PyImport_AddModule("__main__");
277 main_dict = PyModule_GetDict(main_mod);
279 PyRun_SimpleString("from math import *");
280 PyGILState_Release(gstate);
282 FacesWithSizeMap.Clear();
283 FaceId2SizeMap.clear();
284 EdgesWithSizeMap.Clear();
285 EdgeId2SizeMap.clear();
286 VerticesWithSizeMap.Clear();
287 VertexId2SizeMap.clear();
288 FaceId2PythonSmp.clear();
289 EdgeId2PythonSmp.clear();
290 VertexId2PythonSmp.clear();
291 FaceId2AttractorCoords.clear();
292 FaceId2ClassAttractor.clear();
293 FaceIndex2ClassAttractor.clear();
294 FacesWithEnforcedVertices.Clear();
295 FaceId2EnforcedVertexCoords.clear();
296 EnfVertexCoords2ProjVertex.clear();
297 EnfVertexCoords2EnfVertexList.clear();
299 _compute_canceled = false;
302 //=============================================================================
306 //=============================================================================
308 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
313 //=============================================================================
317 //=============================================================================
319 bool BLSURFPlugin_BLSURF::CheckHypothesis
321 const TopoDS_Shape& aShape,
322 SMESH_Hypothesis::Hypothesis_Status& aStatus)
325 _haveViscousLayers = false;
327 list<const SMESHDS_Hypothesis*>::const_iterator itl;
328 const SMESHDS_Hypothesis* theHyp;
330 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape,
331 /*ignoreAuxiliary=*/false);
332 aStatus = SMESH_Hypothesis::HYP_OK;
335 return true; // can work with no hypothesis
338 for ( itl = hyps.begin(); itl != hyps.end() && ( aStatus == HYP_OK ); ++itl )
341 string hypName = theHyp->GetName();
342 if ( hypName == BLSURFPlugin_Hypothesis::GetHypType(true) ||
343 hypName == BLSURFPlugin_Hypothesis::GetHypType(false) )
345 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
347 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
348 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
349 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
350 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
352 else if ( hypName == StdMeshers_ViscousLayers2D::GetHypType() )
354 if ( !_haveViscousLayers )
356 if ( error( StdMeshers_ViscousLayers2D::CheckHypothesis( aMesh, aShape, aStatus )))
357 _haveViscousLayers = true;
362 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
365 return aStatus == SMESH_Hypothesis::HYP_OK;
368 //=============================================================================
370 * Pass parameters to MG-CADSurf
372 //=============================================================================
374 inline std::string val_to_string(double d)
376 std::ostringstream o;
381 inline std::string val_to_string_rel(double d)
383 std::ostringstream o;
389 inline std::string val_to_string(int i)
391 std::ostringstream o;
396 inline std::string val_to_string_rel(int i)
398 std::ostringstream o;
404 double _smp_phy_size;
405 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
406 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
407 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
414 /////////////////////////////////////////////////////////
416 projectionPoint getProjectionPoint(TopoDS_Face& theFace, const gp_Pnt& thePoint)
418 projectionPoint myPoint;
420 if ( theFace.IsNull() )
422 TopoDS_Shape foundFace, myShape = theHelper->GetSubShape();
423 TopTools_MapOfShape checkedFaces;
424 std::map< double, std::pair< TopoDS_Face, gp_Pnt2d > > dist2face;
426 for ( TopExp_Explorer exp ( myShape, TopAbs_FACE ); exp.More(); exp.Next())
428 const TopoDS_Face& face = TopoDS::Face( exp.Current() );
429 if ( !checkedFaces.Add( face )) continue;
431 // check distance to face
432 Handle(ShapeAnalysis_Surface) surface = theHelper->GetSurface( face );
433 gp_Pnt2d uv = surface->ValueOfUV( thePoint, Precision::Confusion());
434 double distance = surface->Gap();
435 if ( distance > Precision::Confusion() )
437 // the face is far, store for future analysis
438 dist2face.insert( std::make_pair( distance, std::make_pair( face, uv )));
442 // check location on the face
443 BRepClass_FaceClassifier FC( face, uv, BRep_Tool::Tolerance( face ));
444 if ( FC.State() == TopAbs_IN )
446 if ( !foundFace.IsNull() )
447 return myPoint; // thePoint seems to be TopAbs_ON
449 myPoint.uv = uv.XY();
450 myPoint.xyz = surface->Value( uv ).XYZ();
453 if ( FC.State() == TopAbs_ON )
457 if ( foundFace.IsNull() )
459 // find the closest face
460 std::map< double, std::pair< TopoDS_Face, gp_Pnt2d > >::iterator d2f = dist2face.begin();
461 for ( ; d2f != dist2face.end(); ++d2f )
463 const TopoDS_Face& face = d2f->second.first;
464 const gp_Pnt2d & uv = d2f->second.second;
465 BRepClass_FaceClassifier FC( face, uv, Precision::Confusion());
466 if ( FC.State() == TopAbs_IN )
469 myPoint.uv = uv.XY();
470 myPoint.xyz = theHelper->GetSurface( face )->Value( uv ).XYZ();
475 // set the resultShape
476 // if ( foundFace.IsNull() )
477 // throw SMESH_ComputeError(COMPERR_BAD_PARMETERS,
478 // "getProjectionPoint: can't find a face by a vertex");
479 theFace = TopoDS::Face( foundFace );
483 Handle(Geom_Surface) surface = BRep_Tool::Surface( theFace );
484 GeomAPI_ProjectPointOnSurf projector( thePoint, surface );
485 if ( !projector.IsDone() || projector.NbPoints()==0 )
486 throw SMESH_ComputeError(COMPERR_BAD_PARMETERS,
487 "getProjectionPoint: can't project a vertex to a face");
489 Quantity_Parameter u,v;
490 projector.LowerDistanceParameters(u,v);
491 myPoint.uv = gp_XY(u,v);
492 gp_Pnt aPnt = projector.NearestPoint();
493 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
495 BRepClass_FaceClassifier FC( theFace, myPoint.uv, Precision::Confusion());
496 if ( FC.State() != TopAbs_IN )
503 /////////////////////////////////////////////////////////
504 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
507 if ( !entry.empty() )
509 GEOM::GEOM_Object_var aGeomObj;
510 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
511 if (!aSObj->_is_nil()) {
512 CORBA::Object_var obj = aSObj->GetObject();
513 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
516 if ( !aGeomObj->_is_nil() )
517 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
522 void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
524 BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
526 // Find the face and get the (u,v) values of the enforced vertex on the face
527 projectionPoint myPoint = getProjectionPoint(faceShape,aPnt);
528 if ( faceShape.IsNull() )
531 enf_coords.push_back(aPnt.X());
532 enf_coords.push_back(aPnt.Y());
533 enf_coords.push_back(aPnt.Z());
535 coords.push_back(myPoint.uv.X());
536 coords.push_back(myPoint.uv.Y());
537 coords.push_back(myPoint.xyz.X());
538 coords.push_back(myPoint.xyz.Y());
539 coords.push_back(myPoint.xyz.Z());
541 s_coords.push_back(myPoint.xyz.X());
542 s_coords.push_back(myPoint.xyz.Y());
543 s_coords.push_back(myPoint.xyz.Z());
545 // Save pair projected vertex / enf vertex
546 EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
547 pair<BLSURFPlugin_Hypothesis::TEnfVertexList::iterator,bool> ret;
548 BLSURFPlugin_Hypothesis::TEnfVertexList::iterator it;
549 ret = EnfVertexCoords2EnfVertexList[s_coords].insert(enfVertex);
550 if (ret.second == false) {
552 (*it)->grpName = enfVertex->grpName;
556 if (! FacesWithEnforcedVertices.Contains(faceShape)) {
557 key = FacesWithEnforcedVertices.Add(faceShape);
560 key = FacesWithEnforcedVertices.FindIndex(faceShape);
563 // If a node is already created by an attractor, do not create enforced vertex
564 int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
565 bool sameAttractor = false;
566 if (attractorKey >= 0)
567 if (FaceId2AttractorCoords.count(attractorKey) > 0)
568 if (FaceId2AttractorCoords[attractorKey] == coords)
569 sameAttractor = true;
571 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
573 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
576 if (! sameAttractor) {
577 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList ens;
579 FaceId2EnforcedVertexCoords[key] = ens;
584 /////////////////////////////////////////////////////////
585 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
587 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
590 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
592 for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
593 enfVertex = *enfVertexListIt;
594 // Case of manual coords
595 if (enfVertex->coords.size() != 0) {
596 aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
597 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
600 // Case of geom vertex coords
601 if (enfVertex->geomEntry != "") {
602 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
603 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
604 if (GeomType == TopAbs_VERTEX)
606 enfVertex->vertex = TopoDS::Vertex( GeomShape );
607 aPnt = BRep_Tool::Pnt( enfVertex->vertex );
608 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
611 if (GeomType == TopAbs_COMPOUND)
613 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next())
614 if (it.Value().ShapeType() == TopAbs_VERTEX)
616 enfVertex->vertex = TopoDS::Vertex( it.Value() );
617 aPnt = BRep_Tool::Pnt( enfVertex->vertex );
618 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
625 /////////////////////////////////////////////////////////
626 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction, double defaultSize)
628 double xa, ya, za; // Coordinates of attractor point
629 double a, b; // Attractor parameter
631 bool createNode=false; // To create a node on attractor projection
633 const char *sep = ";";
634 // atIt->second has the following pattern:
635 // ATTRACTOR(xa;ya;za;a;b;True|False;d)
637 // xa;ya;za : coordinates of attractor
638 // a : desired size on attractor
639 // b : distance of influence of attractor
640 // d : distance until which the size remains constant
642 // We search the parameters in the string
644 pos1 = AttractorFunction.find(sep);
645 if (pos1!=string::npos)
646 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
648 pos2 = AttractorFunction.find(sep, pos1+1);
649 if (pos2!=string::npos) {
650 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
654 pos2 = AttractorFunction.find(sep, pos1+1);
655 if (pos2!=string::npos) {
656 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
660 pos2 = AttractorFunction.find(sep, pos1+1);
661 if (pos2!=string::npos) {
662 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
666 pos2 = AttractorFunction.find(sep, pos1+1);
667 if (pos2!=string::npos) {
668 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
672 pos2 = AttractorFunction.find(sep, pos1+1);
673 if (pos2!=string::npos) {
674 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
675 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
679 pos2 = AttractorFunction.find(")");
680 if (pos2!=string::npos) {
681 d = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
684 // Get the (u,v) values of the attractor on the face
685 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
686 gp_XY uvPoint = myPoint.uv;
687 gp_XYZ xyzPoint = myPoint.xyz;
688 Standard_Real u0 = uvPoint.X();
689 Standard_Real v0 = uvPoint.Y();
690 Standard_Real x0 = xyzPoint.X();
691 Standard_Real y0 = xyzPoint.Y();
692 Standard_Real z0 = xyzPoint.Z();
693 std::vector<double> coords;
694 coords.push_back(u0);
695 coords.push_back(v0);
696 coords.push_back(x0);
697 coords.push_back(y0);
698 coords.push_back(z0);
699 // We construct the python function
700 ostringstream attractorFunctionStream;
701 attractorFunctionStream << "def f(u,v): return ";
702 attractorFunctionStream << defaultSize << "-(" << defaultSize <<"-" << a << ")";
703 //attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
704 // rnc: make possible to keep the size constant until
705 // a defined distance. Distance is expressed as the positiv part
706 // of r-d where r is the distance to (u0,v0)
707 attractorFunctionStream << "*exp(-(0.5*(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"+abs(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"))/(" << b << "))**2)";
710 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
711 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
714 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
716 FaceId2SizeMap[key] =attractorFunctionStream.str();
718 FaceId2AttractorCoords[key] = coords;
720 // // Test for new attractors
721 // gp_Pnt myP(xyzPoint);
722 // TopoDS_Vertex myV = BRepBuilderAPI_MakeVertex(myP);
723 // BLSURFPlugin_Attractor myAttractor(TopoDS::Face(GeomShape),myV,200);
724 // myAttractor.SetParameters(a, defaultSize, b, d);
725 // myAttractor.SetType(1);
726 // FaceId2ClassAttractor[key] = myAttractor;
727 // if(FaceId2ClassAttractor[key].GetFace().IsNull()){
731 // One sub-shape to get ids from
732 BLSURFPlugin_BLSURF::TListOfIDs _getSubShapeIDsInMainShape(const TopoDS_Shape& theMainShape,
733 const TopoDS_Shape& theSubShape,
734 TopAbs_ShapeEnum theShapeType)
736 BLSURFPlugin_BLSURF::TListOfIDs face_ids;
738 TopTools_MapOfShape subShapes;
739 TopTools_IndexedMapOfShape anIndices;
740 TopExp::MapShapes(theMainShape, theShapeType, anIndices);
742 for (TopExp_Explorer face_iter(theSubShape,theShapeType);face_iter.More();face_iter.Next())
744 if ( subShapes.Add( face_iter.Current() )) // issue 23416
746 int face_id = anIndices.FindIndex( face_iter.Current() );
748 throw SALOME_Exception( "Periodicity: sub_shape not found in main_shape");
749 face_ids.push_back( face_id );
755 BLSURFPlugin_BLSURF::TListOfIDs _getSubShapeIDsInMainShape(SMESH_Mesh* theMesh,
756 TopoDS_Shape theSubShape,
757 TopAbs_ShapeEnum theShapeType)
759 BLSURFPlugin_BLSURF::TListOfIDs face_ids;
761 for (TopExp_Explorer face_iter(theSubShape,theShapeType);face_iter.More();face_iter.Next())
763 int face_id = theMesh->GetMeshDS()->ShapeToIndex(face_iter.Current());
765 throw SALOME_Exception ( "Periodicity: sub_shape not found in main_shape");
766 face_ids.push_back(face_id);
772 void BLSURFPlugin_BLSURF::addCoordsFromVertices(const std::vector<std::string> &theVerticesEntries, std::vector<double> &theVerticesCoords)
774 for (std::vector<std::string>::const_iterator it = theVerticesEntries.begin(); it != theVerticesEntries.end(); it++)
776 BLSURFPlugin_Hypothesis::TEntry theVertexEntry = *it;
777 addCoordsFromVertex(theVertexEntry, theVerticesCoords);
782 void BLSURFPlugin_BLSURF::addCoordsFromVertex(BLSURFPlugin_Hypothesis::TEntry theVertexEntry, std::vector<double> &theVerticesCoords)
784 if (theVertexEntry!="")
786 TopoDS_Shape aShape = entryToShape(theVertexEntry);
788 gp_Pnt aPnt = BRep_Tool::Pnt( TopoDS::Vertex( aShape ) );
789 double theX, theY, theZ;
794 theVerticesCoords.push_back(theX);
795 theVerticesCoords.push_back(theY);
796 theVerticesCoords.push_back(theZ);
800 /////////////////////////////////////////////////////////
801 void BLSURFPlugin_BLSURF::createPreCadFacesPeriodicity(TopoDS_Shape theGeomShape, const BLSURFPlugin_Hypothesis::TPreCadPeriodicity &preCadPeriodicity)
803 TopoDS_Shape geomShape1 = entryToShape(preCadPeriodicity.shape1Entry);
804 TopoDS_Shape geomShape2 = entryToShape(preCadPeriodicity.shape2Entry);
806 TListOfIDs theFace1_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape1, TopAbs_FACE);
807 TListOfIDs theFace2_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape2, TopAbs_FACE);
809 TPreCadPeriodicityIDs preCadFacesPeriodicityIDs;
810 preCadFacesPeriodicityIDs.shape1IDs = theFace1_ids;
811 preCadFacesPeriodicityIDs.shape2IDs = theFace2_ids;
813 addCoordsFromVertices(preCadPeriodicity.theSourceVerticesEntries, preCadFacesPeriodicityIDs.theSourceVerticesCoords);
814 addCoordsFromVertices(preCadPeriodicity.theTargetVerticesEntries, preCadFacesPeriodicityIDs.theTargetVerticesCoords);
816 _preCadFacesIDsPeriodicityVector.push_back(preCadFacesPeriodicityIDs);
819 /////////////////////////////////////////////////////////
820 void BLSURFPlugin_BLSURF::createPreCadEdgesPeriodicity(TopoDS_Shape theGeomShape, const BLSURFPlugin_Hypothesis::TPreCadPeriodicity &preCadPeriodicity)
822 TopoDS_Shape geomShape1 = entryToShape(preCadPeriodicity.shape1Entry);
823 TopoDS_Shape geomShape2 = entryToShape(preCadPeriodicity.shape2Entry);
825 TListOfIDs theEdge1_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape1, TopAbs_EDGE);
826 TListOfIDs theEdge2_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape2, TopAbs_EDGE);
828 TPreCadPeriodicityIDs preCadEdgesPeriodicityIDs;
829 preCadEdgesPeriodicityIDs.shape1IDs = theEdge1_ids;
830 preCadEdgesPeriodicityIDs.shape2IDs = theEdge2_ids;
832 addCoordsFromVertices(preCadPeriodicity.theSourceVerticesEntries, preCadEdgesPeriodicityIDs.theSourceVerticesCoords);
833 addCoordsFromVertices(preCadPeriodicity.theTargetVerticesEntries, preCadEdgesPeriodicityIDs.theTargetVerticesCoords);
835 _preCadEdgesIDsPeriodicityVector.push_back(preCadEdgesPeriodicityIDs);
839 /////////////////////////////////////////////////////////
841 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
842 cadsurf_session_t * css,
843 const TopoDS_Shape& theGeomShape)
846 // Clear map so that it is not stored in the algorithm with old enforced vertices in it
847 FacesWithSizeMap.Clear();
848 FaceId2SizeMap.clear();
849 EdgesWithSizeMap.Clear();
850 EdgeId2SizeMap.clear();
851 VerticesWithSizeMap.Clear();
852 VertexId2SizeMap.clear();
853 FaceId2PythonSmp.clear();
854 EdgeId2PythonSmp.clear();
855 VertexId2PythonSmp.clear();
856 FaceId2AttractorCoords.clear();
857 FaceId2ClassAttractor.clear();
858 FaceIndex2ClassAttractor.clear();
859 FacesWithEnforcedVertices.Clear();
860 FaceId2EnforcedVertexCoords.clear();
861 EnfVertexCoords2ProjVertex.clear();
862 EnfVertexCoords2EnfVertexList.clear();
864 double diagonal = SMESH_Mesh::GetShapeDiagonalSize( theGeomShape );
865 double bbSegmentation = _gen->GetBoundaryBoxSegmentation();
866 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
867 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
868 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize(diagonal, bbSegmentation);
869 bool _phySizeRel = BLSURFPlugin_Hypothesis::GetDefaultPhySizeRel();
870 double _minSize = BLSURFPlugin_Hypothesis::GetDefaultMinSize(diagonal);
871 bool _minSizeRel = BLSURFPlugin_Hypothesis::GetDefaultMinSizeRel();
872 double _maxSize = BLSURFPlugin_Hypothesis::GetDefaultMaxSize(diagonal);
873 bool _maxSizeRel = BLSURFPlugin_Hypothesis::GetDefaultMaxSizeRel();
874 double _use_gradation = BLSURFPlugin_Hypothesis::GetDefaultUseGradation();
875 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
876 double _use_volume_gradation = BLSURFPlugin_Hypothesis::GetDefaultUseVolumeGradation();
877 double _volume_gradation = BLSURFPlugin_Hypothesis::GetDefaultVolumeGradation();
878 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
879 double _angleMesh = BLSURFPlugin_Hypothesis::GetDefaultAngleMesh();
880 double _chordalError = BLSURFPlugin_Hypothesis::GetDefaultChordalError(diagonal);
881 bool _anisotropic = BLSURFPlugin_Hypothesis::GetDefaultAnisotropic();
882 double _anisotropicRatio = BLSURFPlugin_Hypothesis::GetDefaultAnisotropicRatio();
883 bool _removeTinyEdges = BLSURFPlugin_Hypothesis::GetDefaultRemoveTinyEdges();
884 double _tinyEdgeLength = BLSURFPlugin_Hypothesis::GetDefaultTinyEdgeLength(diagonal);
885 bool _optimiseTinyEdges = BLSURFPlugin_Hypothesis::GetDefaultOptimiseTinyEdges();
886 double _tinyEdgeOptimisLength = BLSURFPlugin_Hypothesis::GetDefaultTinyEdgeOptimisationLength(diagonal);
887 bool _correctSurfaceIntersec= BLSURFPlugin_Hypothesis::GetDefaultCorrectSurfaceIntersection();
888 double _corrSurfaceIntersCost = BLSURFPlugin_Hypothesis::GetDefaultCorrectSurfaceIntersectionMaxCost();
889 bool _badElementRemoval = BLSURFPlugin_Hypothesis::GetDefaultBadElementRemoval();
890 double _badElementAspectRatio = BLSURFPlugin_Hypothesis::GetDefaultBadElementAspectRatio();
891 bool _optimizeMesh = BLSURFPlugin_Hypothesis::GetDefaultOptimizeMesh();
892 bool _quadraticMesh = BLSURFPlugin_Hypothesis::GetDefaultQuadraticMesh();
893 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
894 //int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
897 //int _precadMergeEdges = BLSURFPlugin_Hypothesis::GetDefaultPreCADMergeEdges();
898 //int _precadRemoveDuplicateCADFaces = BLSURFPlugin_Hypothesis::GetDefaultPreCADRemoveDuplicateCADFaces();
899 //int _precadProcess3DTopology = BLSURFPlugin_Hypothesis::GetDefaultPreCADProcess3DTopology();
900 //int _precadDiscardInput = BLSURFPlugin_Hypothesis::GetDefaultPreCADDiscardInput();
902 const BLSURFPlugin_Hypothesis::TPreCadPeriodicityVector preCadFacesPeriodicityVector = BLSURFPlugin_Hypothesis::GetPreCadFacesPeriodicityVector(hyp);
905 _physicalMesh = (int) hyp->GetPhysicalMesh();
906 _geometricMesh = (int) hyp->GetGeometricMesh();
907 if (hyp->GetPhySize() > 0) {
908 _phySize = hyp->GetPhySize();
909 // if user size is not explicitly specified, "relative" flag is ignored
910 _phySizeRel = hyp->IsPhySizeRel();
912 if (hyp->GetMinSize() > 0) {
913 _minSize = hyp->GetMinSize();
914 // if min size is not explicitly specified, "relative" flag is ignored
915 _minSizeRel = hyp->IsMinSizeRel();
917 if (hyp->GetMaxSize() > 0) {
918 _maxSize = hyp->GetMaxSize();
919 // if max size is not explicitly specified, "relative" flag is ignored
920 _maxSizeRel = hyp->IsMaxSizeRel();
922 _use_gradation = hyp->GetUseGradation();
923 if (hyp->GetGradation() > 0 && _use_gradation)
924 _gradation = hyp->GetGradation();
925 _use_volume_gradation = hyp->GetUseVolumeGradation();
926 if (hyp->GetVolumeGradation() > 0 && _use_volume_gradation )
927 _volume_gradation = hyp->GetVolumeGradation();
928 _quadAllowed = hyp->GetQuadAllowed();
929 if (hyp->GetAngleMesh() > 0)
930 _angleMesh = hyp->GetAngleMesh();
931 if (hyp->GetChordalError() > 0)
932 _chordalError = hyp->GetChordalError();
933 _anisotropic = hyp->GetAnisotropic();
934 if (hyp->GetAnisotropicRatio() >= 0)
935 _anisotropicRatio = hyp->GetAnisotropicRatio();
936 _removeTinyEdges = hyp->GetRemoveTinyEdges();
937 if (hyp->GetTinyEdgeLength() > 0)
938 _tinyEdgeLength = hyp->GetTinyEdgeLength();
939 _optimiseTinyEdges = hyp->GetOptimiseTinyEdges();
940 if (hyp->GetTinyEdgeOptimisationLength() > 0)
941 _tinyEdgeOptimisLength = hyp->GetTinyEdgeOptimisationLength();
942 _correctSurfaceIntersec = hyp->GetCorrectSurfaceIntersection();
943 if (hyp->GetCorrectSurfaceIntersectionMaxCost() > 0)
944 _corrSurfaceIntersCost = hyp->GetCorrectSurfaceIntersectionMaxCost();
945 _badElementRemoval = hyp->GetBadElementRemoval();
946 if (hyp->GetBadElementAspectRatio() >= 0)
947 _badElementAspectRatio = hyp->GetBadElementAspectRatio();
948 _optimizeMesh = hyp->GetOptimizeMesh();
949 _quadraticMesh = hyp->GetQuadraticMesh();
950 _verb = hyp->GetVerbosity();
951 //_topology = (int) hyp->GetTopology();
953 //_precadMergeEdges = hyp->GetPreCADMergeEdges();
954 //_precadRemoveDuplicateCADFaces = hyp->GetPreCADRemoveDuplicateCADFaces();
955 //_precadProcess3DTopology = hyp->GetPreCADProcess3DTopology();
956 //_precadDiscardInput = hyp->GetPreCADDiscardInput();
958 const BLSURFPlugin_Hypothesis::TOptionValues& opts = hyp->GetOptionValues();
959 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
960 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt ){
961 MESSAGE("OptionValue: " << opIt->first.c_str() << ", value: " << opIt->second.c_str());
962 if ( !opIt->second.empty() ) {
963 // With MeshGems 2.4-5, there are issues with periodicity and multithread
964 // => As a temporary workaround, we enforce to use only one thread if periodicity is used.
965 if (opIt->first == "max_number_of_threads" && opIt->second != "1" && ! preCadFacesPeriodicityVector.empty()){
966 std::cout << "INFO: Disabling multithread to avoid periodicity issues" << std::endl;
967 set_param(css, opIt->first.c_str(), "1");
970 set_param(css, opIt->first.c_str(), opIt->second.c_str());
974 const BLSURFPlugin_Hypothesis::TOptionValues& custom_opts = hyp->GetCustomOptionValues();
975 for ( opIt = custom_opts.begin(); opIt != custom_opts.end(); ++opIt )
976 if ( !opIt->second.empty() ) {
977 set_param(css, opIt->first.c_str(), opIt->second.c_str());
980 const BLSURFPlugin_Hypothesis::TOptionValues& preCADopts = hyp->GetPreCADOptionValues();
981 for ( opIt = preCADopts.begin(); opIt != preCADopts.end(); ++opIt )
982 if ( !opIt->second.empty() ) {
983 set_param(css, opIt->first.c_str(), opIt->second.c_str());
987 if ( BLSURFPlugin_Hypothesis::HasPreCADOptions( hyp ))
989 cadsurf_set_param(css, "use_precad", "yes" ); // for old versions
991 // PreProcessor (formerly PreCAD) -- commented params are preCADoptions (since 0023307)
992 //set_param(css, "merge_edges", _precadMergeEdges ? "yes" : "no");
993 //set_param(css, "remove_duplicate_cad_faces", _precadRemoveDuplicateCADFaces ? "yes" : "no");
994 //set_param(css, "process_3d_topology", _precadProcess3DTopology ? "1" : "0");
995 //set_param(css, "discard_input_topology", _precadDiscardInput ? "1" : "0");
996 //set_param(css, "max_number_of_points_per_patch", "1000000");
998 bool useGradation = false;
999 switch (_physicalMesh)
1001 case BLSURFPlugin_Hypothesis::PhysicalGlobalSize:
1002 set_param(css, "physical_size_mode", "global");
1003 set_param(css, "global_physical_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1005 case BLSURFPlugin_Hypothesis::PhysicalLocalSize:
1006 set_param(css, "physical_size_mode", "local");
1007 set_param(css, "global_physical_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1008 useGradation = true;
1011 set_param(css, "physical_size_mode", "none");
1014 switch (_geometricMesh)
1016 case BLSURFPlugin_Hypothesis::GeometricalGlobalSize:
1017 set_param(css, "geometric_size_mode", "global");
1018 set_param(css, "geometric_approximation", val_to_string(_angleMesh).c_str());
1019 set_param(css, "chordal_error", val_to_string(_chordalError).c_str());
1020 useGradation = true;
1022 case BLSURFPlugin_Hypothesis::GeometricalLocalSize:
1023 set_param(css, "geometric_size_mode", "local");
1024 set_param(css, "geometric_approximation", val_to_string(_angleMesh).c_str());
1025 set_param(css, "chordal_error", val_to_string(_chordalError).c_str());
1026 useGradation = true;
1029 set_param(css, "geometric_size_mode", "none");
1032 if ( hyp && hyp->GetPhySize() > 0 ) {
1033 // user size is explicitly specified via hypothesis parameters
1034 // min and max sizes should be compared with explicitly specified user size
1035 // - compute absolute min size
1036 double mins = _minSizeRel ? _minSize * diagonal : _minSize;
1037 // - min size should not be greater than user size
1038 if ( _phySize < mins )
1039 set_param(css, "min_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1041 set_param(css, "min_size", _minSizeRel ? val_to_string_rel(_minSize).c_str() : val_to_string(_minSize).c_str());
1042 // - compute absolute max size
1043 double maxs = _maxSizeRel ? _maxSize * diagonal : _maxSize;
1044 // - max size should not be less than user size
1045 if ( _phySize > maxs )
1046 set_param(css, "max_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1048 set_param(css, "max_size", _maxSizeRel ? val_to_string_rel(_maxSize).c_str() : val_to_string(_maxSize).c_str());
1051 // user size is not explicitly specified
1052 // - if minsize is not explicitly specified, we pass default value computed automatically, in this case "relative" flag is ignored
1053 set_param(css, "min_size", _minSizeRel ? val_to_string_rel(_minSize).c_str() : val_to_string(_minSize).c_str());
1054 // - if maxsize is not explicitly specified, we pass default value computed automatically, in this case "relative" flag is ignored
1055 set_param(css, "max_size", _maxSizeRel ? val_to_string_rel(_maxSize).c_str() : val_to_string(_maxSize).c_str());
1057 // anisotropic and quadrangle mesh requires disabling gradation
1058 if ( _anisotropic && _quadAllowed )
1059 useGradation = false; // limitation of V1.3
1060 if ( useGradation && _use_gradation )
1061 set_param(css, "gradation", val_to_string(_gradation).c_str());
1062 if ( useGradation && _use_volume_gradation )
1063 set_param(css, "volume_gradation", val_to_string(_volume_gradation).c_str());
1064 set_param(css, "element_generation", _quadAllowed ? "quad_dominant" : "triangle");
1067 set_param(css, "metric", _anisotropic ? "anisotropic" : "isotropic");
1069 set_param(css, "anisotropic_ratio", val_to_string(_anisotropicRatio).c_str());
1070 set_param(css, "remove_tiny_edges", _removeTinyEdges ? "1" : "0");
1071 if ( _removeTinyEdges )
1072 set_param(css, "tiny_edge_length", val_to_string(_tinyEdgeLength).c_str());
1073 set_param(css, "optimise_tiny_edges", _optimiseTinyEdges ? "1" : "0");
1074 if ( _optimiseTinyEdges )
1075 set_param(css, "tiny_edge_optimisation_length", val_to_string(_tinyEdgeOptimisLength).c_str());
1076 set_param(css, "correct_surface_intersections", _correctSurfaceIntersec ? "1" : "0");
1077 if ( _correctSurfaceIntersec )
1078 set_param(css, "surface_intersections_processing_max_cost", val_to_string(_corrSurfaceIntersCost ).c_str());
1079 set_param(css, "force_bad_surface_element_removal", _badElementRemoval ? "1" : "0");
1080 if ( _badElementRemoval )
1081 set_param(css, "bad_surface_element_aspect_ratio", val_to_string(_badElementAspectRatio).c_str());
1082 set_param(css, "optimisation", _optimizeMesh ? "yes" : "no");
1083 set_param(css, "element_order", _quadraticMesh ? "quadratic" : "linear");
1084 set_param(css, "verbose", val_to_string(_verb).c_str());
1086 _smp_phy_size = _phySizeRel ? _phySize*diagonal : _phySize;
1088 std::cout << "_smp_phy_size = " << _smp_phy_size << std::endl;
1090 if (_physicalMesh == BLSURFPlugin_Hypothesis::PhysicalLocalSize)
1092 TopoDS_Shape GeomShape;
1093 TopoDS_Shape AttShape;
1094 TopAbs_ShapeEnum GeomType;
1096 // Standard Size Maps
1098 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
1099 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
1100 for ( ; smIt != sizeMaps.end(); ++smIt ) {
1101 if ( !smIt->second.empty() ) {
1102 GeomShape = entryToShape(smIt->first);
1103 GeomType = GeomShape.ShapeType();
1106 if (GeomType == TopAbs_COMPOUND) {
1107 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1109 if (it.Value().ShapeType() == TopAbs_FACE){
1110 HasSizeMapOnFace = true;
1111 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
1112 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
1115 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
1117 FaceId2SizeMap[key] = smIt->second;
1120 if (it.Value().ShapeType() == TopAbs_EDGE){
1121 HasSizeMapOnEdge = true;
1122 HasSizeMapOnFace = true;
1123 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
1124 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
1127 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
1129 EdgeId2SizeMap[key] = smIt->second;
1131 // Group of vertices
1132 if (it.Value().ShapeType() == TopAbs_VERTEX){
1133 HasSizeMapOnVertex = true;
1134 HasSizeMapOnEdge = true;
1135 HasSizeMapOnFace = true;
1136 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
1137 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
1140 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
1142 VertexId2SizeMap[key] = smIt->second;
1147 if (GeomType == TopAbs_FACE){
1148 HasSizeMapOnFace = true;
1149 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
1150 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
1153 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
1155 FaceId2SizeMap[key] = smIt->second;
1158 if (GeomType == TopAbs_EDGE){
1159 HasSizeMapOnEdge = true;
1160 HasSizeMapOnFace = true;
1161 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
1162 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
1165 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
1167 EdgeId2SizeMap[key] = smIt->second;
1170 if (GeomType == TopAbs_VERTEX){
1171 HasSizeMapOnVertex = true;
1172 HasSizeMapOnEdge = true;
1173 HasSizeMapOnFace = true;
1174 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
1175 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
1178 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
1180 VertexId2SizeMap[key] = smIt->second;
1188 // TODO appeler le constructeur des attracteurs directement ici
1189 // if ( !_phySizeRel ) {
1190 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
1191 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
1192 for ( ; atIt != attractors.end(); ++atIt ) {
1193 if ( !atIt->second.empty() ) {
1194 GeomShape = entryToShape(atIt->first);
1195 GeomType = GeomShape.ShapeType();
1197 if (GeomType == TopAbs_COMPOUND){
1198 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1199 if (it.Value().ShapeType() == TopAbs_FACE){
1200 HasSizeMapOnFace = true;
1201 createAttractorOnFace(it.Value(), atIt->second, _phySizeRel ? _phySize*diagonal : _phySize);
1206 if (GeomType == TopAbs_FACE){
1207 HasSizeMapOnFace = true;
1208 createAttractorOnFace(GeomShape, atIt->second, _phySizeRel ? _phySize*diagonal : _phySize);
1211 if (GeomType == TopAbs_EDGE){
1212 HasSizeMapOnEdge = true;
1213 HasSizeMapOnFace = true;
1214 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
1216 if (GeomType == TopAbs_VERTEX){
1217 HasSizeMapOnVertex = true;
1218 HasSizeMapOnEdge = true;
1219 HasSizeMapOnFace = true;
1220 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
1228 // temporary commented out for testing
1230 // - Fill in the BLSURFPlugin_Hypothesis::TAttractorMap map in the hypothesis
1231 // - Uncomment and complete this part to construct the attractors from the attractor shape and the passed parameters on each face of the map
1232 // - To do this use the public methodss: SetParameters(several double parameters) and SetType(int type)
1234 // - Construct the attractors with an empty dist. map in the hypothesis
1235 // - 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()
1236 // -> define a bool _mapbuilt in the class that is set to false by default and set to true when calling _buildmap() OK
1238 theNbAttractors = 0;
1239 const BLSURFPlugin_Hypothesis::TAttractorMap class_attractors = BLSURFPlugin_Hypothesis::GetClassAttractorEntries(hyp);
1241 BLSURFPlugin_Hypothesis::TAttractorMap::const_iterator AtIt = class_attractors.begin();
1242 for ( ; AtIt != class_attractors.end(); ++AtIt ) {
1243 if ( !AtIt->second->Empty() ) {
1244 GeomShape = entryToShape(AtIt->first);
1245 if ( !SMESH_MesherHelper::IsSubShape( GeomShape, theGeomShape ))
1247 AttShape = AtIt->second->GetAttractorShape();
1248 GeomType = GeomShape.ShapeType();
1250 // if (GeomType == TopAbs_COMPOUND){
1251 // for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1252 // if (it.Value().ShapeType() == TopAbs_FACE){
1253 // HasAttractorOnFace = true;
1254 // myAttractor = BLSURFPluginAttractor(GeomShape, AttShape);
1259 if (GeomType == TopAbs_FACE
1260 && (AttShape.ShapeType() == TopAbs_VERTEX
1261 || AttShape.ShapeType() == TopAbs_EDGE
1262 || AttShape.ShapeType() == TopAbs_WIRE
1263 || AttShape.ShapeType() == TopAbs_COMPOUND) ){
1264 HasSizeMapOnFace = true;
1266 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
1268 FaceId2ClassAttractor[key].push_back( AtIt->second );
1272 MESSAGE("Wrong shape type !!")
1280 // Enforced Vertices
1282 const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
1283 BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
1284 for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
1285 if ( !enfIt->second.empty() ) {
1286 GeomShape = entryToShape(enfIt->first);
1287 if ( GeomShape.IsNull() )
1289 createEnforcedVertexOnFace( GeomShape, enfIt->second );
1292 else if ( GeomShape.ShapeType() == TopAbs_COMPOUND)
1294 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1295 if (it.Value().ShapeType() == TopAbs_FACE){
1296 HasSizeMapOnFace = true;
1297 createEnforcedVertexOnFace(it.Value(), enfIt->second);
1301 else if ( GeomShape.ShapeType() == TopAbs_FACE)
1303 HasSizeMapOnFace = true;
1304 createEnforcedVertexOnFace(GeomShape, enfIt->second);
1309 // Internal vertices
1310 bool useInternalVertexAllFaces = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFaces(hyp);
1311 if (useInternalVertexAllFaces) {
1312 std::string grpName = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFacesGroup(hyp);
1314 TopExp_Explorer exp (theGeomShape, TopAbs_FACE);
1315 for (; exp.More(); exp.Next()){
1316 TopExp_Explorer exp_face (exp.Current(), TopAbs_VERTEX, TopAbs_EDGE);
1317 for (; exp_face.More(); exp_face.Next())
1319 // Get coords of vertex
1320 // Check if current coords is already in enfVertexList
1321 // If coords not in enfVertexList, add new enfVertex
1322 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(exp_face.Current()));
1323 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex = new BLSURFPlugin_Hypothesis::TEnfVertex();
1324 enfVertex->coords.push_back(aPnt.X());
1325 enfVertex->coords.push_back(aPnt.Y());
1326 enfVertex->coords.push_back(aPnt.Z());
1327 enfVertex->name = "";
1328 enfVertex->faceEntries.clear();
1329 enfVertex->geomEntry = "";
1330 enfVertex->grpName = grpName;
1331 enfVertex->vertex = TopoDS::Vertex( exp_face.Current() );
1332 _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()), aPnt, enfVertex);
1333 HasSizeMapOnFace = true;
1338 cadsurf_set_sizemap_iso_cad_face(css, size_on_surface, &_smp_phy_size);
1340 if (HasSizeMapOnEdge){
1341 cadsurf_set_sizemap_iso_cad_edge(css, size_on_edge, &_smp_phy_size);
1343 if (HasSizeMapOnVertex){
1344 cadsurf_set_sizemap_iso_cad_point(css, size_on_vertex, &_smp_phy_size);
1351 _preCadFacesIDsPeriodicityVector.clear();
1352 _preCadEdgesIDsPeriodicityVector.clear();
1354 for (std::size_t i = 0; i<preCadFacesPeriodicityVector.size(); i++){
1355 createPreCadFacesPeriodicity(theGeomShape, preCadFacesPeriodicityVector[i]);
1358 const BLSURFPlugin_Hypothesis::TPreCadPeriodicityVector preCadEdgesPeriodicityVector = BLSURFPlugin_Hypothesis::GetPreCadEdgesPeriodicityVector(hyp);
1360 for (std::size_t i = 0; i<preCadEdgesPeriodicityVector.size(); i++){
1361 createPreCadEdgesPeriodicity(theGeomShape, preCadEdgesPeriodicityVector[i]);
1365 //================================================================================
1367 * \brief Throws an exception if a parameter name is wrong
1369 //================================================================================
1371 void BLSURFPlugin_BLSURF::set_param(cadsurf_session_t *css,
1372 const char * option_name,
1373 const char * option_value)
1375 status_t status = cadsurf_set_param(css, option_name, option_value );
1377 if ( _hypothesis && _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
1378 cout << option_name << " = " << option_value << endl;
1380 if ( status != MESHGEMS_STATUS_OK )
1382 if ( status == MESHGEMS_STATUS_UNKNOWN_PARAMETER ) {
1383 throw SALOME_Exception
1384 ( SMESH_Comment("Invalid name of CADSURF parameter: ") << option_name );
1386 else if ( status == MESHGEMS_STATUS_NOLICENSE )
1387 throw SALOME_Exception
1388 ( "No valid license available" );
1390 throw SALOME_Exception
1391 ( SMESH_Comment("Either wrong name or unacceptable value of CADSURF parameter '")
1392 << option_name << "': " << option_value);
1398 // --------------------------------------------------------------------------
1400 * \brief Class correctly terminating usage of MG-CADSurf library at destruction
1402 struct BLSURF_Cleaner
1405 cadsurf_session_t* _css;
1409 BLSURF_Cleaner(context_t * ctx,
1410 cadsurf_session_t* css=0,
1421 Clean( /*exceptContext=*/false );
1423 void Clean(const bool exceptContext)
1427 cadsurf_session_delete(_css); _css = 0;
1429 // #if BLSURF_VERSION_LONG >= "3.1.1"
1430 // // if(geo_sizemap_e)
1431 // // distene_sizemap_delete(geo_sizemap_e);
1432 // // if(geo_sizemap_f)
1433 // // distene_sizemap_delete(geo_sizemap_f);
1434 // if(iso_sizemap_p)
1435 // distene_sizemap_delete(iso_sizemap_p);
1436 // if(iso_sizemap_e)
1437 // distene_sizemap_delete(iso_sizemap_e);
1438 // if(iso_sizemap_f)
1439 // distene_sizemap_delete(iso_sizemap_f);
1441 // // if(clean_geo_sizemap_e)
1442 // // distene_sizemap_delete(clean_geo_sizemap_e);
1443 // // if(clean_geo_sizemap_f)
1444 // // distene_sizemap_delete(clean_geo_sizemap_f);
1445 // if(clean_iso_sizemap_p)
1446 // distene_sizemap_delete(clean_iso_sizemap_p);
1447 // if(clean_iso_sizemap_e)
1448 // distene_sizemap_delete(clean_iso_sizemap_e);
1449 // if(clean_iso_sizemap_f)
1450 // distene_sizemap_delete(clean_iso_sizemap_f);
1453 cad_delete(_cad); _cad = 0;
1454 dcad_delete(_dcad); _dcad = 0;
1455 if ( !exceptContext )
1457 context_delete(_ctx); _ctx = 0;
1463 // --------------------------------------------------------------------------
1464 // comparator to sort nodes and sub-meshes
1465 struct ShapeTypeCompare
1467 // sort nodes by position in the following order:
1468 // SMDS_TOP_FACE=2, SMDS_TOP_EDGE=1, SMDS_TOP_VERTEX=0, SMDS_TOP_3DSPACE=3
1469 bool operator()( const SMDS_MeshNode* n1, const SMDS_MeshNode* n2 ) const
1471 // NEW ORDER: nodes earlier added to sub-mesh are considered "less"
1472 return n1->getIdInShape() < n2->getIdInShape();
1473 // SMDS_TypeOfPosition pos1 = n1->GetPosition()->GetTypeOfPosition();
1474 // SMDS_TypeOfPosition pos2 = n2->GetPosition()->GetTypeOfPosition();
1475 // if ( pos1 == pos2 ) return 0;
1476 // if ( pos1 < pos2 || pos1 == SMDS_TOP_3DSPACE ) return 1;
1479 // sort sub-meshes in order: EDGE, VERTEX
1480 bool operator()( const SMESHDS_SubMesh* s1, const SMESHDS_SubMesh* s2 ) const
1482 int isVertex1 = ( s1 && s1->NbElements() == 0 );
1483 int isVertex2 = ( s2 && s2->NbElements() == 0 );
1484 if ( isVertex1 == isVertex2 )
1486 return isVertex1 < isVertex2;
1490 //================================================================================
1492 * \brief Fills groups of nodes to be merged
1494 //================================================================================
1496 void getNodeGroupsToMerge( const SMESHDS_SubMesh* smDS,
1497 const TopoDS_Shape& shape,
1498 SMESH_MeshEditor::TListOfListOfNodes& nodeGroupsToMerge)
1500 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1501 switch ( shape.ShapeType() )
1503 case TopAbs_VERTEX: {
1504 std::list< const SMDS_MeshNode* > nodes;
1505 while ( nIt->more() )
1506 nodes.push_back( nIt->next() );
1507 if ( nodes.size() > 1 )
1508 nodeGroupsToMerge.push_back( nodes );
1512 std::multimap< double, const SMDS_MeshNode* > u2node;
1513 const SMDS_EdgePosition* ePos;
1514 while ( nIt->more() )
1516 const SMDS_MeshNode* n = nIt->next();
1517 if (( ePos = dynamic_cast< const SMDS_EdgePosition* >( n->GetPosition() )))
1518 u2node.insert( make_pair( ePos->GetUParameter(), n ));
1520 if ( u2node.size() < 2 ) return;
1522 //double tol = (( u2node.rbegin()->first - u2node.begin()->first ) / 20.) / u2node.size();
1524 BRep_Tool::Range( TopoDS::Edge( shape ), f,l );
1525 double tol = (( l - f ) / 20.) / u2node.size();
1527 std::multimap< double, const SMDS_MeshNode* >::iterator un2, un1;
1528 for ( un2 = u2node.begin(), un1 = un2++; un2 != u2node.end(); un1 = un2++ )
1530 if (( un2->first - un1->first ) <= tol )
1532 std::list< const SMDS_MeshNode* > nodes;
1533 nodes.push_back( un1->second );
1534 while (( un2->first - un1->first ) <= tol )
1536 nodes.push_back( un2->second );
1537 if ( ++un2 == u2node.end()) {
1542 // make nodes created on the boundary of viscous layer replace nodes created
1543 // by MG-CADSurf as their SMDS_Position is more correct
1544 nodes.sort( ShapeTypeCompare() );
1545 nodeGroupsToMerge.push_back( nodes );
1552 // SMESH_MeshEditor::TListOfListOfNodes::const_iterator nll = nodeGroupsToMerge.begin();
1553 // for ( ; nll != nodeGroupsToMerge.end(); ++nll )
1555 // cout << "Merge ";
1556 // const std::list< const SMDS_MeshNode* >& nl = *nll;
1557 // std::list< const SMDS_MeshNode* >::const_iterator nIt = nl.begin();
1558 // for ( ; nIt != nl.end(); ++nIt )
1559 // cout << (*nIt) << " ";
1565 //================================================================================
1567 * \brief A temporary mesh used to compute mesh on a proxy FACE
1569 //================================================================================
1571 struct TmpMesh: public SMESH_Mesh
1573 typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
1574 TN2NMap _tmp2origNN;
1575 TopoDS_Face _proxyFace;
1579 _myMeshDS = new SMESHDS_Mesh( _id, true );
1581 //--------------------------------------------------------------------------------
1583 * \brief Creates a FACE bound by viscous layers and mesh each its EDGE with 1 segment
1585 //--------------------------------------------------------------------------------
1587 const TopoDS_Face& makeProxyFace( SMESH_ProxyMesh::Ptr& viscousMesh,
1588 const TopoDS_Face& origFace)
1590 SMESH_Mesh* origMesh = viscousMesh->GetMesh();
1592 SMESH_MesherHelper helper( *origMesh );
1593 helper.SetSubShape( origFace );
1594 const bool hasSeam = helper.HasRealSeam();
1596 // get data of nodes on inner boundary of viscous layers
1598 TSideVector wireVec = StdMeshers_FaceSide::GetFaceWires(origFace, *origMesh,
1599 /*skipMediumNodes = */true,
1600 err, &helper, viscousMesh );
1601 if ( err && err->IsKO() )
1602 throw *err.get(); // it should be caught at SMESH_subMesh
1604 // proxy nodes and corresponding tmp VERTEXes
1605 std::vector<const SMDS_MeshNode*> origNodes;
1606 std::vector<TopoDS_Vertex> tmpVertex;
1608 // create a proxy FACE
1609 TopoDS_Face origFaceCopy = TopoDS::Face( origFace.EmptyCopied() );
1610 BRepBuilderAPI_MakeFace newFace( origFaceCopy );
1611 bool hasPCurves = false;
1612 for ( size_t iW = 0; iW != wireVec.size(); ++iW )
1614 StdMeshers_FaceSidePtr& wireData = wireVec[iW];
1615 const UVPtStructVec& wirePoints = wireData->GetUVPtStruct();
1616 if ( wirePoints.size() < 3 )
1619 BRepBuilderAPI_MakePolygon polygon;
1620 const size_t i0 = tmpVertex.size();
1621 for ( size_t iN = 0; iN < wirePoints.size(); ++iN )
1623 polygon.Add( SMESH_TNodeXYZ( wirePoints[ iN ].node ));
1624 origNodes.push_back( wirePoints[ iN ].node );
1625 tmpVertex.push_back( polygon.LastVertex() );
1627 // check presence of a pcurve
1628 checkPCurve( polygon, origFaceCopy, hasPCurves, &wirePoints[ iN-1 ] );
1630 tmpVertex[ i0 ] = polygon.FirstVertex(); // polygon.LastVertex()==NULL for 1 vertex in wire
1632 if ( !polygon.IsDone() )
1633 throw SALOME_Exception("BLSURFPlugin_BLSURF: BRepBuilderAPI_MakePolygon failed");
1634 TopoDS_Wire wire = polygon;
1636 wire = updateSeam( wire, origNodes );
1637 newFace.Add( wire );
1639 _proxyFace = newFace;
1641 // set a new shape to mesh
1642 TopoDS_Compound auxCompoundToMesh;
1643 BRep_Builder shapeBuilder;
1644 shapeBuilder.MakeCompound( auxCompoundToMesh );
1645 shapeBuilder.Add( auxCompoundToMesh, _proxyFace );
1646 shapeBuilder.Add( auxCompoundToMesh, origMesh->GetShapeToMesh() );
1648 ShapeToMesh( auxCompoundToMesh );
1651 // Make input mesh for MG-CADSurf: segments on EDGE's of newFace
1653 // make nodes and fill in _tmp2origNN
1655 SMESHDS_Mesh* tmpMeshDS = GetMeshDS();
1656 for ( size_t i = 0; i < origNodes.size(); ++i )
1658 GetSubMesh( tmpVertex[i] )->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1659 if ( const SMDS_MeshNode* tmpN = SMESH_Algo::VertexNode( tmpVertex[i], tmpMeshDS ))
1660 _tmp2origNN.insert( _tmp2origNN.end(), make_pair( tmpN, origNodes[i] ));
1661 // else -- it can be a seam vertex replaced by updateSeam()
1662 // throw SALOME_Exception("BLSURFPlugin_BLSURF: a proxy vertex not meshed");
1666 TopoDS_Vertex v1, v2;
1667 for ( TopExp_Explorer edge( _proxyFace, TopAbs_EDGE ); edge.More(); edge.Next() )
1669 const TopoDS_Edge& E = TopoDS::Edge( edge.Current() );
1670 TopExp::Vertices( E, v1, v2 );
1671 const SMDS_MeshNode* n1 = SMESH_Algo::VertexNode( v1, tmpMeshDS );
1672 const SMDS_MeshNode* n2 = SMESH_Algo::VertexNode( v2, tmpMeshDS );
1674 if ( SMDS_MeshElement* seg = tmpMeshDS->AddEdge( n1, n2 ))
1675 tmpMeshDS->SetMeshElementOnShape( seg, E );
1681 //--------------------------------------------------------------------------------
1683 * \brief Add pcurve to the last edge of a wire
1685 //--------------------------------------------------------------------------------
1687 void checkPCurve( BRepBuilderAPI_MakePolygon& wire,
1688 const TopoDS_Face& face,
1690 const uvPtStruct * wirePoints )
1694 TopoDS_Edge edge = wire.Edge();
1695 if ( edge.IsNull() ) return;
1697 if ( BRep_Tool::CurveOnSurface(edge, face, f, l))
1702 gp_XY p1 = wirePoints[ 0 ].UV(), p2 = wirePoints[ 1 ].UV();
1703 Handle(Geom2d_Line) pcurve = new Geom2d_Line( p1, gp_Dir2d( p2 - p1 ));
1704 BRep_Builder().UpdateEdge( edge, Handle(Geom_Curve)(), Precision::Confusion() );
1705 BRep_Builder().UpdateEdge( edge, pcurve, face, Precision::Confusion() );
1706 BRep_Builder().Range( edge, 0, ( p2 - p1 ).Modulus() );
1707 // cout << "n1 = mesh.AddNode( " << p1.X()*10 << ", " << p1.Y() << ", 0 )" << endl
1708 // << "n2 = mesh.AddNode( " << p2.X()*10 << ", " << p2.Y() << ", 0 )" << endl
1709 // << "mesh.AddEdge( [ n1, n2 ] )" << endl;
1712 //--------------------------------------------------------------------------------
1714 * \brief Replace coincident EDGEs with reversed copies.
1716 //--------------------------------------------------------------------------------
1718 TopoDS_Wire updateSeam( const TopoDS_Wire& wire,
1719 const std::vector<const SMDS_MeshNode*>& nodesOfVertices )
1721 BRepBuilderAPI_MakeWire newWire;
1723 typedef NCollection_DataMap<SMESH_TLink, TopoDS_Edge, SMESH_TLink > TSeg2EdgeMap;
1724 TSeg2EdgeMap seg2EdgeMap;
1726 TopoDS_Iterator edgeIt( wire );
1727 for ( int iSeg = 1; edgeIt.More(); edgeIt.Next(), ++iSeg )
1729 SMESH_TLink link( nodesOfVertices[ iSeg-1 ], nodesOfVertices[ iSeg ]);
1730 TopoDS_Edge edge( TopoDS::Edge( edgeIt.Value() ));
1732 TopoDS_Edge* edgeInMap = seg2EdgeMap.Bound( link, edge );
1733 bool isSeam = ( *edgeInMap != edge );
1736 edgeInMap->Reverse();
1739 newWire.Add( edge );
1744 //--------------------------------------------------------------------------------
1746 * \brief Fill in the origMesh with faces computed by MG-CADSurf in this tmp mesh
1748 //--------------------------------------------------------------------------------
1750 void FillInOrigMesh( SMESH_Mesh& origMesh,
1751 const TopoDS_Face& origFace )
1753 SMESH_MesherHelper helper( origMesh );
1754 helper.SetSubShape( origFace );
1755 helper.SetElementsOnShape( true );
1757 SMESH_MesherHelper tmpHelper( *this );
1758 tmpHelper.SetSubShape( _proxyFace );
1760 // iterate over tmp faces and copy them in origMesh
1761 const SMDS_MeshNode* nodes[27];
1762 const SMDS_MeshNode* nullNode = 0;
1764 SMDS_FaceIteratorPtr fIt = GetMeshDS()->facesIterator(/*idInceasingOrder=*/true);
1765 while ( fIt->more() )
1767 const SMDS_MeshElement* f = fIt->next();
1768 SMDS_ElemIteratorPtr nIt = f->nodesIterator();
1770 for ( ; nIt->more(); ++nbN )
1772 const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nIt->next() );
1773 TN2NMap::iterator n2nIt =
1774 _tmp2origNN.insert( _tmp2origNN.end(), make_pair( n, nullNode ));
1775 if ( !n2nIt->second ) {
1777 gp_XY uv = tmpHelper.GetNodeUV( _proxyFace, n );
1778 n2nIt->second = helper.AddNode( xyz[0], xyz[1], xyz[2], uv.X(), uv.Y() );
1780 nodes[ nbN ] = n2nIt->second;
1783 case 3: helper.AddFace( nodes[0], nodes[1], nodes[2] ); break;
1784 // case 6: helper.AddFace( nodes[0], nodes[1], nodes[2],
1785 // nodes[3], nodes[4], nodes[5]); break;
1786 case 4: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3] ); break;
1787 // case 9: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3],
1788 // nodes[4], nodes[5], nodes[6], nodes[7], nodes[8]); break;
1789 // case 8: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3],
1790 // nodes[4], nodes[5], nodes[6], nodes[7]); break;
1797 * \brief A structure holding an error description and a verbisity level
1799 struct message_cb_user_data
1801 std::string * _error;
1808 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
1809 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1810 real *duu, real *duv, real *dvv, void *user_data);
1811 status_t message_cb(message_t *msg, void *user_data);
1812 status_t interrupt_cb(integer *interrupt_status, void *user_data);
1814 //=============================================================================
1818 //=============================================================================
1820 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
1822 // Fix problem with locales
1823 Kernel_Utils::Localizer aLocalizer;
1825 this->SMESH_Algo::_progress = 1e-3; // prevent progress advancment while computing attractors
1827 bool viscousLayersMade =
1828 ( aShape.ShapeType() == TopAbs_FACE &&
1829 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( aShape ), aMesh ));
1831 if ( !viscousLayersMade )
1832 if ( !compute( aMesh, aShape, /*allowSubMeshClearing=*/true ))
1835 if ( _haveViscousLayers || viscousLayersMade )
1837 // Compute viscous layers
1839 TopTools_MapOfShape map;
1840 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1842 const TopoDS_Face& F = TopoDS::Face(face_iter.Current());
1843 if ( !map.Add( F )) continue;
1844 SMESH_ProxyMesh::Ptr viscousMesh = StdMeshers_ViscousLayers2D::Compute( aMesh, F );
1846 return false; // error in StdMeshers_ViscousLayers2D::Compute()
1848 // Compute MG-CADSurf mesh on viscous layers
1850 if ( viscousMesh->NbProxySubMeshes() > 0 )
1853 const TopoDS_Face& proxyFace = tmpMesh.makeProxyFace( viscousMesh, F );
1854 if ( !compute( tmpMesh, proxyFace, /*allowSubMeshClearing=*/false ))
1856 tmpMesh.FillInOrigMesh( aMesh, F );
1860 // Re-compute MG-CADSurf mesh on the rest faces if the mesh was cleared
1862 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1864 const TopoDS_Face& F = TopoDS::Face(face_iter.Current());
1865 SMESH_subMesh* fSM = aMesh.GetSubMesh( F );
1866 if ( fSM->IsMeshComputed() ) continue;
1868 if ( !compute( aMesh, aShape, /*allowSubMeshClearing=*/true ))
1876 //=============================================================================
1880 //=============================================================================
1882 bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh& aMesh,
1883 const TopoDS_Shape& aShape,
1884 bool allowSubMeshClearing)
1886 /* create a distene context (generic object) */
1887 status_t status = STATUS_ERROR;
1889 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1890 SMESH_MesherHelper helper( aMesh ), helperWithShape( aMesh );
1891 myHelper = theHelper = & helperWithShape;
1892 // do not call helper.IsQuadraticSubMesh() because sub-meshes
1893 // may be cleaned and helper.myTLinkNodeMap gets invalid in such a case
1894 bool haveQuadraticSubMesh = helperWithShape.IsQuadraticSubMesh( aShape );
1895 bool quadraticSubMeshAndViscousLayer = false;
1896 bool needMerge = false;
1897 typedef set< SMESHDS_SubMesh*, ShapeTypeCompare > TSubMeshSet;
1898 TSubMeshSet edgeSubmeshes;
1899 TSubMeshSet& mergeSubmeshes = edgeSubmeshes;
1901 TopTools_IndexedMapOfShape pmap, emap, fmap;
1903 TopTools_IndexedDataMapOfShapeListOfShape e2ffmap;
1904 TopExp::MapShapesAndAncestors( aShape, TopAbs_EDGE, TopAbs_FACE, e2ffmap );
1906 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1908 feclearexcept( FE_ALL_EXCEPT );
1909 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1912 context_t *ctx = context_new();
1914 /* Set the message callback in the working context */
1915 message_cb_user_data mcud;
1916 mcud._error = & this->SMESH_Algo::_comment;
1917 mcud._progress = & this->SMESH_Algo::_progress;
1919 _hypothesis ? _hypothesis->GetVerbosity() : BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
1920 context_set_message_callback(ctx, message_cb, &mcud);
1922 /* set the interruption callback */
1923 _compute_canceled = false;
1924 context_set_interrupt_callback(ctx, interrupt_cb, this);
1926 /* create the CAD object we will work on. It is associated to the context ctx. */
1927 cad_t *c = cad_new(ctx);
1928 dcad_t *dcad = dcad_new(c);
1930 // To enable multithreading
1931 cad_set_thread_safety(c, 1);
1933 /* Now fill the CAD object with data from your CAD
1934 * environement. This is the most complex part of a successfull
1940 cadsurf_session_t *css = cadsurf_session_new(ctx);
1942 // an object that correctly deletes all cadsurf objects at destruction
1943 BLSURF_Cleaner cleaner( ctx,css,c,dcad );
1945 SetParameters(_hypothesis, css, aShape);
1947 haveQuadraticSubMesh = haveQuadraticSubMesh || (_hypothesis != NULL && _hypothesis->GetQuadraticMesh());
1948 helper.SetIsQuadratic( haveQuadraticSubMesh );
1950 // To remove as soon as quadratic mesh is allowed - BEGIN
1951 // GDD: Viscous layer is not allowed with quadratic mesh
1952 if (_haveViscousLayers && haveQuadraticSubMesh ) {
1953 quadraticSubMeshAndViscousLayer = true;
1954 _haveViscousLayers = !haveQuadraticSubMesh;
1955 _comment += "Warning: Viscous layer is not possible with a quadratic mesh, it is ignored.";
1956 error(COMPERR_WARNING, _comment);
1958 // To remove as soon as quadratic mesh is allowed - END
1960 // needed to prevent the opencascade memory managmement from freeing things
1961 vector<Handle(Geom2d_Curve)> curves;
1962 vector<Handle(Geom_Surface)> surfaces;
1966 FaceId2PythonSmp.clear();
1967 EdgeId2PythonSmp.clear();
1968 VertexId2PythonSmp.clear();
1970 /****************************************************************************************
1972 *****************************************************************************************/
1974 string bad_end = "return";
1976 TopTools_IndexedMapOfShape _map;
1977 TopExp::MapShapes(aShape,TopAbs_VERTEX,_map);
1978 int ienf = _map.Extent();
1980 assert(Py_IsInitialized());
1981 PyGILState_STATE gstate;
1983 string theSizeMapStr;
1985 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1987 TopoDS_Face f = TopoDS::Face(face_iter.Current());
1989 SMESH_subMesh* fSM = aMesh.GetSubMesh( f );
1990 if ( !fSM->IsEmpty() ) continue; // skip already meshed FACE with viscous layers
1992 // make INTERNAL face oriented FORWARD (issue 0020993)
1993 if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
1994 f.Orientation(TopAbs_FORWARD);
1996 iface = fmap.Add(f);
1998 surfaces.push_back(BRep_Tool::Surface(f));
2000 /* create an object representing the face for cadsurf */
2001 /* where face_id is an integer identifying the face.
2002 * surf_function is the function that defines the surface
2003 * (For this face, it will be called by cadsurf with your_face_object_ptr
2004 * as last parameter.
2006 #if OCC_VERSION_MAJOR < 7
2007 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
2009 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back().get());
2012 /* by default a face has no tag (color).
2013 The following call sets it to the same value as the Geom module ID : */
2014 int faceTag = meshDS->ShapeToIndex(f);
2015 faceTag = BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis );
2016 cad_face_set_tag(fce, faceTag);
2018 /* Set face orientation (optional if you want a well oriented output mesh)*/
2019 if(f.Orientation() != TopAbs_FORWARD)
2020 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
2022 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
2024 if (HasSizeMapOnFace /*&& !use_precad*/) //22903: use_precad seems not to interfere
2026 // -----------------
2028 // -----------------
2029 faceKey = FacesWithSizeMap.FindIndex(f);
2032 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end())
2034 theSizeMapStr = FaceId2SizeMap[faceKey];
2035 // check if function ends with "return"
2036 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2038 // Expr To Python function, verification is performed at validation in GUI
2039 gstate = PyGILState_Ensure();
2040 PyObject * obj = NULL;
2041 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2043 PyObject * func = NULL;
2044 func = PyObject_GetAttrString(main_mod, "f");
2045 FaceId2PythonSmp[iface]=func;
2046 FaceId2SizeMap.erase(faceKey);
2047 PyGILState_Release(gstate);
2050 // Specific size map = Attractor
2051 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
2053 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
2054 if (attractor_iter->first == faceKey)
2056 double xyzCoords[3] = {attractor_iter->second[2],
2057 attractor_iter->second[3],
2058 attractor_iter->second[4]};
2060 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
2061 BRepClass_FaceClassifier scl(f,P,1e-7);
2062 scl.Perform(f, P, 1e-7);
2063 TopAbs_State result = scl.State();
2064 if ( result == TopAbs_OUT )
2065 MESSAGE("Point is out of face: node is not created");
2066 if ( result == TopAbs_UNKNOWN )
2067 MESSAGE("Point position on face is unknown: node is not created");
2068 if ( result == TopAbs_ON )
2069 MESSAGE("Point is on border of face: node is not created");
2070 if ( result == TopAbs_IN )
2072 // Point is inside face and not on border
2073 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
2075 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
2076 cad_point_set_tag(point_p, ienf);
2078 FaceId2AttractorCoords.erase(faceKey);
2082 // -----------------
2084 // -----------------
2085 TId2ClsAttractorVec::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
2086 if (clAttractor_iter != FaceId2ClassAttractor.end()){
2087 std::vector< BLSURFPlugin_Attractor* > & attVec = clAttractor_iter->second;
2088 for ( size_t i = 0; i < attVec.size(); ++i )
2089 if ( !attVec[i]->IsMapBuilt() ) {
2090 std::cout<<"Compute " << theNbAttractors-- << "-th attractor" <<std::endl;
2091 attVec[i]->BuildMap();
2093 FaceIndex2ClassAttractor[iface].swap( attVec );
2094 FaceId2ClassAttractor.erase(clAttractor_iter);
2096 } // if (HasSizeMapOnFace && !use_precad)
2098 // ------------------
2099 // Enforced Vertices
2100 // ------------------
2101 faceKey = FacesWithEnforcedVertices.FindIndex(f);
2102 std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
2103 if (evmIt != FaceId2EnforcedVertexCoords.end())
2105 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl = evmIt->second;
2106 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
2107 for (; evlIt != evl.end(); ++evlIt)
2109 double uvCoords[2] = { evlIt->at(0), evlIt->at(1) };
2111 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
2113 BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
2114 xyzCoords.push_back(evlIt->at(2));
2115 xyzCoords.push_back(evlIt->at(3));
2116 xyzCoords.push_back(evlIt->at(4));
2117 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(xyzCoords);
2118 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end() &&
2119 !enfCoordsIt->second.empty() )
2121 // to merge nodes of an INTERNAL vertex belonging to several faces
2122 TopoDS_Vertex v = (*enfCoordsIt->second.begin() )->vertex;
2123 if ( v.IsNull() ) v = (*enfCoordsIt->second.rbegin())->vertex;
2124 if ( !v.IsNull() && meshDS->ShapeToIndex( v ) > 0 )
2126 tag = pmap.Add( v );
2127 SMESH_subMesh* vSM = aMesh.GetSubMesh( v );
2128 vSM->ComputeStateEngine( SMESH_subMesh::COMPUTE );
2129 mergeSubmeshes.insert( vSM->GetSubMeshDS() );
2130 // //if ( tag != pmap.Extent() )
2131 // needMerge = true;
2134 if ( tag == 0 ) tag = ienf;
2135 cad_point_set_tag(point_p, tag);
2137 FaceId2EnforcedVertexCoords.erase(faceKey);
2141 /****************************************************************************************
2143 now create the edges associated to this face
2144 *****************************************************************************************/
2146 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next())
2148 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
2149 int ic = emap.FindIndex(e);
2154 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
2156 if (HasSizeMapOnEdge){
2157 edgeKey = EdgesWithSizeMap.FindIndex(e);
2158 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end())
2160 theSizeMapStr = EdgeId2SizeMap[edgeKey];
2161 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2163 // Expr To Python function, verification is performed at validation in GUI
2164 gstate = PyGILState_Ensure();
2165 PyObject * obj = NULL;
2166 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2168 PyObject * func = NULL;
2169 func = PyObject_GetAttrString(main_mod, "f");
2170 EdgeId2PythonSmp[ic]=func;
2171 EdgeId2SizeMap.erase(edgeKey);
2172 PyGILState_Release(gstate);
2175 /* data of nodes existing on the edge */
2176 StdMeshers_FaceSidePtr nodeData;
2177 SMESH_subMesh* sm = aMesh.GetSubMesh( e );
2178 if ( !sm->IsEmpty() )
2180 // SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
2181 // /*complexFirst=*/false);
2182 // while ( subsmIt->more() )
2183 // edgeSubmeshes.insert( subsmIt->next()->GetSubMeshDS() );
2184 edgeSubmeshes.insert( sm->GetSubMeshDS() );
2186 nodeData.reset( new StdMeshers_FaceSide( f, e, &aMesh, /*isForwrd = */true,
2187 /*ignoreMedium=*/haveQuadraticSubMesh));
2188 if ( nodeData->MissVertexNode() )
2189 return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
2191 const std::vector<UVPtStruct>& nodeDataVec = nodeData->GetUVPtStruct();
2192 if ( !nodeDataVec.empty() )
2194 if ( Abs( nodeDataVec[0].param - tmin ) > Abs( nodeDataVec.back().param - tmin ))
2196 nodeData->Reverse();
2197 nodeData->GetUVPtStruct(); // nodeData recomputes nodeDataVec
2199 // tmin and tmax can change in case of viscous layer on an adjacent edge
2200 tmin = nodeDataVec.front().param;
2201 tmax = nodeDataVec.back().param;
2205 cout << "---------------- Invalid nodeData" << endl;
2210 /* attach the edge to the current cadsurf face */
2211 #if OCC_VERSION_MAJOR < 7
2212 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
2214 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back().get());
2217 /* by default an edge has no tag (color).
2218 The following call sets it to the same value as the edge_id : */
2219 // IMP23368. Do not set tag to an EDGE shared by FACEs of a hyper-patch
2220 bool isInHyperPatch = false;
2222 std::set< int > faceTags, faceIDs;
2223 TopTools_ListIteratorOfListOfShape fIt( e2ffmap.FindFromKey( e ));
2224 for ( ; fIt.More(); fIt.Next() )
2226 int faceTag = meshDS->ShapeToIndex( fIt.Value() );
2227 if ( !faceIDs.insert( faceTag ).second )
2228 continue; // a face encounters twice for a seam edge
2229 int hpTag = BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis );
2230 if ( !faceTags.insert( hpTag ).second )
2232 isInHyperPatch = true;
2237 if ( !isInHyperPatch )
2238 cad_edge_set_tag(edg, ic);
2240 /* by default, an edge does not necessalry appear in the resulting mesh,
2241 unless the following property is set :
2243 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
2245 /* by default an edge is a boundary edge */
2246 if (e.Orientation() == TopAbs_INTERNAL)
2247 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
2249 // pass existing nodes of sub-meshes to MG-CADSurf
2252 const std::vector<UVPtStruct>& nodeDataVec = nodeData->GetUVPtStruct();
2253 const int nbNodes = nodeDataVec.size();
2255 dcad_edge_discretization_t *dedge;
2256 dcad_get_edge_discretization(dcad, edg, &dedge);
2257 dcad_edge_discretization_set_vertex_count( dedge, nbNodes );
2259 // cout << endl << " EDGE " << ic << endl;
2260 // cout << "tmin = "<<tmin << ", tmax = "<< tmax << endl;
2261 for ( int iN = 0; iN < nbNodes; ++iN )
2263 const UVPtStruct& nData = nodeDataVec[ iN ];
2264 double t = nData.param;
2265 real uv[2] = { nData.u, nData.v };
2266 SMESH_TNodeXYZ nXYZ( nData.node );
2267 // cout << "\tt = " << t
2268 // << "\t uv = ( " << uv[0] << ","<< uv[1] << " ) "
2269 // << "\t u = " << nData.param
2270 // << "\t ID = " << nData.node->GetID() << endl;
2271 dcad_edge_discretization_set_vertex_coordinates( dedge, iN+1, t, uv, nXYZ._xyz );
2273 dcad_edge_discretization_set_property(dedge, DISTENE_DCAD_PROPERTY_REQUIRED);
2276 /****************************************************************************************
2278 *****************************************************************************************/
2282 gp_Pnt2d e0 = curves.back()->Value(tmin);
2283 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
2284 Standard_Real d1=0,d2=0;
2287 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
2288 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
2292 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
2295 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
2297 *ip = pmap.FindIndex(v);
2300 // SMESH_subMesh* sm = aMesh.GetSubMesh(v);
2301 // if ( sm->IsMeshComputed() )
2302 // edgeSubmeshes.insert( sm->GetSubMeshDS() );
2305 // std::string aFileName = "fmap_vertex_";
2306 // aFileName.append(val_to_string(*ip));
2307 // aFileName.append(".brep");
2308 // BRepTools::Write(v,aFileName.c_str());
2310 if (HasSizeMapOnVertex){
2311 vertexKey = VerticesWithSizeMap.FindIndex(v);
2312 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
2313 theSizeMapStr = VertexId2SizeMap[vertexKey];
2314 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2316 // Expr To Python function, verification is performed at validation in GUI
2317 gstate = PyGILState_Ensure();
2318 PyObject * obj = NULL;
2319 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2321 PyObject * func = NULL;
2322 func = PyObject_GetAttrString(main_mod, "f");
2323 VertexId2PythonSmp[*ip]=func;
2324 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
2325 PyGILState_Release(gstate);
2330 // should not happen
2331 MESSAGE("An edge does not have 2 extremities.");
2334 // This defines the curves extremity connectivity
2335 cad_edge_set_extremities(edg, ip1, ip2);
2336 /* set the tag (color) to the same value as the extremity id : */
2337 cad_edge_set_extremities_tag(edg, ip1, ip2);
2340 cad_edge_set_extremities(edg, ip2, ip1);
2341 cad_edge_set_extremities_tag(edg, ip2, ip1);
2347 // Clear mesh from already meshed edges if possible else
2348 // remember that merge is needed
2349 TSubMeshSet::iterator smIt = edgeSubmeshes.begin();
2350 for ( ; smIt != edgeSubmeshes.end(); ++smIt ) // loop on already meshed EDGEs
2352 SMESHDS_SubMesh* smDS = *smIt;
2353 if ( !smDS ) continue;
2354 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2357 const SMDS_MeshNode* n = nIt->next();
2358 if ( n->NbInverseElements( SMDSAbs_Face ) > 0 )
2360 needMerge = true; // to correctly sew with viscous mesh
2361 // add existing medium nodes to helper
2362 if ( aMesh.NbEdges( ORDER_QUADRATIC ) > 0 )
2364 SMDS_ElemIteratorPtr edgeIt = smDS->GetElements();
2365 while ( edgeIt->more() )
2366 helper.AddTLinks( static_cast<const SMDS_MeshEdge*>(edgeIt->next()));
2371 if ( allowSubMeshClearing )
2373 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
2374 while ( eIt->more() ) meshDS->RemoveFreeElement( eIt->next(), 0 );
2375 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2376 while ( nIt->more() ) meshDS->RemoveFreeNode( nIt->next(), 0 );
2385 ///////////////////////
2387 ///////////////////////
2389 if (! _preCadFacesIDsPeriodicityVector.empty())
2391 for (std::size_t i=0; i < _preCadFacesIDsPeriodicityVector.size(); i++){
2392 std::vector<int> theFace1_ids = _preCadFacesIDsPeriodicityVector[i].shape1IDs;
2393 std::vector<int> theFace2_ids = _preCadFacesIDsPeriodicityVector[i].shape2IDs;
2394 int* theFace1_ids_c = &theFace1_ids[0];
2395 int* theFace2_ids_c = &theFace2_ids[0];
2396 std::ostringstream o;
2397 o << "_preCadFacesIDsPeriodicityVector[" << i << "] = [";
2398 for (std::size_t j=0; j < theFace1_ids.size(); j++)
2399 o << theFace1_ids[j] << ", ";
2401 for (std::size_t j=0; j < theFace2_ids.size(); j++)
2402 o << theFace2_ids[j] << ", ";
2404 // if ( _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
2405 // cout << o.str() << endl;
2406 if (_preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords.empty())
2408 // If no source points, call periodicity without transformation function
2409 meshgems_cad_periodicity_transformation_t periodicity_transformation = NULL;
2410 status = cad_add_face_multiple_periodicity_with_transformation_function(c, theFace1_ids_c, theFace1_ids.size(),
2411 theFace2_ids_c, theFace2_ids.size(), periodicity_transformation, NULL);
2412 if(status != STATUS_OK)
2413 cout << "cad_add_face_multiple_periodicity_with_transformation_function failed with error code " << status << "\n";
2417 // get the transformation vertices
2418 double* theSourceVerticesCoords_c = &_preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords[0];
2419 double* theTargetVerticesCoords_c = &_preCadFacesIDsPeriodicityVector[i].theTargetVerticesCoords[0];
2420 int nbSourceVertices = _preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords.size()/3;
2421 int nbTargetVertices = _preCadFacesIDsPeriodicityVector[i].theTargetVerticesCoords.size()/3;
2423 status = cad_add_face_multiple_periodicity_with_transformation_function_by_points(c, theFace1_ids_c, theFace1_ids.size(),
2424 theFace2_ids_c, theFace2_ids.size(), theSourceVerticesCoords_c, nbSourceVertices, theTargetVerticesCoords_c, nbTargetVertices);
2425 if(status != STATUS_OK)
2426 cout << "cad_add_face_multiple_periodicity_with_transformation_function_by_points failed with error code " << status << "\n";
2431 if (! _preCadEdgesIDsPeriodicityVector.empty())
2433 for (std::size_t i=0; i < _preCadEdgesIDsPeriodicityVector.size(); i++){
2434 std::vector<int> theEdge1_ids = _preCadEdgesIDsPeriodicityVector[i].shape1IDs;
2435 std::vector<int> theEdge2_ids = _preCadEdgesIDsPeriodicityVector[i].shape2IDs;
2436 // Use the address of the first element of the vector to initialize the array
2437 int* theEdge1_ids_c = &theEdge1_ids[0];
2438 int* theEdge2_ids_c = &theEdge2_ids[0];
2440 std::ostringstream o;
2441 o << "_preCadEdgesIDsPeriodicityVector[" << i << "] = [";
2442 for (std::size_t j=0; j < theEdge1_ids.size(); j++)
2443 o << theEdge1_ids[j] << ", ";
2445 for (std::size_t j=0; j < theEdge2_ids.size(); j++)
2446 o << theEdge2_ids[j] << ", ";
2448 // if ( _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
2449 // cout << o.str() << endl;
2451 if (_preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords.empty())
2453 // If no source points, call periodicity without transformation function
2454 meshgems_cad_periodicity_transformation_t periodicity_transformation = NULL;
2455 status = cad_add_edge_multiple_periodicity_with_transformation_function(c, theEdge1_ids_c, theEdge1_ids.size(),
2456 theEdge2_ids_c, theEdge2_ids.size(), periodicity_transformation, NULL);
2457 if(status != STATUS_OK)
2458 cout << "cad_add_edge_multiple_periodicity_with_transformation_function failed with error code " << status << "\n";
2462 // get the transformation vertices
2463 double* theSourceVerticesCoords_c = &_preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords[0];
2464 double* theTargetVerticesCoords_c = &_preCadEdgesIDsPeriodicityVector[i].theTargetVerticesCoords[0];
2465 int nbSourceVertices = _preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords.size()/3;
2466 int nbTargetVertices = _preCadEdgesIDsPeriodicityVector[i].theTargetVerticesCoords.size()/3;
2468 status = cad_add_edge_multiple_periodicity_with_transformation_function_by_points(c, theEdge1_ids_c, theEdge1_ids.size(),
2469 theEdge2_ids_c, theEdge2_ids.size(), theSourceVerticesCoords_c, nbSourceVertices, theTargetVerticesCoords_c, nbTargetVertices);
2470 if(status != STATUS_OK)
2471 cout << "cad_add_edge_multiple_periodicity_with_transformation_function_by_points failed with error code " << status << "\n";
2477 // TODO: be able to use a mesh in input.
2478 // See imsh usage in Products/templates/mg-cadsurf_template_common.cpp
2479 // => cadsurf_set_mesh
2481 // Use the original dcad
2482 cadsurf_set_dcad(css, dcad);
2484 // Use the original cad
2485 cadsurf_set_cad(css, c);
2487 std::cout << std::endl;
2488 std::cout << "Beginning of Surface Mesh generation" << std::endl;
2489 std::cout << std::endl;
2494 status = cadsurf_compute_mesh(css);
2497 catch ( std::exception& exc ) {
2498 _comment += exc.what();
2500 catch (Standard_Failure& ex) {
2501 _comment += ex.DynamicType()->Name();
2502 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
2504 _comment += ex.GetMessageString();
2508 if ( _comment.empty() )
2509 _comment = "Exception in cadsurf_compute_mesh()";
2512 std::cout << std::endl;
2513 std::cout << "End of Surface Mesh generation" << std::endl;
2514 std::cout << std::endl;
2517 cadsurf_get_mesh(css, &msh);
2519 /* release the mesh object */
2520 cadsurf_regain_mesh(css, msh);
2521 return error(_comment);
2524 std::string GMFFileName = BLSURFPlugin_Hypothesis::GetDefaultGMFFile();
2526 GMFFileName = _hypothesis->GetGMFFile();
2527 if (GMFFileName != "") {
2528 bool asciiFound = (GMFFileName.find(".mesh", GMFFileName.length()-5) != std::string::npos);
2529 bool binaryFound = (GMFFileName.find(".meshb",GMFFileName.length()-6) != std::string::npos);
2530 if (!asciiFound && !binaryFound)
2531 GMFFileName.append(".mesh");
2532 mesh_write_mesh(msh, GMFFileName.c_str());
2535 /* retrieve mesh data (see meshgems/mesh.h) */
2536 integer nv, ne, nt, nq, vtx[4], tag, nb_tag;
2537 integer *evedg, *evtri, *evquad, *tags_buff, type;
2540 mesh_get_vertex_count(msh, &nv);
2541 mesh_get_edge_count(msh, &ne);
2542 mesh_get_triangle_count(msh, &nt);
2543 mesh_get_quadrangle_count(msh, &nq);
2545 evedg = (integer *)mesh_calloc_generic_buffer(msh);
2546 evtri = (integer *)mesh_calloc_generic_buffer(msh);
2547 evquad = (integer *)mesh_calloc_generic_buffer(msh);
2548 tags_buff = (integer*)mesh_calloc_generic_buffer(msh);
2550 std::vector<const SMDS_MeshNode*> nodes(nv+1);
2551 std::vector<bool> tags(nv+1);
2553 /* enumerated vertices */
2554 for(int iv=1;iv<=nv;iv++) {
2555 mesh_get_vertex_coordinates(msh, iv, xyz);
2556 mesh_get_vertex_tag(msh, iv, &tag);
2557 // Issue 0020656. Use vertex coordinates
2559 if ( tag > 0 && tag <= pmap.Extent() ) {
2560 TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
2561 double tol = BRep_Tool::Tolerance( v );
2562 gp_Pnt p = BRep_Tool::Pnt( v );
2563 if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
2564 xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
2566 tag = 0; // enforced or attracted vertex
2567 nodes[iv] = SMESH_Algo::VertexNode( v, meshDS );
2570 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
2572 // Create group of enforced vertices if requested
2573 BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
2575 projVertex.push_back((double)xyz[0]);
2576 projVertex.push_back((double)xyz[1]);
2577 projVertex.push_back((double)xyz[2]);
2578 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
2579 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end())
2581 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
2582 BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
2583 for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
2584 currentEnfVertex = (*enfListIt);
2585 if (currentEnfVertex->grpName != "") {
2586 bool groupDone = false;
2587 SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
2588 while (grIt->more()) {
2589 SMESH_Group * group = grIt->next();
2590 if ( !group ) continue;
2591 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
2592 if ( !groupDS ) continue;
2593 if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
2594 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
2595 aGroupDS->SMDSGroup().Add(nodes[iv]);
2596 // How can I inform the hypothesis ?
2597 // _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
2605 SMESH_Group* aGroup = aMesh.AddGroup(SMDSAbs_Node, currentEnfVertex->grpName.c_str(), groupId);
2606 aGroup->SetName( currentEnfVertex->grpName.c_str() );
2607 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
2608 aGroupDS->SMDSGroup().Add(nodes[iv]);
2612 throw SALOME_Exception(LOCALIZED("An enforced vertex node was not added to a group"));
2615 MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
2619 // internal points are tagged to zero
2620 if(tag > 0 && tag <= pmap.Extent() ){
2621 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
2628 /* enumerate edges */
2629 for(int it=1;it<=ne;it++) {
2631 mesh_get_edge_vertices(msh, it, vtx);
2632 mesh_get_edge_extra_vertices(msh, it, &type, evedg);
2633 mesh_get_edge_tag(msh, it, &tag);
2635 // If PreCAD performed some cleaning operations (remove tiny edges,
2636 // merge edges ...) an output tag can indeed represent several original tags.
2637 // Get the initial tags corresponding to the output tag and redefine the tag as
2638 // the last of the two initial tags (else the output tag is out of emap and hasn't any meaning)
2639 mesh_get_composite_tag_definition(msh, tag, &nb_tag, tags_buff);
2641 tag=tags_buff[nb_tag-1];
2642 if ( tag < 1 || tag > emap.Extent() )
2644 std::cerr << "MG-CADSurf BUG:::: Edge tag " << tag
2645 << " does not point to a CAD edge (nb edges " << emap.Extent() << ")" << std::endl;
2649 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
2650 tags[vtx[0]] = false;
2653 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
2654 tags[vtx[1]] = false;
2656 if (type == MESHGEMS_MESH_ELEMENT_TYPE_EDGE3) {
2658 if (tags[evedg[0]]) {
2659 Set_NodeOnEdge(meshDS, nodes[evedg[0]], emap(tag));
2660 tags[evedg[0]] = false;
2662 edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]], nodes[evedg[0]]);
2665 edg = helper.AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
2667 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
2670 /* enumerate triangles */
2671 for(int it=1;it<=nt;it++) {
2673 mesh_get_triangle_vertices(msh, it, vtx);
2674 mesh_get_triangle_extra_vertices(msh, it, &type, evtri);
2675 mesh_get_triangle_tag(msh, it, &tag);
2677 meshDS->SetNodeOnFace(nodes[vtx[0]], tag);
2678 tags[vtx[0]] = false;
2681 meshDS->SetNodeOnFace(nodes[vtx[1]], tag);
2682 tags[vtx[1]] = false;
2685 meshDS->SetNodeOnFace(nodes[vtx[2]], tag);
2686 tags[vtx[2]] = false;
2688 if (type == MESHGEMS_MESH_ELEMENT_TYPE_TRIA6) {
2689 // QUADRATIC TRIANGLE
2690 if (tags[evtri[0]]) {
2691 meshDS->SetNodeOnFace(nodes[evtri[0]], tag);
2692 tags[evtri[0]] = false;
2694 if (tags[evtri[1]]) {
2695 meshDS->SetNodeOnFace(nodes[evtri[1]], tag);
2696 tags[evtri[1]] = false;
2698 if (tags[evtri[2]]) {
2699 meshDS->SetNodeOnFace(nodes[evtri[2]], tag);
2700 tags[evtri[2]] = false;
2702 tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]],
2703 nodes[evtri[0]], nodes[evtri[1]], nodes[evtri[2]]);
2706 tri = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
2708 meshDS->SetMeshElementOnShape(tri, tag);
2711 /* enumerate quadrangles */
2712 for(int it=1;it<=nq;it++) {
2713 SMDS_MeshFace* quad;
2714 mesh_get_quadrangle_vertices(msh, it, vtx);
2715 mesh_get_quadrangle_extra_vertices(msh, it, &type, evquad);
2716 mesh_get_quadrangle_tag(msh, it, &tag);
2718 meshDS->SetNodeOnFace(nodes[vtx[0]], tag);
2719 tags[vtx[0]] = false;
2722 meshDS->SetNodeOnFace(nodes[vtx[1]], tag);
2723 tags[vtx[1]] = false;
2726 meshDS->SetNodeOnFace(nodes[vtx[2]], tag);
2727 tags[vtx[2]] = false;
2730 meshDS->SetNodeOnFace(nodes[vtx[3]], tag);
2731 tags[vtx[3]] = false;
2733 if (type == MESHGEMS_MESH_ELEMENT_TYPE_QUAD9) {
2734 // QUADRATIC QUADRANGLE
2735 std::cout << "This is a quadratic quadrangle" << std::endl;
2736 if (tags[evquad[0]]) {
2737 meshDS->SetNodeOnFace(nodes[evquad[0]], tag);
2738 tags[evquad[0]] = false;
2740 if (tags[evquad[1]]) {
2741 meshDS->SetNodeOnFace(nodes[evquad[1]], tag);
2742 tags[evquad[1]] = false;
2744 if (tags[evquad[2]]) {
2745 meshDS->SetNodeOnFace(nodes[evquad[2]], tag);
2746 tags[evquad[2]] = false;
2748 if (tags[evquad[3]]) {
2749 meshDS->SetNodeOnFace(nodes[evquad[3]], tag);
2750 tags[evquad[3]] = false;
2752 if (tags[evquad[4]]) {
2753 meshDS->SetNodeOnFace(nodes[evquad[4]], tag);
2754 tags[evquad[4]] = false;
2756 quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]],
2757 nodes[evquad[0]], nodes[evquad[1]], nodes[evquad[2]], nodes[evquad[3]],
2761 quad = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
2763 meshDS->SetMeshElementOnShape(quad, tag);
2766 /* release the mesh object, the rest is released by cleaner */
2767 cadsurf_regain_mesh(css, msh);
2770 // Remove free nodes that can appear e.g. if "remove tiny edges"(IPAL53235)
2771 for(int iv=1;iv<=nv;iv++)
2772 if ( nodes[iv] && nodes[iv]->NbInverseElements() == 0 )
2773 meshDS->RemoveFreeNode( nodes[iv], 0, /*fromGroups=*/false );
2776 if ( needMerge ) // sew mesh computed by MG-CADSurf with pre-existing mesh
2778 SMESH_MeshEditor editor( &aMesh );
2779 SMESH_MeshEditor::TListOfListOfNodes nodeGroupsToMerge;
2780 TIDSortedElemSet segementsOnEdge;
2781 TSubMeshSet::iterator smIt;
2782 SMESHDS_SubMesh* smDS;
2784 // merge nodes on EDGE's with ones computed by MG-CADSurf
2785 for ( smIt = mergeSubmeshes.begin(); smIt != mergeSubmeshes.end(); ++smIt )
2787 if (! (smDS = *smIt) ) continue;
2788 getNodeGroupsToMerge( smDS, meshDS->IndexToShape((*smIt)->GetID()), nodeGroupsToMerge );
2790 SMDS_ElemIteratorPtr segIt = smDS->GetElements();
2791 while ( segIt->more() )
2792 segementsOnEdge.insert( segIt->next() );
2795 editor.MergeNodes( nodeGroupsToMerge );
2798 SMESH_MeshEditor::TListOfListOfElementsID equalSegments;
2799 editor.FindEqualElements( segementsOnEdge, equalSegments );
2800 editor.MergeElements( equalSegments );
2802 // remove excess segments created on the boundary of viscous layers
2803 const SMDS_TypeOfPosition onFace = SMDS_TOP_FACE;
2804 for ( int i = 1; i <= emap.Extent(); ++i )
2806 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( emap( i )))
2808 SMDS_ElemIteratorPtr segIt = smDS->GetElements();
2809 while ( segIt->more() )
2811 const SMDS_MeshElement* seg = segIt->next();
2812 if ( seg->GetNode(0)->GetPosition()->GetTypeOfPosition() == onFace ||
2813 seg->GetNode(1)->GetPosition()->GetTypeOfPosition() == onFace )
2814 meshDS->RemoveFreeElement( seg, smDS );
2821 // SetIsAlwaysComputed( true ) to sub-meshes of EDGEs w/o mesh
2822 for (int i = 1; i <= emap.Extent(); i++)
2823 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
2824 sm->SetIsAlwaysComputed( true );
2825 for (int i = 1; i <= pmap.Extent(); i++)
2826 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( pmap( i )))
2827 if ( !sm->IsMeshComputed() )
2828 sm->SetIsAlwaysComputed( true );
2830 // Set error to FACE's w/o elements
2831 SMESH_ComputeErrorName err = COMPERR_ALGO_FAILED;
2832 if ( _comment.empty() && status == STATUS_OK )
2834 err = COMPERR_WARNING;
2835 _comment = "No mesh elements assigned to a face";
2837 bool badFaceFound = false;
2838 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
2840 TopoDS_Face f = TopoDS::Face(face_iter.Current());
2841 SMESH_subMesh* sm = aMesh.GetSubMesh( f );
2842 if ( !sm->GetSubMeshDS() || sm->GetSubMeshDS()->NbElements() == 0 )
2844 int faceTag = sm->GetId();
2845 if ( faceTag != BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis ))
2847 // triangles are assigned to the first face of hyper-patch
2848 sm->SetIsAlwaysComputed( true );
2852 sm->GetComputeError().reset( new SMESH_ComputeError( err, _comment, this ));
2853 badFaceFound = true;
2857 if ( err == COMPERR_WARNING )
2861 if ( status != STATUS_OK && !badFaceFound ) {
2865 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
2867 if ( oldFEFlags > 0 )
2868 feenableexcept( oldFEFlags );
2869 feclearexcept( FE_ALL_EXCEPT );
2873 std::cout << "FacesWithSizeMap" << std::endl;
2874 FacesWithSizeMap.Statistics(std::cout);
2875 std::cout << "EdgesWithSizeMap" << std::endl;
2876 EdgesWithSizeMap.Statistics(std::cout);
2877 std::cout << "VerticesWithSizeMap" << std::endl;
2878 VerticesWithSizeMap.Statistics(std::cout);
2879 std::cout << "FacesWithEnforcedVertices" << std::endl;
2880 FacesWithEnforcedVertices.Statistics(std::cout);
2883 return ( status == STATUS_OK && !quadraticSubMeshAndViscousLayer );
2886 //================================================================================
2888 * \brief Compute a mesh basing on discrete CAD description
2890 //================================================================================
2892 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
2894 if ( aMesh.NbFaces() == 0 )
2895 return error( COMPERR_BAD_INPUT_MESH, "2D elements are missing" );
2897 context_t *ctx = context_new();
2898 if (!ctx) return error("Pb in context_new()");
2900 BLSURF_Cleaner cleaner( ctx );
2902 message_cb_user_data mcud;
2903 mcud._error = & this->SMESH_Algo::_comment;
2904 mcud._progress = & this->SMESH_Algo::_progress;
2906 _hypothesis ? _hypothesis->GetVerbosity() : BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
2907 meshgems_status_t ret = context_set_message_callback(ctx, message_cb, &mcud);
2908 if (ret != STATUS_OK) return error("Pb. in context_set_message_callback() ");
2910 cadsurf_session_t * css = cadsurf_session_new(ctx);
2911 if(!css) return error( "Pb. in cadsurf_session_new() " );
2915 // Fill an input mesh
2917 mesh_t * msh = meshgems_mesh_new_in_memory( ctx );
2918 if ( !msh ) return error("Pb. in meshgems_mesh_new_in_memory()");
2920 // mark nodes used by 2D elements
2921 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
2922 SMDS_NodeIteratorPtr nodeIt = meshDS->nodesIterator();
2923 while ( nodeIt->more() )
2925 const SMDS_MeshNode* n = nodeIt->next();
2926 n->setIsMarked( n->NbInverseElements( SMDSAbs_Face ));
2928 meshgems_mesh_set_vertex_count( msh, meshDS->NbNodes() );
2930 // set node coordinates
2931 if ( meshDS->NbNodes() != meshDS->MaxNodeID() )
2933 meshDS->compactMesh();
2935 SMESH_TNodeXYZ nXYZ;
2936 nodeIt = meshDS->nodesIterator();
2938 for ( i = 1; nodeIt->more(); ++i )
2940 nXYZ.Set( nodeIt->next() );
2941 meshgems_mesh_set_vertex_coordinates( msh, i, nXYZ._xyz );
2944 // set nodes of faces
2945 meshgems_mesh_set_triangle_count ( msh, meshDS->GetMeshInfo().NbTriangles() );
2946 meshgems_mesh_set_quadrangle_count( msh, meshDS->GetMeshInfo().NbQuadrangles() );
2947 meshgems_integer nodeIDs[4];
2948 meshgems_integer iT = 1, iQ = 1;
2949 SMDS_FaceIteratorPtr faceIt = meshDS->facesIterator();
2950 while ( faceIt->more() )
2952 const SMDS_MeshElement* face = faceIt->next();
2953 meshgems_integer nbNodes = face->NbCornerNodes();
2954 if ( nbNodes > 4 || face->IsPoly() ) continue;
2956 for ( i = 0; i < nbNodes; ++i )
2957 nodeIDs[i] = face->GetNode( i )->GetID();
2959 meshgems_mesh_set_triangle_vertices ( msh, iT++, nodeIDs );
2961 meshgems_mesh_set_quadrangle_vertices( msh, iQ++, nodeIDs );
2964 ret = cadsurf_set_mesh(css, msh);
2965 if ( ret != STATUS_OK ) return error("Pb in cadsurf_set_mesh()");
2970 SetParameters(_hypothesis, css, aMesh.GetShapeToMesh() );
2972 ret = cadsurf_compute_mesh(css);
2973 if ( ret != STATUS_OK ) return false;
2976 cadsurf_get_mesh(css, &omsh);
2977 if ( !omsh ) return error( "Pb. in cadsurf_get_mesh()" );
2980 // Update SALOME mesh
2982 // remove quadrangles and triangles
2983 for ( faceIt = meshDS->facesIterator(); faceIt->more(); )
2985 const SMDS_MeshElement* face = faceIt->next();
2986 if ( !face->IsPoly() )
2987 meshDS->RemoveFreeElement( face, /*sm=*/0, /*fromGroups=*/true );
2989 // remove edges that bound the just removed faces
2990 for ( SMDS_EdgeIteratorPtr edgeIt = meshDS->edgesIterator(); edgeIt->more(); )
2992 const SMDS_MeshElement* edge = edgeIt->next();
2993 const SMDS_MeshNode* n0 = edge->GetNode(0);
2994 const SMDS_MeshNode* n1 = edge->GetNode(1);
2995 if ( n0->isMarked() &&
2997 n0->NbInverseElements( SMDSAbs_Volume ) == 0 &&
2998 n1->NbInverseElements( SMDSAbs_Volume ) == 0 )
2999 meshDS->RemoveFreeElement( edge, /*sm=*/0, /*fromGroups=*/true );
3001 // remove nodes that just became free
3002 for ( nodeIt = meshDS->nodesIterator(); nodeIt->more(); )
3004 const SMDS_MeshNode* n = nodeIt->next();
3005 if ( n->isMarked() && n->NbInverseElements() == 0 )
3006 meshDS->RemoveFreeNode( n, /*sm=*/0, /*fromGroups=*/true );
3010 meshgems_integer nbvtx = 0, nodeID;
3011 meshgems_mesh_get_vertex_count( omsh, &nbvtx );
3012 meshgems_real xyz[3];
3013 for ( i = 1; i <= nbvtx; ++i )
3015 meshgems_mesh_get_vertex_coordinates( omsh, i, xyz );
3016 SMDS_MeshNode* n = meshDS->AddNode( xyz[0], xyz[1], xyz[2] );
3017 nodeID = n->GetID();
3018 meshgems_mesh_set_vertex_tag( omsh, i, &nodeID ); // save mapping of IDs in MG and SALOME meshes
3022 meshgems_integer nbtri = 0;
3023 meshgems_mesh_get_triangle_count( omsh, &nbtri );
3024 const SMDS_MeshNode* nodes[3];
3025 for ( i = 1; i <= nbtri; ++i )
3027 meshgems_mesh_get_triangle_vertices( omsh, i, nodeIDs );
3028 for ( int j = 0; j < 3; ++j )
3030 meshgems_mesh_get_vertex_tag( omsh, nodeIDs[j], &nodeID );
3031 nodes[j] = meshDS->FindNode( nodeID );
3033 meshDS->AddFace( nodes[0], nodes[1], nodes[2] );
3036 cadsurf_regain_mesh(css, omsh);
3038 // as we don't assign the new triangles to a shape (the pseudo-shape),
3039 // we mark the shape as always computed to avoid the error messages
3040 // that no elements assigned to the shape
3041 aMesh.GetSubMesh( aHelper->GetSubShape() )->SetIsAlwaysComputed( true );
3046 //================================================================================
3048 * \brief Terminates computation
3050 //================================================================================
3052 void BLSURFPlugin_BLSURF::CancelCompute()
3054 _compute_canceled = true;
3057 //=============================================================================
3061 //=============================================================================
3063 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS,
3064 const SMDS_MeshNode* node,
3065 const TopoDS_Shape& ed)
3067 const TopoDS_Edge edge = TopoDS::Edge(ed);
3069 gp_Pnt pnt(node->X(), node->Y(), node->Z());
3071 Standard_Real p0 = 0.0;
3072 Standard_Real p1 = 1.0;
3073 TopLoc_Location loc;
3074 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
3075 if ( curve.IsNull() )
3077 // issue 22499. Node at a sphere apex
3078 meshDS->SetNodeOnEdge(node, edge, p0);
3082 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
3083 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
3086 if ( proj.NbPoints() > 0 )
3088 pa = (double)proj.LowerDistanceParameter();
3089 // Issue 0020656. Move node if it is too far from edge
3090 gp_Pnt curve_pnt = curve->Value( pa );
3091 double dist2 = pnt.SquareDistance( curve_pnt );
3092 double tol = BRep_Tool::Tolerance( edge );
3093 if ( 1e-14 < dist2 && dist2 <= 1000*tol ) // large enough and within tolerance
3095 curve_pnt.Transform( loc );
3096 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
3100 meshDS->SetNodeOnEdge(node, edge, pa);
3103 /* Curve definition function See cad_curv_t in file meshgems/cad.h for
3105 * NOTE : if when your CAD systems evaluates second
3106 * order derivatives it also computes first order derivatives and
3107 * function evaluation, you can optimize this example by making only
3108 * one CAD call and filling the necessary uv, dt, dtt arrays.
3110 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
3112 /* t is given. It contains the t (time) 1D parametric coordintaes
3113 of the point PreCAD/MG-CADSurf is querying on the curve */
3115 /* user_data identifies the edge PreCAD/MG-CADSurf is querying
3116 * (see cad_edge_new later in this example) */
3117 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
3120 /* MG-CADSurf is querying the function evaluation */
3123 uv[0]=P.X(); uv[1]=P.Y();
3127 /* query for the first order derivatives */
3130 dt[0]=V1.X(); dt[1]=V1.Y();
3134 /* query for the second order derivatives */
3137 dtt[0]=V2.X(); dtt[1]=V2.Y();
3143 /* Surface definition function.
3144 * See cad_surf_t in file meshgems/cad.h for more information.
3145 * NOTE : if when your CAD systems evaluates second order derivatives it also
3146 * computes first order derivatives and function evaluation, you can optimize
3147 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
3150 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
3151 real *duu, real *duv, real *dvv, void *user_data)
3153 /* uv[2] is given. It contains the u,v coordinates of the point
3154 * PreCAD/MG-CADSurf is querying on the surface */
3156 /* user_data identifies the face PreCAD/MG-CADSurf is querying (see
3157 * cad_face_new later in this example)*/
3158 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
3162 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
3163 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
3170 geometry->D1(uv[0],uv[1],P,D1U,D1V);
3171 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
3172 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
3175 if(duu && duv && dvv){
3179 gp_Vec D2U,D2V,D2UV;
3181 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
3182 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
3183 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
3184 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
3191 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
3193 TId2ClsAttractorVec::iterator f2attVec;
3194 if (FaceId2PythonSmp.count(face_id) != 0){
3195 assert(Py_IsInitialized());
3196 PyGILState_STATE gstate;
3197 gstate = PyGILState_Ensure();
3198 PyObject* pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
3200 if ( pyresult != NULL) {
3201 result = PyFloat_AsDouble(pyresult);
3202 Py_DECREF(pyresult);
3207 string err_description="";
3208 PyObject* new_stderr = newPyStdOut(err_description);
3209 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3210 Py_INCREF(old_stderr);
3211 PySys_SetObject((char*)"stderr", new_stderr);
3213 PySys_SetObject((char*)"stderr", old_stderr);
3214 Py_DECREF(new_stderr);
3215 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
3216 result = *((real*)user_data);
3219 PyGILState_Release(gstate);
3221 else if (( f2attVec = FaceIndex2ClassAttractor.find(face_id)) != FaceIndex2ClassAttractor.end() && !f2attVec->second.empty())
3225 std::vector< BLSURFPlugin_Attractor* > & attVec = f2attVec->second;
3226 for ( size_t i = 0; i < attVec.size(); ++i )
3228 //result += attVec[i]->GetSize(uv[0],uv[1]);
3229 result = Min( result, attVec[i]->GetSize(uv[0],uv[1]));
3231 //*size = result / attVec.size(); // mean of sizes defined by all attractors
3235 *size = *((real*)user_data);
3237 // std::cout << "Size_on_surface sur la face " << face_id << " donne une size de: " << *size << std::endl;
3241 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
3243 if (EdgeId2PythonSmp.count(edge_id) != 0){
3244 assert(Py_IsInitialized());
3245 PyGILState_STATE gstate;
3246 gstate = PyGILState_Ensure();
3247 PyObject* pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
3249 if ( pyresult != NULL) {
3250 result = PyFloat_AsDouble(pyresult);
3251 Py_DECREF(pyresult);
3256 string err_description="";
3257 PyObject* new_stderr = newPyStdOut(err_description);
3258 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3259 Py_INCREF(old_stderr);
3260 PySys_SetObject((char*)"stderr", new_stderr);
3262 PySys_SetObject((char*)"stderr", old_stderr);
3263 Py_DECREF(new_stderr);
3264 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
3265 result = *((real*)user_data);
3268 PyGILState_Release(gstate);
3271 *size = *((real*)user_data);
3276 status_t size_on_vertex(integer point_id, real *size, void *user_data)
3278 if (VertexId2PythonSmp.count(point_id) != 0){
3279 assert(Py_IsInitialized());
3280 PyGILState_STATE gstate;
3281 gstate = PyGILState_Ensure();
3282 PyObject* pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
3284 if ( pyresult != NULL) {
3285 result = PyFloat_AsDouble(pyresult);
3286 Py_DECREF(pyresult);
3291 string err_description="";
3292 PyObject* new_stderr = newPyStdOut(err_description);
3293 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3294 Py_INCREF(old_stderr);
3295 PySys_SetObject((char*)"stderr", new_stderr);
3297 PySys_SetObject((char*)"stderr", old_stderr);
3298 Py_DECREF(new_stderr);
3299 MESSAGE("Can't evaluate f()" << " error is " << err_description);
3300 result = *((real*)user_data);
3303 PyGILState_Release(gstate);
3306 *size = *((real*)user_data);
3312 * The following function will be called for PreCAD/MG-CADSurf message
3313 * printing. See context_set_message_callback (later in this
3314 * template) for how to set user_data.
3316 status_t message_cb(message_t *msg, void *user_data)
3318 integer errnumber = 0;
3320 message_get_number(msg, &errnumber);
3321 message_get_description(msg, &desc);
3323 message_cb_user_data * mcud = (message_cb_user_data*)user_data;
3324 // Get all the error message and some warning messages related to license and periodicity
3325 if ( errnumber < 0 ||
3326 err.find("license" ) != string::npos ||
3327 err.find("periodicity") != string::npos )
3329 // remove ^A from the tail
3330 int len = strlen( desc );
3331 while (len > 0 && desc[len-1] != '\n')
3333 mcud->_error->append( desc, len );
3336 if ( errnumber == 3009001 )
3337 * mcud->_progress = atof( desc + 11 ) / 100.;
3338 if ( mcud->_verbosity > 0 )
3339 std::cout << desc << std::endl;
3344 /* This is the interrupt callback. PreCAD/MG-CADSurf will call this
3345 * function regularily. See the file meshgems/interrupt.h
3347 status_t interrupt_cb(integer *interrupt_status, void *user_data)
3349 integer you_want_to_continue = 1;
3350 BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
3351 you_want_to_continue = !tmp->computeCanceled();
3353 if(you_want_to_continue)
3355 *interrupt_status = INTERRUPT_CONTINUE;
3358 else /* you want to stop MG-CADSurf */
3360 *interrupt_status = INTERRUPT_STOP;
3361 return STATUS_ERROR;
3365 //=============================================================================
3369 //=============================================================================
3370 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
3371 const TopoDS_Shape& aShape,
3372 MapShapeNbElems& aResMap)
3374 double diagonal = aMesh.GetShapeDiagonalSize();
3375 double bbSegmentation = _gen->GetBoundaryBoxSegmentation();
3376 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
3377 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize(diagonal, bbSegmentation);
3378 bool _phySizeRel = BLSURFPlugin_Hypothesis::GetDefaultPhySizeRel();
3379 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
3380 double _angleMesh = BLSURFPlugin_Hypothesis::GetDefaultAngleMesh();
3381 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
3383 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
3384 _phySizeRel = _hypothesis->IsPhySizeRel();
3385 if ( _hypothesis->GetPhySize() > 0)
3386 _phySize = _phySizeRel ? diagonal*_hypothesis->GetPhySize() : _hypothesis->GetPhySize();
3387 //_geometricMesh = (int) hyp->GetGeometricMesh();
3388 if (_hypothesis->GetAngleMesh() > 0)
3389 _angleMesh = _hypothesis->GetAngleMesh();
3390 _quadAllowed = _hypothesis->GetQuadAllowed();
3392 //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
3393 // GetDefaultPhySize() sometimes leads to computation failure
3394 _phySize = aMesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
3397 bool IsQuadratic = _quadraticMesh;
3402 TopTools_DataMapOfShapeInteger EdgesMap;
3403 double fullLen = 0.0;
3404 double fullNbSeg = 0;
3405 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
3406 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3407 if( EdgesMap.IsBound(E) )
3409 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
3410 double aLen = SMESH_Algo::EdgeLength(E);
3413 if(_physicalMesh==1) {
3414 nb1d = (int)( aLen/_phySize + 1 );
3419 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
3420 double fullAng = 0.0;
3421 double dp = (l-f)/200;
3426 for(int j=2; j<=200; j++) {
3429 fullAng += fabs(V1.Angle(V2));
3433 nb1d = (int)( fullAng/_angleMesh + 1 );
3436 std::vector<int> aVec(SMDSEntity_Last);
3437 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
3438 if( IsQuadratic > 0 ) {
3439 aVec[SMDSEntity_Node] = 2*nb1d - 1;
3440 aVec[SMDSEntity_Quad_Edge] = nb1d;
3443 aVec[SMDSEntity_Node] = nb1d - 1;
3444 aVec[SMDSEntity_Edge] = nb1d;
3446 aResMap.insert(std::make_pair(sm,aVec));
3447 EdgesMap.Bind(E,nb1d);
3449 double ELen = fullLen/fullNbSeg;
3453 // try to evaluate as in MEFISTO
3454 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
3455 TopoDS_Face F = TopoDS::Face( exp.Current() );
3456 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
3458 BRepGProp::SurfaceProperties(F,G);
3459 double anArea = G.Mass();
3461 std::vector<int> nb1dVec;
3462 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
3463 int nbSeg = EdgesMap.Find(exp1.Current());
3465 nb1dVec.push_back( nbSeg );
3468 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
3469 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
3472 if ( nb1dVec.size() == 4 ) // quadrangle geom face
3474 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
3476 nbNodes = (n1 + 1) * (n2 + 1);
3481 nbTria = nbQuad = nbTria / 3 + 1;
3484 std::vector<int> aVec(SMDSEntity_Last,0);
3486 int nb1d_in = (nbTria*3 - nb1d) / 2;
3487 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3488 aVec[SMDSEntity_Quad_Triangle] = nbTria;
3489 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
3492 aVec[SMDSEntity_Node] = nbNodes;
3493 aVec[SMDSEntity_Triangle] = nbTria;
3494 aVec[SMDSEntity_Quadrangle] = nbQuad;
3496 aResMap.insert(std::make_pair(sm,aVec));
3503 BRepGProp::VolumeProperties(aShape,G);
3504 double aVolume = G.Mass();
3505 double tetrVol = 0.1179*ELen*ELen*ELen;
3506 int nbVols = int(aVolume/tetrVol);
3507 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3508 std::vector<int> aVec(SMDSEntity_Last);
3509 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
3511 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3512 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3515 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3516 aVec[SMDSEntity_Tetra] = nbVols;
3518 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
3519 aResMap.insert(std::make_pair(sm,aVec));