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,
247 : SMESH_2D_Algo(hypId, gen)
249 _name = theHasGEOM ? "MG-CADSurf" : "MG-CADSurf_NOGEOM";//"BLSURF";
250 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
251 _compatibleHypothesis.push_back(BLSURFPlugin_Hypothesis::GetHypType(theHasGEOM));
253 _compatibleHypothesis.push_back(StdMeshers_ViscousLayers2D::GetHypType());
254 _requireDiscreteBoundary = false;
255 _onlyUnaryInput = false;
257 _supportSubmeshes = true;
258 _requireShape = theHasGEOM;
260 /* Initialize the Python interpreter */
261 assert(Py_IsInitialized());
262 PyGILState_STATE gstate;
263 gstate = PyGILState_Ensure();
266 main_mod = PyImport_AddModule("__main__");
269 main_dict = PyModule_GetDict(main_mod);
271 PyRun_SimpleString("from math import *");
272 PyGILState_Release(gstate);
274 FacesWithSizeMap.Clear();
275 FaceId2SizeMap.clear();
276 EdgesWithSizeMap.Clear();
277 EdgeId2SizeMap.clear();
278 VerticesWithSizeMap.Clear();
279 VertexId2SizeMap.clear();
280 FaceId2PythonSmp.clear();
281 EdgeId2PythonSmp.clear();
282 VertexId2PythonSmp.clear();
283 FaceId2AttractorCoords.clear();
284 FaceId2ClassAttractor.clear();
285 FaceIndex2ClassAttractor.clear();
286 FacesWithEnforcedVertices.Clear();
287 FaceId2EnforcedVertexCoords.clear();
288 EnfVertexCoords2ProjVertex.clear();
289 EnfVertexCoords2EnfVertexList.clear();
291 _compute_canceled = false;
294 //=============================================================================
298 //=============================================================================
300 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
305 //=============================================================================
309 //=============================================================================
311 bool BLSURFPlugin_BLSURF::CheckHypothesis
313 const TopoDS_Shape& aShape,
314 SMESH_Hypothesis::Hypothesis_Status& aStatus)
317 _haveViscousLayers = false;
319 list<const SMESHDS_Hypothesis*>::const_iterator itl;
320 const SMESHDS_Hypothesis* theHyp;
322 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape,
323 /*ignoreAuxiliary=*/false);
324 aStatus = SMESH_Hypothesis::HYP_OK;
327 return true; // can work with no hypothesis
330 for ( itl = hyps.begin(); itl != hyps.end() && ( aStatus == HYP_OK ); ++itl )
333 string hypName = theHyp->GetName();
334 if ( hypName == BLSURFPlugin_Hypothesis::GetHypType(true) ||
335 hypName == BLSURFPlugin_Hypothesis::GetHypType(false) )
337 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
339 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
340 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
341 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
342 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
344 else if ( hypName == StdMeshers_ViscousLayers2D::GetHypType() )
346 if ( !_haveViscousLayers )
348 if ( error( StdMeshers_ViscousLayers2D::CheckHypothesis( aMesh, aShape, aStatus )))
349 _haveViscousLayers = true;
354 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
357 return aStatus == SMESH_Hypothesis::HYP_OK;
360 //=============================================================================
362 * Pass parameters to MG-CADSurf
364 //=============================================================================
366 inline std::string val_to_string(double d)
368 std::ostringstream o;
373 inline std::string val_to_string_rel(double d)
375 std::ostringstream o;
381 inline std::string val_to_string(int i)
383 std::ostringstream o;
388 inline std::string val_to_string_rel(int i)
390 std::ostringstream o;
396 double _smp_phy_size;
397 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
398 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
399 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
406 /////////////////////////////////////////////////////////
408 projectionPoint getProjectionPoint(TopoDS_Face& theFace, const gp_Pnt& thePoint)
410 projectionPoint myPoint;
412 if ( theFace.IsNull() )
414 TopoDS_Shape foundFace, myShape = theHelper->GetSubShape();
415 TopTools_MapOfShape checkedFaces;
416 std::map< double, std::pair< TopoDS_Face, gp_Pnt2d > > dist2face;
418 for ( TopExp_Explorer exp ( myShape, TopAbs_FACE ); exp.More(); exp.Next())
420 const TopoDS_Face& face = TopoDS::Face( exp.Current() );
421 if ( !checkedFaces.Add( face )) continue;
423 // check distance to face
424 Handle(ShapeAnalysis_Surface) surface = theHelper->GetSurface( face );
425 gp_Pnt2d uv = surface->ValueOfUV( thePoint, Precision::Confusion());
426 double distance = surface->Gap();
427 if ( distance > Precision::Confusion() )
429 // the face is far, store for future analysis
430 dist2face.insert( std::make_pair( distance, std::make_pair( face, uv )));
434 // check location on the face
435 BRepClass_FaceClassifier FC( face, uv, BRep_Tool::Tolerance( face ));
436 if ( FC.State() == TopAbs_IN )
438 if ( !foundFace.IsNull() )
439 return myPoint; // thePoint seems to be TopAbs_ON
441 myPoint.uv = uv.XY();
442 myPoint.xyz = surface->Value( uv ).XYZ();
445 if ( FC.State() == TopAbs_ON )
449 if ( foundFace.IsNull() )
451 // find the closest face
452 std::map< double, std::pair< TopoDS_Face, gp_Pnt2d > >::iterator d2f = dist2face.begin();
453 for ( ; d2f != dist2face.end(); ++d2f )
455 const TopoDS_Face& face = d2f->second.first;
456 const gp_Pnt2d & uv = d2f->second.second;
457 BRepClass_FaceClassifier FC( face, uv, Precision::Confusion());
458 if ( FC.State() == TopAbs_IN )
461 myPoint.uv = uv.XY();
462 myPoint.xyz = theHelper->GetSurface( face )->Value( uv ).XYZ();
467 // set the resultShape
468 // if ( foundFace.IsNull() )
469 // throw SMESH_ComputeError(COMPERR_BAD_PARMETERS,
470 // "getProjectionPoint: can't find a face by a vertex");
471 theFace = TopoDS::Face( foundFace );
475 Handle(Geom_Surface) surface = BRep_Tool::Surface( theFace );
476 GeomAPI_ProjectPointOnSurf projector( thePoint, surface );
477 if ( !projector.IsDone() || projector.NbPoints()==0 )
478 throw SMESH_ComputeError(COMPERR_BAD_PARMETERS,
479 "getProjectionPoint: can't project a vertex to a face");
481 Quantity_Parameter u,v;
482 projector.LowerDistanceParameters(u,v);
483 myPoint.uv = gp_XY(u,v);
484 gp_Pnt aPnt = projector.NearestPoint();
485 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
487 BRepClass_FaceClassifier FC( theFace, myPoint.uv, Precision::Confusion());
488 if ( FC.State() != TopAbs_IN )
495 /////////////////////////////////////////////////////////
496 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
498 GEOM::GEOM_Object_var aGeomObj;
499 TopoDS_Shape S = TopoDS_Shape();
500 SALOMEDS::SObject_var aSObj = SMESH_Gen_i::getStudyServant()->FindObjectID( entry.c_str() );
501 if (!aSObj->_is_nil()) {
502 CORBA::Object_var obj = aSObj->GetObject();
503 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
506 if ( !aGeomObj->_is_nil() )
507 S = SMESH_Gen_i::GetSMESHGen()->GeomObjectToShape( aGeomObj.in() );
511 void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
513 BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
515 // Find the face and get the (u,v) values of the enforced vertex on the face
516 projectionPoint myPoint = getProjectionPoint(faceShape,aPnt);
517 if ( faceShape.IsNull() )
520 enf_coords.push_back(aPnt.X());
521 enf_coords.push_back(aPnt.Y());
522 enf_coords.push_back(aPnt.Z());
524 coords.push_back(myPoint.uv.X());
525 coords.push_back(myPoint.uv.Y());
526 coords.push_back(myPoint.xyz.X());
527 coords.push_back(myPoint.xyz.Y());
528 coords.push_back(myPoint.xyz.Z());
530 s_coords.push_back(myPoint.xyz.X());
531 s_coords.push_back(myPoint.xyz.Y());
532 s_coords.push_back(myPoint.xyz.Z());
534 // Save pair projected vertex / enf vertex
535 EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
536 pair<BLSURFPlugin_Hypothesis::TEnfVertexList::iterator,bool> ret;
537 BLSURFPlugin_Hypothesis::TEnfVertexList::iterator it;
538 ret = EnfVertexCoords2EnfVertexList[s_coords].insert(enfVertex);
539 if (ret.second == false) {
541 (*it)->grpName = enfVertex->grpName;
545 if (! FacesWithEnforcedVertices.Contains(faceShape)) {
546 key = FacesWithEnforcedVertices.Add(faceShape);
549 key = FacesWithEnforcedVertices.FindIndex(faceShape);
552 // If a node is already created by an attractor, do not create enforced vertex
553 int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
554 bool sameAttractor = false;
555 if (attractorKey >= 0)
556 if (FaceId2AttractorCoords.count(attractorKey) > 0)
557 if (FaceId2AttractorCoords[attractorKey] == coords)
558 sameAttractor = true;
560 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
562 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
565 if (! sameAttractor) {
566 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList ens;
568 FaceId2EnforcedVertexCoords[key] = ens;
573 /////////////////////////////////////////////////////////
574 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
576 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
579 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
581 for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
582 enfVertex = *enfVertexListIt;
583 // Case of manual coords
584 if (enfVertex->coords.size() != 0) {
585 aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
586 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
589 // Case of geom vertex coords
590 if (enfVertex->geomEntry != "") {
591 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
592 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
593 if (GeomType == TopAbs_VERTEX)
595 enfVertex->vertex = TopoDS::Vertex( GeomShape );
596 aPnt = BRep_Tool::Pnt( enfVertex->vertex );
597 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
600 if (GeomType == TopAbs_COMPOUND)
602 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next())
603 if (it.Value().ShapeType() == TopAbs_VERTEX)
605 enfVertex->vertex = TopoDS::Vertex( it.Value() );
606 aPnt = BRep_Tool::Pnt( enfVertex->vertex );
607 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
614 /////////////////////////////////////////////////////////
615 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction, double defaultSize)
617 double xa, ya, za; // Coordinates of attractor point
618 double a, b; // Attractor parameter
620 bool createNode=false; // To create a node on attractor projection
622 const char *sep = ";";
623 // atIt->second has the following pattern:
624 // ATTRACTOR(xa;ya;za;a;b;True|False;d)
626 // xa;ya;za : coordinates of attractor
627 // a : desired size on attractor
628 // b : distance of influence of attractor
629 // d : distance until which the size remains constant
631 // We search the parameters in the string
633 pos1 = AttractorFunction.find(sep);
634 if (pos1!=string::npos)
635 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
637 pos2 = AttractorFunction.find(sep, pos1+1);
638 if (pos2!=string::npos) {
639 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
643 pos2 = AttractorFunction.find(sep, pos1+1);
644 if (pos2!=string::npos) {
645 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
649 pos2 = AttractorFunction.find(sep, pos1+1);
650 if (pos2!=string::npos) {
651 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
655 pos2 = AttractorFunction.find(sep, pos1+1);
656 if (pos2!=string::npos) {
657 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
661 pos2 = AttractorFunction.find(sep, pos1+1);
662 if (pos2!=string::npos) {
663 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
664 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
668 pos2 = AttractorFunction.find(")");
669 if (pos2!=string::npos) {
670 d = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
673 // Get the (u,v) values of the attractor on the face
674 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
675 gp_XY uvPoint = myPoint.uv;
676 gp_XYZ xyzPoint = myPoint.xyz;
677 Standard_Real u0 = uvPoint.X();
678 Standard_Real v0 = uvPoint.Y();
679 Standard_Real x0 = xyzPoint.X();
680 Standard_Real y0 = xyzPoint.Y();
681 Standard_Real z0 = xyzPoint.Z();
682 std::vector<double> coords;
683 coords.push_back(u0);
684 coords.push_back(v0);
685 coords.push_back(x0);
686 coords.push_back(y0);
687 coords.push_back(z0);
688 // We construct the python function
689 ostringstream attractorFunctionStream;
690 attractorFunctionStream << "def f(u,v): return ";
691 attractorFunctionStream << defaultSize << "-(" << defaultSize <<"-" << a << ")";
692 //attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
693 // rnc: make possible to keep the size constant until
694 // a defined distance. Distance is expressed as the positiv part
695 // of r-d where r is the distance to (u0,v0)
696 attractorFunctionStream << "*exp(-(0.5*(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"+abs(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"))/(" << b << "))**2)";
699 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
700 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
703 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
705 FaceId2SizeMap[key] =attractorFunctionStream.str();
707 FaceId2AttractorCoords[key] = coords;
709 // // Test for new attractors
710 // gp_Pnt myP(xyzPoint);
711 // TopoDS_Vertex myV = BRepBuilderAPI_MakeVertex(myP);
712 // BLSURFPlugin_Attractor myAttractor(TopoDS::Face(GeomShape),myV,200);
713 // myAttractor.SetParameters(a, defaultSize, b, d);
714 // myAttractor.SetType(1);
715 // FaceId2ClassAttractor[key] = myAttractor;
716 // if(FaceId2ClassAttractor[key].GetFace().IsNull()){
720 // One sub-shape to get ids from
721 BLSURFPlugin_BLSURF::TListOfIDs _getSubShapeIDsInMainShape(const TopoDS_Shape& theMainShape,
722 const TopoDS_Shape& theSubShape,
723 TopAbs_ShapeEnum theShapeType)
725 BLSURFPlugin_BLSURF::TListOfIDs face_ids;
727 TopTools_MapOfShape subShapes;
728 TopTools_IndexedMapOfShape anIndices;
729 TopExp::MapShapes(theMainShape, theShapeType, anIndices);
731 for (TopExp_Explorer face_iter(theSubShape,theShapeType);face_iter.More();face_iter.Next())
733 if ( subShapes.Add( face_iter.Current() )) // issue 23416
735 int face_id = anIndices.FindIndex( face_iter.Current() );
737 throw SALOME_Exception( "Periodicity: sub_shape not found in main_shape");
738 face_ids.push_back( face_id );
744 BLSURFPlugin_BLSURF::TListOfIDs _getSubShapeIDsInMainShape(SMESH_Mesh* theMesh,
745 TopoDS_Shape theSubShape,
746 TopAbs_ShapeEnum theShapeType)
748 BLSURFPlugin_BLSURF::TListOfIDs face_ids;
750 for (TopExp_Explorer face_iter(theSubShape,theShapeType);face_iter.More();face_iter.Next())
752 int face_id = theMesh->GetMeshDS()->ShapeToIndex(face_iter.Current());
754 throw SALOME_Exception ( "Periodicity: sub_shape not found in main_shape");
755 face_ids.push_back(face_id);
761 void BLSURFPlugin_BLSURF::addCoordsFromVertices(const std::vector<std::string> &theVerticesEntries, std::vector<double> &theVerticesCoords)
763 for (std::vector<std::string>::const_iterator it = theVerticesEntries.begin(); it != theVerticesEntries.end(); it++)
765 BLSURFPlugin_Hypothesis::TEntry theVertexEntry = *it;
766 addCoordsFromVertex(theVertexEntry, theVerticesCoords);
771 void BLSURFPlugin_BLSURF::addCoordsFromVertex(BLSURFPlugin_Hypothesis::TEntry theVertexEntry, std::vector<double> &theVerticesCoords)
773 if (theVertexEntry!="")
775 TopoDS_Shape aShape = entryToShape(theVertexEntry);
777 gp_Pnt aPnt = BRep_Tool::Pnt( TopoDS::Vertex( aShape ) );
778 double theX, theY, theZ;
783 theVerticesCoords.push_back(theX);
784 theVerticesCoords.push_back(theY);
785 theVerticesCoords.push_back(theZ);
789 /////////////////////////////////////////////////////////
790 void BLSURFPlugin_BLSURF::createPreCadFacesPeriodicity(TopoDS_Shape theGeomShape, const BLSURFPlugin_Hypothesis::TPreCadPeriodicity &preCadPeriodicity)
792 TopoDS_Shape geomShape1 = entryToShape(preCadPeriodicity.shape1Entry);
793 TopoDS_Shape geomShape2 = entryToShape(preCadPeriodicity.shape2Entry);
795 TListOfIDs theFace1_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape1, TopAbs_FACE);
796 TListOfIDs theFace2_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape2, TopAbs_FACE);
798 TPreCadPeriodicityIDs preCadFacesPeriodicityIDs;
799 preCadFacesPeriodicityIDs.shape1IDs = theFace1_ids;
800 preCadFacesPeriodicityIDs.shape2IDs = theFace2_ids;
802 addCoordsFromVertices(preCadPeriodicity.theSourceVerticesEntries, preCadFacesPeriodicityIDs.theSourceVerticesCoords);
803 addCoordsFromVertices(preCadPeriodicity.theTargetVerticesEntries, preCadFacesPeriodicityIDs.theTargetVerticesCoords);
805 _preCadFacesIDsPeriodicityVector.push_back(preCadFacesPeriodicityIDs);
808 /////////////////////////////////////////////////////////
809 void BLSURFPlugin_BLSURF::createPreCadEdgesPeriodicity(TopoDS_Shape theGeomShape, const BLSURFPlugin_Hypothesis::TPreCadPeriodicity &preCadPeriodicity)
811 TopoDS_Shape geomShape1 = entryToShape(preCadPeriodicity.shape1Entry);
812 TopoDS_Shape geomShape2 = entryToShape(preCadPeriodicity.shape2Entry);
814 TListOfIDs theEdge1_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape1, TopAbs_EDGE);
815 TListOfIDs theEdge2_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape2, TopAbs_EDGE);
817 TPreCadPeriodicityIDs preCadEdgesPeriodicityIDs;
818 preCadEdgesPeriodicityIDs.shape1IDs = theEdge1_ids;
819 preCadEdgesPeriodicityIDs.shape2IDs = theEdge2_ids;
821 addCoordsFromVertices(preCadPeriodicity.theSourceVerticesEntries, preCadEdgesPeriodicityIDs.theSourceVerticesCoords);
822 addCoordsFromVertices(preCadPeriodicity.theTargetVerticesEntries, preCadEdgesPeriodicityIDs.theTargetVerticesCoords);
824 _preCadEdgesIDsPeriodicityVector.push_back(preCadEdgesPeriodicityIDs);
828 /////////////////////////////////////////////////////////
830 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
831 cadsurf_session_t * css,
832 const TopoDS_Shape& theGeomShape)
835 // Clear map so that it is not stored in the algorithm with old enforced vertices in it
836 FacesWithSizeMap.Clear();
837 FaceId2SizeMap.clear();
838 EdgesWithSizeMap.Clear();
839 EdgeId2SizeMap.clear();
840 VerticesWithSizeMap.Clear();
841 VertexId2SizeMap.clear();
842 FaceId2PythonSmp.clear();
843 EdgeId2PythonSmp.clear();
844 VertexId2PythonSmp.clear();
845 FaceId2AttractorCoords.clear();
846 FaceId2ClassAttractor.clear();
847 FaceIndex2ClassAttractor.clear();
848 FacesWithEnforcedVertices.Clear();
849 FaceId2EnforcedVertexCoords.clear();
850 EnfVertexCoords2ProjVertex.clear();
851 EnfVertexCoords2EnfVertexList.clear();
853 double diagonal = SMESH_Mesh::GetShapeDiagonalSize( theGeomShape );
854 double bbSegmentation = _gen->GetBoundaryBoxSegmentation();
855 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
856 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
857 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize(diagonal, bbSegmentation);
858 bool _phySizeRel = BLSURFPlugin_Hypothesis::GetDefaultPhySizeRel();
859 double _minSize = BLSURFPlugin_Hypothesis::GetDefaultMinSize(diagonal);
860 bool _minSizeRel = BLSURFPlugin_Hypothesis::GetDefaultMinSizeRel();
861 double _maxSize = BLSURFPlugin_Hypothesis::GetDefaultMaxSize(diagonal);
862 bool _maxSizeRel = BLSURFPlugin_Hypothesis::GetDefaultMaxSizeRel();
863 double _use_gradation = BLSURFPlugin_Hypothesis::GetDefaultUseGradation();
864 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
865 double _use_volume_gradation = BLSURFPlugin_Hypothesis::GetDefaultUseVolumeGradation();
866 double _volume_gradation = BLSURFPlugin_Hypothesis::GetDefaultVolumeGradation();
867 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
868 double _angleMesh = BLSURFPlugin_Hypothesis::GetDefaultAngleMesh();
869 double _chordalError = BLSURFPlugin_Hypothesis::GetDefaultChordalError(diagonal);
870 bool _anisotropic = BLSURFPlugin_Hypothesis::GetDefaultAnisotropic();
871 double _anisotropicRatio = BLSURFPlugin_Hypothesis::GetDefaultAnisotropicRatio();
872 bool _removeTinyEdges = BLSURFPlugin_Hypothesis::GetDefaultRemoveTinyEdges();
873 double _tinyEdgeLength = BLSURFPlugin_Hypothesis::GetDefaultTinyEdgeLength(diagonal);
874 bool _optimiseTinyEdges = BLSURFPlugin_Hypothesis::GetDefaultOptimiseTinyEdges();
875 double _tinyEdgeOptimisLength = BLSURFPlugin_Hypothesis::GetDefaultTinyEdgeOptimisationLength(diagonal);
876 bool _correctSurfaceIntersec= BLSURFPlugin_Hypothesis::GetDefaultCorrectSurfaceIntersection();
877 double _corrSurfaceIntersCost = BLSURFPlugin_Hypothesis::GetDefaultCorrectSurfaceIntersectionMaxCost();
878 bool _badElementRemoval = BLSURFPlugin_Hypothesis::GetDefaultBadElementRemoval();
879 double _badElementAspectRatio = BLSURFPlugin_Hypothesis::GetDefaultBadElementAspectRatio();
880 bool _optimizeMesh = BLSURFPlugin_Hypothesis::GetDefaultOptimizeMesh();
881 bool _quadraticMesh = BLSURFPlugin_Hypothesis::GetDefaultQuadraticMesh();
882 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
883 //int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
886 //int _precadMergeEdges = BLSURFPlugin_Hypothesis::GetDefaultPreCADMergeEdges();
887 //int _precadRemoveDuplicateCADFaces = BLSURFPlugin_Hypothesis::GetDefaultPreCADRemoveDuplicateCADFaces();
888 //int _precadProcess3DTopology = BLSURFPlugin_Hypothesis::GetDefaultPreCADProcess3DTopology();
889 //int _precadDiscardInput = BLSURFPlugin_Hypothesis::GetDefaultPreCADDiscardInput();
891 const BLSURFPlugin_Hypothesis::TPreCadPeriodicityVector preCadFacesPeriodicityVector = BLSURFPlugin_Hypothesis::GetPreCadFacesPeriodicityVector(hyp);
894 _physicalMesh = (int) hyp->GetPhysicalMesh();
895 _geometricMesh = (int) hyp->GetGeometricMesh();
896 if (hyp->GetPhySize() > 0) {
897 _phySize = hyp->GetPhySize();
898 // if user size is not explicitly specified, "relative" flag is ignored
899 _phySizeRel = hyp->IsPhySizeRel();
901 if (hyp->GetMinSize() > 0) {
902 _minSize = hyp->GetMinSize();
903 // if min size is not explicitly specified, "relative" flag is ignored
904 _minSizeRel = hyp->IsMinSizeRel();
906 if (hyp->GetMaxSize() > 0) {
907 _maxSize = hyp->GetMaxSize();
908 // if max size is not explicitly specified, "relative" flag is ignored
909 _maxSizeRel = hyp->IsMaxSizeRel();
911 _use_gradation = hyp->GetUseGradation();
912 if (hyp->GetGradation() > 0 && _use_gradation)
913 _gradation = hyp->GetGradation();
914 _use_volume_gradation = hyp->GetUseVolumeGradation();
915 if (hyp->GetVolumeGradation() > 0 && _use_volume_gradation )
916 _volume_gradation = hyp->GetVolumeGradation();
917 _quadAllowed = hyp->GetQuadAllowed();
918 if (hyp->GetAngleMesh() > 0)
919 _angleMesh = hyp->GetAngleMesh();
920 if (hyp->GetChordalError() > 0)
921 _chordalError = hyp->GetChordalError();
922 _anisotropic = hyp->GetAnisotropic();
923 if (hyp->GetAnisotropicRatio() >= 0)
924 _anisotropicRatio = hyp->GetAnisotropicRatio();
925 _removeTinyEdges = hyp->GetRemoveTinyEdges();
926 if (hyp->GetTinyEdgeLength() > 0)
927 _tinyEdgeLength = hyp->GetTinyEdgeLength();
928 _optimiseTinyEdges = hyp->GetOptimiseTinyEdges();
929 if (hyp->GetTinyEdgeOptimisationLength() > 0)
930 _tinyEdgeOptimisLength = hyp->GetTinyEdgeOptimisationLength();
931 _correctSurfaceIntersec = hyp->GetCorrectSurfaceIntersection();
932 if (hyp->GetCorrectSurfaceIntersectionMaxCost() > 0)
933 _corrSurfaceIntersCost = hyp->GetCorrectSurfaceIntersectionMaxCost();
934 _badElementRemoval = hyp->GetBadElementRemoval();
935 if (hyp->GetBadElementAspectRatio() >= 0)
936 _badElementAspectRatio = hyp->GetBadElementAspectRatio();
937 _optimizeMesh = hyp->GetOptimizeMesh();
938 _quadraticMesh = hyp->GetQuadraticMesh();
939 _verb = hyp->GetVerbosity();
940 //_topology = (int) hyp->GetTopology();
942 //_precadMergeEdges = hyp->GetPreCADMergeEdges();
943 //_precadRemoveDuplicateCADFaces = hyp->GetPreCADRemoveDuplicateCADFaces();
944 //_precadProcess3DTopology = hyp->GetPreCADProcess3DTopology();
945 //_precadDiscardInput = hyp->GetPreCADDiscardInput();
947 const BLSURFPlugin_Hypothesis::TOptionValues& opts = hyp->GetOptionValues();
948 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
949 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt ){
950 MESSAGE("OptionValue: " << opIt->first.c_str() << ", value: " << opIt->second.c_str());
951 if ( !opIt->second.empty() ) {
952 // With MeshGems 2.4-5, there are issues with periodicity and multithread
953 // => As a temporary workaround, we enforce to use only one thread if periodicity is used.
954 if (opIt->first == "max_number_of_threads" && opIt->second != "1" && ! preCadFacesPeriodicityVector.empty()){
955 std::cout << "INFO: Disabling multithread to avoid periodicity issues" << std::endl;
956 set_param(css, opIt->first.c_str(), "1");
959 set_param(css, opIt->first.c_str(), opIt->second.c_str());
963 const BLSURFPlugin_Hypothesis::TOptionValues& custom_opts = hyp->GetCustomOptionValues();
964 for ( opIt = custom_opts.begin(); opIt != custom_opts.end(); ++opIt )
965 if ( !opIt->second.empty() ) {
966 set_param(css, opIt->first.c_str(), opIt->second.c_str());
969 const BLSURFPlugin_Hypothesis::TOptionValues& preCADopts = hyp->GetPreCADOptionValues();
970 for ( opIt = preCADopts.begin(); opIt != preCADopts.end(); ++opIt )
971 if ( !opIt->second.empty() ) {
972 set_param(css, opIt->first.c_str(), opIt->second.c_str());
976 if ( BLSURFPlugin_Hypothesis::HasPreCADOptions( hyp ))
978 cadsurf_set_param(css, "use_precad", "yes" ); // for old versions
980 // PreProcessor (formerly PreCAD) -- commented params are preCADoptions (since 0023307)
981 //set_param(css, "merge_edges", _precadMergeEdges ? "yes" : "no");
982 //set_param(css, "remove_duplicate_cad_faces", _precadRemoveDuplicateCADFaces ? "yes" : "no");
983 //set_param(css, "process_3d_topology", _precadProcess3DTopology ? "1" : "0");
984 //set_param(css, "discard_input_topology", _precadDiscardInput ? "1" : "0");
985 //set_param(css, "max_number_of_points_per_patch", "1000000");
987 bool useGradation = false;
988 switch (_physicalMesh)
990 case BLSURFPlugin_Hypothesis::PhysicalGlobalSize:
991 set_param(css, "physical_size_mode", "global");
992 set_param(css, "global_physical_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
994 case BLSURFPlugin_Hypothesis::PhysicalLocalSize:
995 set_param(css, "physical_size_mode", "local");
996 set_param(css, "global_physical_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1000 set_param(css, "physical_size_mode", "none");
1003 switch (_geometricMesh)
1005 case BLSURFPlugin_Hypothesis::GeometricalGlobalSize:
1006 set_param(css, "geometric_size_mode", "global");
1007 set_param(css, "geometric_approximation", val_to_string(_angleMesh).c_str());
1008 set_param(css, "chordal_error", val_to_string(_chordalError).c_str());
1009 useGradation = true;
1011 case BLSURFPlugin_Hypothesis::GeometricalLocalSize:
1012 set_param(css, "geometric_size_mode", "local");
1013 set_param(css, "geometric_approximation", val_to_string(_angleMesh).c_str());
1014 set_param(css, "chordal_error", val_to_string(_chordalError).c_str());
1015 useGradation = true;
1018 set_param(css, "geometric_size_mode", "none");
1021 if ( hyp && hyp->GetPhySize() > 0 ) {
1022 // user size is explicitly specified via hypothesis parameters
1023 // min and max sizes should be compared with explicitly specified user size
1024 // - compute absolute min size
1025 double mins = _minSizeRel ? _minSize * diagonal : _minSize;
1026 // - min size should not be greater than user size
1027 if ( _phySize < mins )
1028 set_param(css, "min_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1030 set_param(css, "min_size", _minSizeRel ? val_to_string_rel(_minSize).c_str() : val_to_string(_minSize).c_str());
1031 // - compute absolute max size
1032 double maxs = _maxSizeRel ? _maxSize * diagonal : _maxSize;
1033 // - max size should not be less than user size
1034 if ( _phySize > maxs )
1035 set_param(css, "max_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1037 set_param(css, "max_size", _maxSizeRel ? val_to_string_rel(_maxSize).c_str() : val_to_string(_maxSize).c_str());
1040 // user size is not explicitly specified
1041 // - if minsize is not explicitly specified, we pass default value computed automatically, in this case "relative" flag is ignored
1042 set_param(css, "min_size", _minSizeRel ? val_to_string_rel(_minSize).c_str() : val_to_string(_minSize).c_str());
1043 // - if maxsize is not explicitly specified, we pass default value computed automatically, in this case "relative" flag is ignored
1044 set_param(css, "max_size", _maxSizeRel ? val_to_string_rel(_maxSize).c_str() : val_to_string(_maxSize).c_str());
1046 // anisotropic and quadrangle mesh requires disabling gradation
1047 if ( _anisotropic && _quadAllowed )
1048 useGradation = false; // limitation of V1.3
1049 if ( useGradation && _use_gradation )
1050 set_param(css, "gradation", val_to_string(_gradation).c_str());
1051 if ( useGradation && _use_volume_gradation )
1052 set_param(css, "volume_gradation", val_to_string(_volume_gradation).c_str());
1053 set_param(css, "element_generation", _quadAllowed ? "quad_dominant" : "triangle");
1056 set_param(css, "metric", _anisotropic ? "anisotropic" : "isotropic");
1058 set_param(css, "anisotropic_ratio", val_to_string(_anisotropicRatio).c_str());
1059 set_param(css, "remove_tiny_edges", _removeTinyEdges ? "1" : "0");
1060 if ( _removeTinyEdges )
1061 set_param(css, "tiny_edge_length", val_to_string(_tinyEdgeLength).c_str());
1062 set_param(css, "optimise_tiny_edges", _optimiseTinyEdges ? "1" : "0");
1063 if ( _optimiseTinyEdges )
1064 set_param(css, "tiny_edge_optimisation_length", val_to_string(_tinyEdgeOptimisLength).c_str());
1065 set_param(css, "correct_surface_intersections", _correctSurfaceIntersec ? "1" : "0");
1066 if ( _correctSurfaceIntersec )
1067 set_param(css, "surface_intersections_processing_max_cost", val_to_string(_corrSurfaceIntersCost ).c_str());
1068 set_param(css, "force_bad_surface_element_removal", _badElementRemoval ? "1" : "0");
1069 if ( _badElementRemoval )
1070 set_param(css, "bad_surface_element_aspect_ratio", val_to_string(_badElementAspectRatio).c_str());
1071 set_param(css, "optimisation", _optimizeMesh ? "yes" : "no");
1072 set_param(css, "element_order", _quadraticMesh ? "quadratic" : "linear");
1073 set_param(css, "verbose", val_to_string(_verb).c_str());
1075 _smp_phy_size = _phySizeRel ? _phySize*diagonal : _phySize;
1077 std::cout << "_smp_phy_size = " << _smp_phy_size << std::endl;
1079 if (_physicalMesh == BLSURFPlugin_Hypothesis::PhysicalLocalSize)
1081 TopoDS_Shape GeomShape;
1082 TopoDS_Shape AttShape;
1083 TopAbs_ShapeEnum GeomType;
1085 // Standard Size Maps
1087 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
1088 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
1089 for ( ; smIt != sizeMaps.end(); ++smIt ) {
1090 if ( !smIt->second.empty() ) {
1091 GeomShape = entryToShape(smIt->first);
1092 GeomType = GeomShape.ShapeType();
1095 if (GeomType == TopAbs_COMPOUND) {
1096 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1098 if (it.Value().ShapeType() == TopAbs_FACE){
1099 HasSizeMapOnFace = true;
1100 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
1101 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
1104 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
1106 FaceId2SizeMap[key] = smIt->second;
1109 if (it.Value().ShapeType() == TopAbs_EDGE){
1110 HasSizeMapOnEdge = true;
1111 HasSizeMapOnFace = true;
1112 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
1113 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
1116 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
1118 EdgeId2SizeMap[key] = smIt->second;
1120 // Group of vertices
1121 if (it.Value().ShapeType() == TopAbs_VERTEX){
1122 HasSizeMapOnVertex = true;
1123 HasSizeMapOnEdge = true;
1124 HasSizeMapOnFace = true;
1125 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
1126 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
1129 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
1131 VertexId2SizeMap[key] = smIt->second;
1136 if (GeomType == TopAbs_FACE){
1137 HasSizeMapOnFace = true;
1138 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
1139 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
1142 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
1144 FaceId2SizeMap[key] = smIt->second;
1147 if (GeomType == TopAbs_EDGE){
1148 HasSizeMapOnEdge = true;
1149 HasSizeMapOnFace = true;
1150 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
1151 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
1154 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
1156 EdgeId2SizeMap[key] = smIt->second;
1159 if (GeomType == TopAbs_VERTEX){
1160 HasSizeMapOnVertex = true;
1161 HasSizeMapOnEdge = true;
1162 HasSizeMapOnFace = true;
1163 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
1164 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
1167 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
1169 VertexId2SizeMap[key] = smIt->second;
1177 // TODO appeler le constructeur des attracteurs directement ici
1178 // if ( !_phySizeRel ) {
1179 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
1180 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
1181 for ( ; atIt != attractors.end(); ++atIt ) {
1182 if ( !atIt->second.empty() ) {
1183 GeomShape = entryToShape(atIt->first);
1184 GeomType = GeomShape.ShapeType();
1186 if (GeomType == TopAbs_COMPOUND){
1187 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1188 if (it.Value().ShapeType() == TopAbs_FACE){
1189 HasSizeMapOnFace = true;
1190 createAttractorOnFace(it.Value(), atIt->second, _phySizeRel ? _phySize*diagonal : _phySize);
1195 if (GeomType == TopAbs_FACE){
1196 HasSizeMapOnFace = true;
1197 createAttractorOnFace(GeomShape, atIt->second, _phySizeRel ? _phySize*diagonal : _phySize);
1200 if (GeomType == TopAbs_EDGE){
1201 HasSizeMapOnEdge = true;
1202 HasSizeMapOnFace = true;
1203 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
1205 if (GeomType == TopAbs_VERTEX){
1206 HasSizeMapOnVertex = true;
1207 HasSizeMapOnEdge = true;
1208 HasSizeMapOnFace = true;
1209 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
1217 // temporary commented out for testing
1219 // - Fill in the BLSURFPlugin_Hypothesis::TAttractorMap map in the hypothesis
1220 // - Uncomment and complete this part to construct the attractors from the attractor shape and the passed parameters on each face of the map
1221 // - To do this use the public methodss: SetParameters(several double parameters) and SetType(int type)
1223 // - Construct the attractors with an empty dist. map in the hypothesis
1224 // - 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()
1225 // -> define a bool _mapbuilt in the class that is set to false by default and set to true when calling _buildmap() OK
1227 theNbAttractors = 0;
1228 const BLSURFPlugin_Hypothesis::TAttractorMap class_attractors = BLSURFPlugin_Hypothesis::GetClassAttractorEntries(hyp);
1230 BLSURFPlugin_Hypothesis::TAttractorMap::const_iterator AtIt = class_attractors.begin();
1231 for ( ; AtIt != class_attractors.end(); ++AtIt ) {
1232 if ( !AtIt->second->Empty() ) {
1233 GeomShape = entryToShape(AtIt->first);
1234 if ( !SMESH_MesherHelper::IsSubShape( GeomShape, theGeomShape ))
1236 AttShape = AtIt->second->GetAttractorShape();
1237 GeomType = GeomShape.ShapeType();
1239 // if (GeomType == TopAbs_COMPOUND){
1240 // for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1241 // if (it.Value().ShapeType() == TopAbs_FACE){
1242 // HasAttractorOnFace = true;
1243 // myAttractor = BLSURFPluginAttractor(GeomShape, AttShape);
1248 if (GeomType == TopAbs_FACE
1249 && (AttShape.ShapeType() == TopAbs_VERTEX
1250 || AttShape.ShapeType() == TopAbs_EDGE
1251 || AttShape.ShapeType() == TopAbs_WIRE
1252 || AttShape.ShapeType() == TopAbs_COMPOUND) ){
1253 HasSizeMapOnFace = true;
1255 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
1257 FaceId2ClassAttractor[key].push_back( AtIt->second );
1261 MESSAGE("Wrong shape type !!")
1269 // Enforced Vertices
1271 const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
1272 BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
1273 for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
1274 if ( !enfIt->second.empty() ) {
1275 GeomShape = entryToShape(enfIt->first);
1276 if ( GeomShape.IsNull() )
1278 createEnforcedVertexOnFace( GeomShape, enfIt->second );
1281 else if ( GeomShape.ShapeType() == TopAbs_COMPOUND)
1283 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1284 if (it.Value().ShapeType() == TopAbs_FACE){
1285 HasSizeMapOnFace = true;
1286 createEnforcedVertexOnFace(it.Value(), enfIt->second);
1290 else if ( GeomShape.ShapeType() == TopAbs_FACE)
1292 HasSizeMapOnFace = true;
1293 createEnforcedVertexOnFace(GeomShape, enfIt->second);
1298 // Internal vertices
1299 bool useInternalVertexAllFaces = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFaces(hyp);
1300 if (useInternalVertexAllFaces) {
1301 std::string grpName = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFacesGroup(hyp);
1303 TopExp_Explorer exp (theGeomShape, TopAbs_FACE);
1304 for (; exp.More(); exp.Next()){
1305 TopExp_Explorer exp_face (exp.Current(), TopAbs_VERTEX, TopAbs_EDGE);
1306 for (; exp_face.More(); exp_face.Next())
1308 // Get coords of vertex
1309 // Check if current coords is already in enfVertexList
1310 // If coords not in enfVertexList, add new enfVertex
1311 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(exp_face.Current()));
1312 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex = new BLSURFPlugin_Hypothesis::TEnfVertex();
1313 enfVertex->coords.push_back(aPnt.X());
1314 enfVertex->coords.push_back(aPnt.Y());
1315 enfVertex->coords.push_back(aPnt.Z());
1316 enfVertex->name = "";
1317 enfVertex->faceEntries.clear();
1318 enfVertex->geomEntry = "";
1319 enfVertex->grpName = grpName;
1320 enfVertex->vertex = TopoDS::Vertex( exp_face.Current() );
1321 _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()), aPnt, enfVertex);
1322 HasSizeMapOnFace = true;
1327 cadsurf_set_sizemap_iso_cad_face(css, size_on_surface, &_smp_phy_size);
1329 if (HasSizeMapOnEdge){
1330 cadsurf_set_sizemap_iso_cad_edge(css, size_on_edge, &_smp_phy_size);
1332 if (HasSizeMapOnVertex){
1333 cadsurf_set_sizemap_iso_cad_point(css, size_on_vertex, &_smp_phy_size);
1340 _preCadFacesIDsPeriodicityVector.clear();
1341 _preCadEdgesIDsPeriodicityVector.clear();
1343 for (std::size_t i = 0; i<preCadFacesPeriodicityVector.size(); i++){
1344 createPreCadFacesPeriodicity(theGeomShape, preCadFacesPeriodicityVector[i]);
1347 const BLSURFPlugin_Hypothesis::TPreCadPeriodicityVector preCadEdgesPeriodicityVector = BLSURFPlugin_Hypothesis::GetPreCadEdgesPeriodicityVector(hyp);
1349 for (std::size_t i = 0; i<preCadEdgesPeriodicityVector.size(); i++){
1350 createPreCadEdgesPeriodicity(theGeomShape, preCadEdgesPeriodicityVector[i]);
1354 //================================================================================
1356 * \brief Throws an exception if a parameter name is wrong
1358 //================================================================================
1360 void BLSURFPlugin_BLSURF::set_param(cadsurf_session_t *css,
1361 const char * option_name,
1362 const char * option_value)
1364 status_t status = cadsurf_set_param(css, option_name, option_value );
1366 if ( _hypothesis && _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
1367 cout << option_name << " = " << option_value << endl;
1369 if ( status != MESHGEMS_STATUS_OK )
1371 if ( status == MESHGEMS_STATUS_UNKNOWN_PARAMETER ) {
1372 throw SALOME_Exception
1373 ( SMESH_Comment("Invalid name of CADSURF parameter: ") << option_name );
1375 else if ( status == MESHGEMS_STATUS_NOLICENSE )
1376 throw SALOME_Exception
1377 ( "No valid license available" );
1379 throw SALOME_Exception
1380 ( SMESH_Comment("Either wrong name or unacceptable value of CADSURF parameter '")
1381 << option_name << "': " << option_value);
1387 // --------------------------------------------------------------------------
1389 * \brief Class correctly terminating usage of MG-CADSurf library at destruction
1391 struct BLSURF_Cleaner
1394 cadsurf_session_t* _css;
1398 BLSURF_Cleaner(context_t * ctx,
1399 cadsurf_session_t* css=0,
1410 Clean( /*exceptContext=*/false );
1412 void Clean(const bool exceptContext)
1416 cadsurf_session_delete(_css); _css = 0;
1418 // #if BLSURF_VERSION_LONG >= "3.1.1"
1419 // // if(geo_sizemap_e)
1420 // // distene_sizemap_delete(geo_sizemap_e);
1421 // // if(geo_sizemap_f)
1422 // // distene_sizemap_delete(geo_sizemap_f);
1423 // if(iso_sizemap_p)
1424 // distene_sizemap_delete(iso_sizemap_p);
1425 // if(iso_sizemap_e)
1426 // distene_sizemap_delete(iso_sizemap_e);
1427 // if(iso_sizemap_f)
1428 // distene_sizemap_delete(iso_sizemap_f);
1430 // // if(clean_geo_sizemap_e)
1431 // // distene_sizemap_delete(clean_geo_sizemap_e);
1432 // // if(clean_geo_sizemap_f)
1433 // // distene_sizemap_delete(clean_geo_sizemap_f);
1434 // if(clean_iso_sizemap_p)
1435 // distene_sizemap_delete(clean_iso_sizemap_p);
1436 // if(clean_iso_sizemap_e)
1437 // distene_sizemap_delete(clean_iso_sizemap_e);
1438 // if(clean_iso_sizemap_f)
1439 // distene_sizemap_delete(clean_iso_sizemap_f);
1442 cad_delete(_cad); _cad = 0;
1443 dcad_delete(_dcad); _dcad = 0;
1444 if ( !exceptContext )
1446 context_delete(_ctx); _ctx = 0;
1452 // --------------------------------------------------------------------------
1453 // comparator to sort nodes and sub-meshes
1454 struct ShapeTypeCompare
1456 // sort nodes by position in the following order:
1457 // SMDS_TOP_FACE=2, SMDS_TOP_EDGE=1, SMDS_TOP_VERTEX=0, SMDS_TOP_3DSPACE=3
1458 bool operator()( const SMDS_MeshNode* n1, const SMDS_MeshNode* n2 ) const
1460 // NEW ORDER: nodes earlier added to sub-mesh are considered "less"
1461 return n1->getIdInShape() < n2->getIdInShape();
1462 // SMDS_TypeOfPosition pos1 = n1->GetPosition()->GetTypeOfPosition();
1463 // SMDS_TypeOfPosition pos2 = n2->GetPosition()->GetTypeOfPosition();
1464 // if ( pos1 == pos2 ) return 0;
1465 // if ( pos1 < pos2 || pos1 == SMDS_TOP_3DSPACE ) return 1;
1468 // sort sub-meshes in order: EDGE, VERTEX
1469 bool operator()( const SMESHDS_SubMesh* s1, const SMESHDS_SubMesh* s2 ) const
1471 int isVertex1 = ( s1 && s1->NbElements() == 0 );
1472 int isVertex2 = ( s2 && s2->NbElements() == 0 );
1473 if ( isVertex1 == isVertex2 )
1475 return isVertex1 < isVertex2;
1479 //================================================================================
1481 * \brief Fills groups of nodes to be merged
1483 //================================================================================
1485 void getNodeGroupsToMerge( const SMESHDS_SubMesh* smDS,
1486 const TopoDS_Shape& shape,
1487 SMESH_MeshEditor::TListOfListOfNodes& nodeGroupsToMerge)
1489 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1490 switch ( shape.ShapeType() )
1492 case TopAbs_VERTEX: {
1493 std::list< const SMDS_MeshNode* > nodes;
1494 while ( nIt->more() )
1495 nodes.push_back( nIt->next() );
1496 if ( nodes.size() > 1 )
1497 nodeGroupsToMerge.push_back( nodes );
1501 std::multimap< double, const SMDS_MeshNode* > u2node;
1502 const SMDS_EdgePosition* ePos;
1503 while ( nIt->more() )
1505 const SMDS_MeshNode* n = nIt->next();
1506 if (( ePos = dynamic_cast< const SMDS_EdgePosition* >( n->GetPosition() )))
1507 u2node.insert( make_pair( ePos->GetUParameter(), n ));
1509 if ( u2node.size() < 2 ) return;
1511 //double tol = (( u2node.rbegin()->first - u2node.begin()->first ) / 20.) / u2node.size();
1513 BRep_Tool::Range( TopoDS::Edge( shape ), f,l );
1514 double tol = (( l - f ) / 20.) / u2node.size();
1516 std::multimap< double, const SMDS_MeshNode* >::iterator un2, un1;
1517 for ( un2 = u2node.begin(), un1 = un2++; un2 != u2node.end(); un1 = un2++ )
1519 if (( un2->first - un1->first ) <= tol )
1521 std::list< const SMDS_MeshNode* > nodes;
1522 nodes.push_back( un1->second );
1523 while (( un2->first - un1->first ) <= tol )
1525 nodes.push_back( un2->second );
1526 if ( ++un2 == u2node.end()) {
1531 // make nodes created on the boundary of viscous layer replace nodes created
1532 // by MG-CADSurf as their SMDS_Position is more correct
1533 nodes.sort( ShapeTypeCompare() );
1534 nodeGroupsToMerge.push_back( nodes );
1541 // SMESH_MeshEditor::TListOfListOfNodes::const_iterator nll = nodeGroupsToMerge.begin();
1542 // for ( ; nll != nodeGroupsToMerge.end(); ++nll )
1544 // cout << "Merge ";
1545 // const std::list< const SMDS_MeshNode* >& nl = *nll;
1546 // std::list< const SMDS_MeshNode* >::const_iterator nIt = nl.begin();
1547 // for ( ; nIt != nl.end(); ++nIt )
1548 // cout << (*nIt) << " ";
1554 //================================================================================
1556 * \brief A temporary mesh used to compute mesh on a proxy FACE
1558 //================================================================================
1560 struct TmpMesh: public SMESH_Mesh
1562 typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
1563 TN2NMap _tmp2origNN;
1564 TopoDS_Face _proxyFace;
1568 _myMeshDS = new SMESHDS_Mesh( _id, true );
1570 //--------------------------------------------------------------------------------
1572 * \brief Creates a FACE bound by viscous layers and mesh each its EDGE with 1 segment
1574 //--------------------------------------------------------------------------------
1576 const TopoDS_Face& makeProxyFace( SMESH_ProxyMesh::Ptr& viscousMesh,
1577 const TopoDS_Face& origFace)
1579 SMESH_Mesh* origMesh = viscousMesh->GetMesh();
1581 SMESH_MesherHelper helper( *origMesh );
1582 helper.SetSubShape( origFace );
1583 const bool hasSeam = helper.HasRealSeam();
1585 // get data of nodes on inner boundary of viscous layers
1587 TSideVector wireVec = StdMeshers_FaceSide::GetFaceWires(origFace, *origMesh,
1588 /*skipMediumNodes = */true,
1589 err, &helper, viscousMesh );
1590 if ( err && err->IsKO() )
1591 throw *err.get(); // it should be caught at SMESH_subMesh
1593 // proxy nodes and corresponding tmp VERTEXes
1594 std::vector<const SMDS_MeshNode*> origNodes;
1595 std::vector<TopoDS_Vertex> tmpVertex;
1597 // create a proxy FACE
1598 TopoDS_Face origFaceCopy = TopoDS::Face( origFace.EmptyCopied() );
1599 BRepBuilderAPI_MakeFace newFace( origFaceCopy );
1600 bool hasPCurves = false;
1601 for ( size_t iW = 0; iW != wireVec.size(); ++iW )
1603 StdMeshers_FaceSidePtr& wireData = wireVec[iW];
1604 const UVPtStructVec& wirePoints = wireData->GetUVPtStruct();
1605 if ( wirePoints.size() < 3 )
1608 BRepBuilderAPI_MakePolygon polygon;
1609 const size_t i0 = tmpVertex.size();
1610 for ( size_t iN = 0; iN < wirePoints.size(); ++iN )
1612 polygon.Add( SMESH_TNodeXYZ( wirePoints[ iN ].node ));
1613 origNodes.push_back( wirePoints[ iN ].node );
1614 tmpVertex.push_back( polygon.LastVertex() );
1616 // check presence of a pcurve
1617 checkPCurve( polygon, origFaceCopy, hasPCurves, &wirePoints[ iN-1 ] );
1619 tmpVertex[ i0 ] = polygon.FirstVertex(); // polygon.LastVertex()==NULL for 1 vertex in wire
1621 if ( !polygon.IsDone() )
1622 throw SALOME_Exception("BLSURFPlugin_BLSURF: BRepBuilderAPI_MakePolygon failed");
1623 TopoDS_Wire wire = polygon;
1625 wire = updateSeam( wire, origNodes );
1626 newFace.Add( wire );
1628 _proxyFace = newFace;
1630 // set a new shape to mesh
1631 TopoDS_Compound auxCompoundToMesh;
1632 BRep_Builder shapeBuilder;
1633 shapeBuilder.MakeCompound( auxCompoundToMesh );
1634 shapeBuilder.Add( auxCompoundToMesh, _proxyFace );
1635 shapeBuilder.Add( auxCompoundToMesh, origMesh->GetShapeToMesh() );
1637 ShapeToMesh( auxCompoundToMesh );
1640 // Make input mesh for MG-CADSurf: segments on EDGE's of newFace
1642 // make nodes and fill in _tmp2origNN
1644 SMESHDS_Mesh* tmpMeshDS = GetMeshDS();
1645 for ( size_t i = 0; i < origNodes.size(); ++i )
1647 GetSubMesh( tmpVertex[i] )->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1648 if ( const SMDS_MeshNode* tmpN = SMESH_Algo::VertexNode( tmpVertex[i], tmpMeshDS ))
1649 _tmp2origNN.insert( _tmp2origNN.end(), make_pair( tmpN, origNodes[i] ));
1650 // else -- it can be a seam vertex replaced by updateSeam()
1651 // throw SALOME_Exception("BLSURFPlugin_BLSURF: a proxy vertex not meshed");
1655 TopoDS_Vertex v1, v2;
1656 for ( TopExp_Explorer edge( _proxyFace, TopAbs_EDGE ); edge.More(); edge.Next() )
1658 const TopoDS_Edge& E = TopoDS::Edge( edge.Current() );
1659 TopExp::Vertices( E, v1, v2 );
1660 const SMDS_MeshNode* n1 = SMESH_Algo::VertexNode( v1, tmpMeshDS );
1661 const SMDS_MeshNode* n2 = SMESH_Algo::VertexNode( v2, tmpMeshDS );
1663 if ( SMDS_MeshElement* seg = tmpMeshDS->AddEdge( n1, n2 ))
1664 tmpMeshDS->SetMeshElementOnShape( seg, E );
1670 //--------------------------------------------------------------------------------
1672 * \brief Add pcurve to the last edge of a wire
1674 //--------------------------------------------------------------------------------
1676 void checkPCurve( BRepBuilderAPI_MakePolygon& wire,
1677 const TopoDS_Face& face,
1679 const uvPtStruct * wirePoints )
1683 TopoDS_Edge edge = wire.Edge();
1684 if ( edge.IsNull() ) return;
1686 if ( BRep_Tool::CurveOnSurface(edge, face, f, l))
1691 gp_XY p1 = wirePoints[ 0 ].UV(), p2 = wirePoints[ 1 ].UV();
1692 Handle(Geom2d_Line) pcurve = new Geom2d_Line( p1, gp_Dir2d( p2 - p1 ));
1693 BRep_Builder().UpdateEdge( edge, Handle(Geom_Curve)(), Precision::Confusion() );
1694 BRep_Builder().UpdateEdge( edge, pcurve, face, Precision::Confusion() );
1695 BRep_Builder().Range( edge, 0, ( p2 - p1 ).Modulus() );
1696 // cout << "n1 = mesh.AddNode( " << p1.X()*10 << ", " << p1.Y() << ", 0 )" << endl
1697 // << "n2 = mesh.AddNode( " << p2.X()*10 << ", " << p2.Y() << ", 0 )" << endl
1698 // << "mesh.AddEdge( [ n1, n2 ] )" << endl;
1701 //--------------------------------------------------------------------------------
1703 * \brief Replace coincident EDGEs with reversed copies.
1705 //--------------------------------------------------------------------------------
1707 TopoDS_Wire updateSeam( const TopoDS_Wire& wire,
1708 const std::vector<const SMDS_MeshNode*>& nodesOfVertices )
1710 BRepBuilderAPI_MakeWire newWire;
1712 typedef NCollection_DataMap<SMESH_TLink, TopoDS_Edge, SMESH_TLink > TSeg2EdgeMap;
1713 TSeg2EdgeMap seg2EdgeMap;
1715 TopoDS_Iterator edgeIt( wire );
1716 for ( int iSeg = 1; edgeIt.More(); edgeIt.Next(), ++iSeg )
1718 SMESH_TLink link( nodesOfVertices[ iSeg-1 ], nodesOfVertices[ iSeg ]);
1719 TopoDS_Edge edge( TopoDS::Edge( edgeIt.Value() ));
1721 TopoDS_Edge* edgeInMap = seg2EdgeMap.Bound( link, edge );
1722 bool isSeam = ( *edgeInMap != edge );
1725 edgeInMap->Reverse();
1728 newWire.Add( edge );
1733 //--------------------------------------------------------------------------------
1735 * \brief Fill in the origMesh with faces computed by MG-CADSurf in this tmp mesh
1737 //--------------------------------------------------------------------------------
1739 void FillInOrigMesh( SMESH_Mesh& origMesh,
1740 const TopoDS_Face& origFace )
1742 SMESH_MesherHelper helper( origMesh );
1743 helper.SetSubShape( origFace );
1744 helper.SetElementsOnShape( true );
1746 SMESH_MesherHelper tmpHelper( *this );
1747 tmpHelper.SetSubShape( _proxyFace );
1749 // iterate over tmp faces and copy them in origMesh
1750 const SMDS_MeshNode* nodes[27];
1751 const SMDS_MeshNode* nullNode = 0;
1753 SMDS_FaceIteratorPtr fIt = GetMeshDS()->facesIterator(/*idInceasingOrder=*/true);
1754 while ( fIt->more() )
1756 const SMDS_MeshElement* f = fIt->next();
1757 SMDS_ElemIteratorPtr nIt = f->nodesIterator();
1759 for ( ; nIt->more(); ++nbN )
1761 const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nIt->next() );
1762 TN2NMap::iterator n2nIt =
1763 _tmp2origNN.insert( _tmp2origNN.end(), make_pair( n, nullNode ));
1764 if ( !n2nIt->second ) {
1766 gp_XY uv = tmpHelper.GetNodeUV( _proxyFace, n );
1767 n2nIt->second = helper.AddNode( xyz[0], xyz[1], xyz[2], uv.X(), uv.Y() );
1769 nodes[ nbN ] = n2nIt->second;
1772 case 3: helper.AddFace( nodes[0], nodes[1], nodes[2] ); break;
1773 // case 6: helper.AddFace( nodes[0], nodes[1], nodes[2],
1774 // nodes[3], nodes[4], nodes[5]); break;
1775 case 4: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3] ); break;
1776 // case 9: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3],
1777 // nodes[4], nodes[5], nodes[6], nodes[7], nodes[8]); break;
1778 // case 8: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3],
1779 // nodes[4], nodes[5], nodes[6], nodes[7]); break;
1786 * \brief A structure holding an error description and a verbisity level
1788 struct message_cb_user_data
1790 std::string * _error;
1797 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
1798 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1799 real *duu, real *duv, real *dvv, void *user_data);
1800 status_t message_cb(message_t *msg, void *user_data);
1801 status_t interrupt_cb(integer *interrupt_status, void *user_data);
1803 //=============================================================================
1807 //=============================================================================
1809 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
1811 // Fix problem with locales
1812 Kernel_Utils::Localizer aLocalizer;
1814 this->SMESH_Algo::_progress = 1e-3; // prevent progress advancment while computing attractors
1816 bool viscousLayersMade =
1817 ( aShape.ShapeType() == TopAbs_FACE &&
1818 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( aShape ), aMesh ));
1820 if ( !viscousLayersMade )
1821 if ( !compute( aMesh, aShape, /*allowSubMeshClearing=*/true ))
1824 if ( _haveViscousLayers || viscousLayersMade )
1826 // Compute viscous layers
1828 TopTools_MapOfShape map;
1829 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1831 const TopoDS_Face& F = TopoDS::Face(face_iter.Current());
1832 if ( !map.Add( F )) continue;
1833 SMESH_ProxyMesh::Ptr viscousMesh = StdMeshers_ViscousLayers2D::Compute( aMesh, F );
1835 return false; // error in StdMeshers_ViscousLayers2D::Compute()
1837 // Compute MG-CADSurf mesh on viscous layers
1839 if ( viscousMesh->NbProxySubMeshes() > 0 )
1842 const TopoDS_Face& proxyFace = tmpMesh.makeProxyFace( viscousMesh, F );
1843 if ( !compute( tmpMesh, proxyFace, /*allowSubMeshClearing=*/false ))
1845 tmpMesh.FillInOrigMesh( aMesh, F );
1849 // Re-compute MG-CADSurf mesh on the rest faces if the mesh was cleared
1851 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1853 const TopoDS_Face& F = TopoDS::Face(face_iter.Current());
1854 SMESH_subMesh* fSM = aMesh.GetSubMesh( F );
1855 if ( fSM->IsMeshComputed() ) continue;
1857 if ( !compute( aMesh, aShape, /*allowSubMeshClearing=*/true ))
1865 //=============================================================================
1869 //=============================================================================
1871 bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh& aMesh,
1872 const TopoDS_Shape& aShape,
1873 bool allowSubMeshClearing)
1875 /* create a distene context (generic object) */
1876 status_t status = STATUS_ERROR;
1878 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1879 SMESH_MesherHelper helper( aMesh ), helperWithShape( aMesh );
1880 myHelper = theHelper = & helperWithShape;
1881 // do not call helper.IsQuadraticSubMesh() because sub-meshes
1882 // may be cleaned and helper.myTLinkNodeMap gets invalid in such a case
1883 bool haveQuadraticSubMesh = helperWithShape.IsQuadraticSubMesh( aShape );
1884 bool quadraticSubMeshAndViscousLayer = false;
1885 bool needMerge = false;
1886 typedef set< SMESHDS_SubMesh*, ShapeTypeCompare > TSubMeshSet;
1887 TSubMeshSet edgeSubmeshes;
1888 TSubMeshSet& mergeSubmeshes = edgeSubmeshes;
1890 TopTools_IndexedMapOfShape pmap, emap, fmap;
1892 TopTools_IndexedDataMapOfShapeListOfShape e2ffmap;
1893 TopExp::MapShapesAndAncestors( aShape, TopAbs_EDGE, TopAbs_FACE, e2ffmap );
1895 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1897 feclearexcept( FE_ALL_EXCEPT );
1898 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1901 context_t *ctx = context_new();
1903 /* Set the message callback in the working context */
1904 message_cb_user_data mcud;
1905 mcud._error = & this->SMESH_Algo::_comment;
1906 mcud._progress = & this->SMESH_Algo::_progress;
1908 _hypothesis ? _hypothesis->GetVerbosity() : BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
1909 context_set_message_callback(ctx, message_cb, &mcud);
1911 /* set the interruption callback */
1912 _compute_canceled = false;
1913 context_set_interrupt_callback(ctx, interrupt_cb, this);
1915 /* create the CAD object we will work on. It is associated to the context ctx. */
1916 cad_t *c = cad_new(ctx);
1917 dcad_t *dcad = dcad_new(c);
1919 // To enable multithreading
1920 cad_set_thread_safety(c, 1);
1922 /* Now fill the CAD object with data from your CAD
1923 * environement. This is the most complex part of a successfull
1929 cadsurf_session_t *css = cadsurf_session_new(ctx);
1931 // an object that correctly deletes all cadsurf objects at destruction
1932 BLSURF_Cleaner cleaner( ctx,css,c,dcad );
1934 SetParameters(_hypothesis, css, aShape);
1936 haveQuadraticSubMesh = haveQuadraticSubMesh || (_hypothesis != NULL && _hypothesis->GetQuadraticMesh());
1937 helper.SetIsQuadratic( haveQuadraticSubMesh );
1939 // To remove as soon as quadratic mesh is allowed - BEGIN
1940 // GDD: Viscous layer is not allowed with quadratic mesh
1941 if (_haveViscousLayers && haveQuadraticSubMesh ) {
1942 quadraticSubMeshAndViscousLayer = true;
1943 _haveViscousLayers = !haveQuadraticSubMesh;
1944 _comment += "Warning: Viscous layer is not possible with a quadratic mesh, it is ignored.";
1945 error(COMPERR_WARNING, _comment);
1947 // To remove as soon as quadratic mesh is allowed - END
1949 // needed to prevent the opencascade memory managmement from freeing things
1950 vector<Handle(Geom2d_Curve)> curves;
1951 vector<Handle(Geom_Surface)> surfaces;
1955 FaceId2PythonSmp.clear();
1956 EdgeId2PythonSmp.clear();
1957 VertexId2PythonSmp.clear();
1959 /****************************************************************************************
1961 *****************************************************************************************/
1963 string bad_end = "return";
1965 TopTools_IndexedMapOfShape _map;
1966 TopExp::MapShapes(aShape,TopAbs_VERTEX,_map);
1967 int ienf = _map.Extent();
1969 assert(Py_IsInitialized());
1970 PyGILState_STATE gstate;
1972 string theSizeMapStr;
1974 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1976 TopoDS_Face f = TopoDS::Face(face_iter.Current());
1978 SMESH_subMesh* fSM = aMesh.GetSubMesh( f );
1979 if ( !fSM->IsEmpty() ) continue; // skip already meshed FACE with viscous layers
1981 // make INTERNAL face oriented FORWARD (issue 0020993)
1982 if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
1983 f.Orientation(TopAbs_FORWARD);
1985 iface = fmap.Add(f);
1987 surfaces.push_back(BRep_Tool::Surface(f));
1989 /* create an object representing the face for cadsurf */
1990 /* where face_id is an integer identifying the face.
1991 * surf_function is the function that defines the surface
1992 * (For this face, it will be called by cadsurf with your_face_object_ptr
1993 * as last parameter.
1995 #if OCC_VERSION_MAJOR < 7
1996 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back());
1998 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back().get());
2001 /* by default a face has no tag (color).
2002 The following call sets it to the same value as the Geom module ID : */
2003 int faceTag = meshDS->ShapeToIndex(f);
2004 faceTag = BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis );
2005 cad_face_set_tag(fce, faceTag);
2007 /* Set face orientation (optional if you want a well oriented output mesh)*/
2008 if(f.Orientation() != TopAbs_FORWARD)
2009 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
2011 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
2013 if (HasSizeMapOnFace /*&& !use_precad*/) //22903: use_precad seems not to interfere
2015 // -----------------
2017 // -----------------
2018 faceKey = FacesWithSizeMap.FindIndex(f);
2021 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end())
2023 theSizeMapStr = FaceId2SizeMap[faceKey];
2024 // check if function ends with "return"
2025 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2027 // Expr To Python function, verification is performed at validation in GUI
2028 gstate = PyGILState_Ensure();
2029 PyObject * obj = NULL;
2030 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2032 PyObject * func = NULL;
2033 func = PyObject_GetAttrString(main_mod, "f");
2034 FaceId2PythonSmp[iface]=func;
2035 FaceId2SizeMap.erase(faceKey);
2036 PyGILState_Release(gstate);
2039 // Specific size map = Attractor
2040 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
2042 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
2043 if (attractor_iter->first == faceKey)
2045 double xyzCoords[3] = {attractor_iter->second[2],
2046 attractor_iter->second[3],
2047 attractor_iter->second[4]};
2049 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
2050 BRepClass_FaceClassifier scl(f,P,1e-7);
2051 scl.Perform(f, P, 1e-7);
2052 TopAbs_State result = scl.State();
2053 if ( result == TopAbs_OUT )
2054 MESSAGE("Point is out of face: node is not created");
2055 if ( result == TopAbs_UNKNOWN )
2056 MESSAGE("Point position on face is unknown: node is not created");
2057 if ( result == TopAbs_ON )
2058 MESSAGE("Point is on border of face: node is not created");
2059 if ( result == TopAbs_IN )
2061 // Point is inside face and not on border
2062 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
2064 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
2065 cad_point_set_tag(point_p, ienf);
2067 FaceId2AttractorCoords.erase(faceKey);
2071 // -----------------
2073 // -----------------
2074 TId2ClsAttractorVec::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
2075 if (clAttractor_iter != FaceId2ClassAttractor.end()){
2076 std::vector< BLSURFPlugin_Attractor* > & attVec = clAttractor_iter->second;
2077 for ( size_t i = 0; i < attVec.size(); ++i )
2078 if ( !attVec[i]->IsMapBuilt() ) {
2079 std::cout<<"Compute " << theNbAttractors-- << "-th attractor" <<std::endl;
2080 attVec[i]->BuildMap();
2082 FaceIndex2ClassAttractor[iface].swap( attVec );
2083 FaceId2ClassAttractor.erase(clAttractor_iter);
2085 } // if (HasSizeMapOnFace && !use_precad)
2087 // ------------------
2088 // Enforced Vertices
2089 // ------------------
2090 faceKey = FacesWithEnforcedVertices.FindIndex(f);
2091 std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
2092 if (evmIt != FaceId2EnforcedVertexCoords.end())
2094 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl = evmIt->second;
2095 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
2096 for (; evlIt != evl.end(); ++evlIt)
2098 double uvCoords[2] = { evlIt->at(0), evlIt->at(1) };
2100 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
2102 BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
2103 xyzCoords.push_back(evlIt->at(2));
2104 xyzCoords.push_back(evlIt->at(3));
2105 xyzCoords.push_back(evlIt->at(4));
2106 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(xyzCoords);
2107 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end() &&
2108 !enfCoordsIt->second.empty() )
2110 // to merge nodes of an INTERNAL vertex belonging to several faces
2111 TopoDS_Vertex v = (*enfCoordsIt->second.begin() )->vertex;
2112 if ( v.IsNull() ) v = (*enfCoordsIt->second.rbegin())->vertex;
2113 if ( !v.IsNull() && meshDS->ShapeToIndex( v ) > 0 )
2115 tag = pmap.Add( v );
2116 SMESH_subMesh* vSM = aMesh.GetSubMesh( v );
2117 vSM->ComputeStateEngine( SMESH_subMesh::COMPUTE );
2118 mergeSubmeshes.insert( vSM->GetSubMeshDS() );
2119 // //if ( tag != pmap.Extent() )
2120 // needMerge = true;
2123 if ( tag == 0 ) tag = ienf;
2124 cad_point_set_tag(point_p, tag);
2126 FaceId2EnforcedVertexCoords.erase(faceKey);
2130 /****************************************************************************************
2132 now create the edges associated to this face
2133 *****************************************************************************************/
2135 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next())
2137 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
2138 int ic = emap.FindIndex(e);
2143 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
2145 if (HasSizeMapOnEdge){
2146 edgeKey = EdgesWithSizeMap.FindIndex(e);
2147 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end())
2149 theSizeMapStr = EdgeId2SizeMap[edgeKey];
2150 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2152 // Expr To Python function, verification is performed at validation in GUI
2153 gstate = PyGILState_Ensure();
2154 PyObject * obj = NULL;
2155 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2157 PyObject * func = NULL;
2158 func = PyObject_GetAttrString(main_mod, "f");
2159 EdgeId2PythonSmp[ic]=func;
2160 EdgeId2SizeMap.erase(edgeKey);
2161 PyGILState_Release(gstate);
2164 /* data of nodes existing on the edge */
2165 StdMeshers_FaceSidePtr nodeData;
2166 SMESH_subMesh* sm = aMesh.GetSubMesh( e );
2167 if ( !sm->IsEmpty() )
2169 // SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
2170 // /*complexFirst=*/false);
2171 // while ( subsmIt->more() )
2172 // edgeSubmeshes.insert( subsmIt->next()->GetSubMeshDS() );
2173 edgeSubmeshes.insert( sm->GetSubMeshDS() );
2175 nodeData.reset( new StdMeshers_FaceSide( f, e, &aMesh, /*isForwrd = */true,
2176 /*ignoreMedium=*/haveQuadraticSubMesh));
2177 if ( nodeData->MissVertexNode() )
2178 return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
2180 const std::vector<UVPtStruct>& nodeDataVec = nodeData->GetUVPtStruct();
2181 if ( !nodeDataVec.empty() )
2183 if ( Abs( nodeDataVec[0].param - tmin ) > Abs( nodeDataVec.back().param - tmin ))
2185 nodeData->Reverse();
2186 nodeData->GetUVPtStruct(); // nodeData recomputes nodeDataVec
2188 // tmin and tmax can change in case of viscous layer on an adjacent edge
2189 tmin = nodeDataVec.front().param;
2190 tmax = nodeDataVec.back().param;
2194 cout << "---------------- Invalid nodeData" << endl;
2199 /* attach the edge to the current cadsurf face */
2200 #if OCC_VERSION_MAJOR < 7
2201 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back());
2203 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back().get());
2206 /* by default an edge has no tag (color).
2207 The following call sets it to the same value as the edge_id : */
2208 // IMP23368. Do not set tag to an EDGE shared by FACEs of a hyper-patch
2209 bool isInHyperPatch = false;
2211 std::set< int > faceTags, faceIDs;
2212 TopTools_ListIteratorOfListOfShape fIt( e2ffmap.FindFromKey( e ));
2213 for ( ; fIt.More(); fIt.Next() )
2215 int faceTag = meshDS->ShapeToIndex( fIt.Value() );
2216 if ( !faceIDs.insert( faceTag ).second )
2217 continue; // a face encounters twice for a seam edge
2218 int hpTag = BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis );
2219 if ( !faceTags.insert( hpTag ).second )
2221 isInHyperPatch = true;
2226 if ( !isInHyperPatch )
2227 cad_edge_set_tag(edg, ic);
2229 /* by default, an edge does not necessalry appear in the resulting mesh,
2230 unless the following property is set :
2232 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
2234 /* by default an edge is a boundary edge */
2235 if (e.Orientation() == TopAbs_INTERNAL)
2236 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
2238 // pass existing nodes of sub-meshes to MG-CADSurf
2241 const std::vector<UVPtStruct>& nodeDataVec = nodeData->GetUVPtStruct();
2242 const int nbNodes = nodeDataVec.size();
2244 dcad_edge_discretization_t *dedge;
2245 dcad_get_edge_discretization(dcad, edg, &dedge);
2246 dcad_edge_discretization_set_vertex_count( dedge, nbNodes );
2248 // cout << endl << " EDGE " << ic << endl;
2249 // cout << "tmin = "<<tmin << ", tmax = "<< tmax << endl;
2250 for ( int iN = 0; iN < nbNodes; ++iN )
2252 const UVPtStruct& nData = nodeDataVec[ iN ];
2253 double t = nData.param;
2254 real uv[2] = { nData.u, nData.v };
2255 SMESH_TNodeXYZ nXYZ( nData.node );
2256 // cout << "\tt = " << t
2257 // << "\t uv = ( " << uv[0] << ","<< uv[1] << " ) "
2258 // << "\t u = " << nData.param
2259 // << "\t ID = " << nData.node->GetID() << endl;
2260 dcad_edge_discretization_set_vertex_coordinates( dedge, iN+1, t, uv, nXYZ._xyz );
2262 dcad_edge_discretization_set_property(dedge, DISTENE_DCAD_PROPERTY_REQUIRED);
2265 /****************************************************************************************
2267 *****************************************************************************************/
2271 gp_Pnt2d e0 = curves.back()->Value(tmin);
2272 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
2273 Standard_Real d1=0,d2=0;
2276 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
2277 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
2281 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
2284 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
2286 *ip = pmap.FindIndex(v);
2289 // SMESH_subMesh* sm = aMesh.GetSubMesh(v);
2290 // if ( sm->IsMeshComputed() )
2291 // edgeSubmeshes.insert( sm->GetSubMeshDS() );
2294 // std::string aFileName = "fmap_vertex_";
2295 // aFileName.append(val_to_string(*ip));
2296 // aFileName.append(".brep");
2297 // BRepTools::Write(v,aFileName.c_str());
2299 if (HasSizeMapOnVertex){
2300 vertexKey = VerticesWithSizeMap.FindIndex(v);
2301 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
2302 theSizeMapStr = VertexId2SizeMap[vertexKey];
2303 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2305 // Expr To Python function, verification is performed at validation in GUI
2306 gstate = PyGILState_Ensure();
2307 PyObject * obj = NULL;
2308 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2310 PyObject * func = NULL;
2311 func = PyObject_GetAttrString(main_mod, "f");
2312 VertexId2PythonSmp[*ip]=func;
2313 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
2314 PyGILState_Release(gstate);
2319 // should not happen
2320 MESSAGE("An edge does not have 2 extremities.");
2323 // This defines the curves extremity connectivity
2324 cad_edge_set_extremities(edg, ip1, ip2);
2325 /* set the tag (color) to the same value as the extremity id : */
2326 cad_edge_set_extremities_tag(edg, ip1, ip2);
2329 cad_edge_set_extremities(edg, ip2, ip1);
2330 cad_edge_set_extremities_tag(edg, ip2, ip1);
2336 // Clear mesh from already meshed edges if possible else
2337 // remember that merge is needed
2338 TSubMeshSet::iterator smIt = edgeSubmeshes.begin();
2339 for ( ; smIt != edgeSubmeshes.end(); ++smIt ) // loop on already meshed EDGEs
2341 SMESHDS_SubMesh* smDS = *smIt;
2342 if ( !smDS ) continue;
2343 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2346 const SMDS_MeshNode* n = nIt->next();
2347 if ( n->NbInverseElements( SMDSAbs_Face ) > 0 )
2349 needMerge = true; // to correctly sew with viscous mesh
2350 // add existing medium nodes to helper
2351 if ( aMesh.NbEdges( ORDER_QUADRATIC ) > 0 )
2353 SMDS_ElemIteratorPtr edgeIt = smDS->GetElements();
2354 while ( edgeIt->more() )
2355 helper.AddTLinks( static_cast<const SMDS_MeshEdge*>(edgeIt->next()));
2360 if ( allowSubMeshClearing )
2362 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
2363 while ( eIt->more() ) meshDS->RemoveFreeElement( eIt->next(), 0 );
2364 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2365 while ( nIt->more() ) meshDS->RemoveFreeNode( nIt->next(), 0 );
2374 ///////////////////////
2376 ///////////////////////
2378 if (! _preCadFacesIDsPeriodicityVector.empty())
2380 for (std::size_t i=0; i < _preCadFacesIDsPeriodicityVector.size(); i++){
2381 std::vector<int> theFace1_ids = _preCadFacesIDsPeriodicityVector[i].shape1IDs;
2382 std::vector<int> theFace2_ids = _preCadFacesIDsPeriodicityVector[i].shape2IDs;
2383 int* theFace1_ids_c = &theFace1_ids[0];
2384 int* theFace2_ids_c = &theFace2_ids[0];
2385 std::ostringstream o;
2386 o << "_preCadFacesIDsPeriodicityVector[" << i << "] = [";
2387 for (std::size_t j=0; j < theFace1_ids.size(); j++)
2388 o << theFace1_ids[j] << ", ";
2390 for (std::size_t j=0; j < theFace2_ids.size(); j++)
2391 o << theFace2_ids[j] << ", ";
2393 // if ( _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
2394 // cout << o.str() << endl;
2395 if (_preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords.empty())
2397 // If no source points, call periodicity without transformation function
2398 meshgems_cad_periodicity_transformation_t periodicity_transformation = NULL;
2399 status = cad_add_face_multiple_periodicity_with_transformation_function(c, theFace1_ids_c, theFace1_ids.size(),
2400 theFace2_ids_c, theFace2_ids.size(), periodicity_transformation, NULL);
2401 if(status != STATUS_OK)
2402 cout << "cad_add_face_multiple_periodicity_with_transformation_function failed with error code " << status << "\n";
2406 // get the transformation vertices
2407 double* theSourceVerticesCoords_c = &_preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords[0];
2408 double* theTargetVerticesCoords_c = &_preCadFacesIDsPeriodicityVector[i].theTargetVerticesCoords[0];
2409 int nbSourceVertices = _preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords.size()/3;
2410 int nbTargetVertices = _preCadFacesIDsPeriodicityVector[i].theTargetVerticesCoords.size()/3;
2412 status = cad_add_face_multiple_periodicity_with_transformation_function_by_points(c, theFace1_ids_c, theFace1_ids.size(),
2413 theFace2_ids_c, theFace2_ids.size(), theSourceVerticesCoords_c, nbSourceVertices, theTargetVerticesCoords_c, nbTargetVertices);
2414 if(status != STATUS_OK)
2415 cout << "cad_add_face_multiple_periodicity_with_transformation_function_by_points failed with error code " << status << "\n";
2420 if (! _preCadEdgesIDsPeriodicityVector.empty())
2422 for (std::size_t i=0; i < _preCadEdgesIDsPeriodicityVector.size(); i++){
2423 std::vector<int> theEdge1_ids = _preCadEdgesIDsPeriodicityVector[i].shape1IDs;
2424 std::vector<int> theEdge2_ids = _preCadEdgesIDsPeriodicityVector[i].shape2IDs;
2425 // Use the address of the first element of the vector to initialize the array
2426 int* theEdge1_ids_c = &theEdge1_ids[0];
2427 int* theEdge2_ids_c = &theEdge2_ids[0];
2429 std::ostringstream o;
2430 o << "_preCadEdgesIDsPeriodicityVector[" << i << "] = [";
2431 for (std::size_t j=0; j < theEdge1_ids.size(); j++)
2432 o << theEdge1_ids[j] << ", ";
2434 for (std::size_t j=0; j < theEdge2_ids.size(); j++)
2435 o << theEdge2_ids[j] << ", ";
2437 // if ( _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
2438 // cout << o.str() << endl;
2440 if (_preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords.empty())
2442 // If no source points, call periodicity without transformation function
2443 meshgems_cad_periodicity_transformation_t periodicity_transformation = NULL;
2444 status = cad_add_edge_multiple_periodicity_with_transformation_function(c, theEdge1_ids_c, theEdge1_ids.size(),
2445 theEdge2_ids_c, theEdge2_ids.size(), periodicity_transformation, NULL);
2446 if(status != STATUS_OK)
2447 cout << "cad_add_edge_multiple_periodicity_with_transformation_function failed with error code " << status << "\n";
2451 // get the transformation vertices
2452 double* theSourceVerticesCoords_c = &_preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords[0];
2453 double* theTargetVerticesCoords_c = &_preCadEdgesIDsPeriodicityVector[i].theTargetVerticesCoords[0];
2454 int nbSourceVertices = _preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords.size()/3;
2455 int nbTargetVertices = _preCadEdgesIDsPeriodicityVector[i].theTargetVerticesCoords.size()/3;
2457 status = cad_add_edge_multiple_periodicity_with_transformation_function_by_points(c, theEdge1_ids_c, theEdge1_ids.size(),
2458 theEdge2_ids_c, theEdge2_ids.size(), theSourceVerticesCoords_c, nbSourceVertices, theTargetVerticesCoords_c, nbTargetVertices);
2459 if(status != STATUS_OK)
2460 cout << "cad_add_edge_multiple_periodicity_with_transformation_function_by_points failed with error code " << status << "\n";
2466 // TODO: be able to use a mesh in input.
2467 // See imsh usage in Products/templates/mg-cadsurf_template_common.cpp
2468 // => cadsurf_set_mesh
2470 // Use the original dcad
2471 cadsurf_set_dcad(css, dcad);
2473 // Use the original cad
2474 cadsurf_set_cad(css, c);
2476 std::cout << std::endl;
2477 std::cout << "Beginning of Surface Mesh generation" << std::endl;
2478 std::cout << std::endl;
2483 status = cadsurf_compute_mesh(css);
2486 catch ( std::exception& exc ) {
2487 _comment += exc.what();
2489 catch (Standard_Failure& ex) {
2490 _comment += ex.DynamicType()->Name();
2491 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
2493 _comment += ex.GetMessageString();
2497 if ( _comment.empty() )
2498 _comment = "Exception in cadsurf_compute_mesh()";
2501 std::cout << std::endl;
2502 std::cout << "End of Surface Mesh generation" << std::endl;
2503 std::cout << std::endl;
2506 cadsurf_get_mesh(css, &msh);
2508 /* release the mesh object */
2509 cadsurf_regain_mesh(css, msh);
2510 return error(_comment);
2513 std::string GMFFileName = BLSURFPlugin_Hypothesis::GetDefaultGMFFile();
2515 GMFFileName = _hypothesis->GetGMFFile();
2516 if (GMFFileName != "") {
2517 bool asciiFound = (GMFFileName.find(".mesh", GMFFileName.length()-5) != std::string::npos);
2518 bool binaryFound = (GMFFileName.find(".meshb",GMFFileName.length()-6) != std::string::npos);
2519 if (!asciiFound && !binaryFound)
2520 GMFFileName.append(".mesh");
2521 mesh_write_mesh(msh, GMFFileName.c_str());
2524 /* retrieve mesh data (see meshgems/mesh.h) */
2525 integer nv, ne, nt, nq, vtx[4], tag, nb_tag;
2526 integer *evedg, *evtri, *evquad, *tags_buff, type;
2529 mesh_get_vertex_count(msh, &nv);
2530 mesh_get_edge_count(msh, &ne);
2531 mesh_get_triangle_count(msh, &nt);
2532 mesh_get_quadrangle_count(msh, &nq);
2534 evedg = (integer *)mesh_calloc_generic_buffer(msh);
2535 evtri = (integer *)mesh_calloc_generic_buffer(msh);
2536 evquad = (integer *)mesh_calloc_generic_buffer(msh);
2537 tags_buff = (integer*)mesh_calloc_generic_buffer(msh);
2539 std::vector<const SMDS_MeshNode*> nodes(nv+1);
2540 std::vector<bool> tags(nv+1);
2542 /* enumerated vertices */
2543 for(int iv=1;iv<=nv;iv++) {
2544 mesh_get_vertex_coordinates(msh, iv, xyz);
2545 mesh_get_vertex_tag(msh, iv, &tag);
2546 // Issue 0020656. Use vertex coordinates
2548 if ( tag > 0 && tag <= pmap.Extent() ) {
2549 TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
2550 double tol = BRep_Tool::Tolerance( v );
2551 gp_Pnt p = BRep_Tool::Pnt( v );
2552 if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 2*tol))
2553 xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
2555 tag = 0; // enforced or attracted vertex
2556 nodes[iv] = SMESH_Algo::VertexNode( v, meshDS );
2559 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
2561 // Create group of enforced vertices if requested
2562 BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
2564 projVertex.push_back((double)xyz[0]);
2565 projVertex.push_back((double)xyz[1]);
2566 projVertex.push_back((double)xyz[2]);
2567 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
2568 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end())
2570 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
2571 BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
2572 for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
2573 currentEnfVertex = (*enfListIt);
2574 if (currentEnfVertex->grpName != "") {
2575 bool groupDone = false;
2576 SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
2577 while (grIt->more()) {
2578 SMESH_Group * group = grIt->next();
2579 if ( !group ) continue;
2580 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
2581 if ( !groupDS ) continue;
2582 if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
2583 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
2584 aGroupDS->SMDSGroup().Add(nodes[iv]);
2585 // How can I inform the hypothesis ?
2586 // _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
2594 SMESH_Group* aGroup = aMesh.AddGroup(SMDSAbs_Node, currentEnfVertex->grpName.c_str(), groupId);
2595 aGroup->SetName( currentEnfVertex->grpName.c_str() );
2596 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
2597 aGroupDS->SMDSGroup().Add(nodes[iv]);
2601 throw SALOME_Exception(LOCALIZED("An enforced vertex node was not added to a group"));
2604 MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
2608 // internal points are tagged to zero
2609 if(tag > 0 && tag <= pmap.Extent() ){
2610 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
2617 /* enumerate edges */
2618 for(int it=1;it<=ne;it++) {
2620 mesh_get_edge_vertices(msh, it, vtx);
2621 mesh_get_edge_extra_vertices(msh, it, &type, evedg);
2622 mesh_get_edge_tag(msh, it, &tag);
2624 // If PreCAD performed some cleaning operations (remove tiny edges,
2625 // merge edges ...) an output tag can indeed represent several original tags.
2626 // Get the initial tags corresponding to the output tag and redefine the tag as
2627 // the last of the two initial tags (else the output tag is out of emap and hasn't any meaning)
2628 mesh_get_composite_tag_definition(msh, tag, &nb_tag, tags_buff);
2630 tag=tags_buff[nb_tag-1];
2631 if ( tag < 1 || tag > emap.Extent() )
2633 std::cerr << "MG-CADSurf BUG:::: Edge tag " << tag
2634 << " does not point to a CAD edge (nb edges " << emap.Extent() << ")" << std::endl;
2638 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
2639 tags[vtx[0]] = false;
2642 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
2643 tags[vtx[1]] = false;
2645 if (type == MESHGEMS_MESH_ELEMENT_TYPE_EDGE3) {
2647 if (tags[evedg[0]]) {
2648 Set_NodeOnEdge(meshDS, nodes[evedg[0]], emap(tag));
2649 tags[evedg[0]] = false;
2651 edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]], nodes[evedg[0]]);
2654 edg = helper.AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
2656 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
2659 /* enumerate triangles */
2660 for(int it=1;it<=nt;it++) {
2662 mesh_get_triangle_vertices(msh, it, vtx);
2663 mesh_get_triangle_extra_vertices(msh, it, &type, evtri);
2664 mesh_get_triangle_tag(msh, it, &tag);
2666 meshDS->SetNodeOnFace(nodes[vtx[0]], tag);
2667 tags[vtx[0]] = false;
2670 meshDS->SetNodeOnFace(nodes[vtx[1]], tag);
2671 tags[vtx[1]] = false;
2674 meshDS->SetNodeOnFace(nodes[vtx[2]], tag);
2675 tags[vtx[2]] = false;
2677 if (type == MESHGEMS_MESH_ELEMENT_TYPE_TRIA6) {
2678 // QUADRATIC TRIANGLE
2679 if (tags[evtri[0]]) {
2680 meshDS->SetNodeOnFace(nodes[evtri[0]], tag);
2681 tags[evtri[0]] = false;
2683 if (tags[evtri[1]]) {
2684 meshDS->SetNodeOnFace(nodes[evtri[1]], tag);
2685 tags[evtri[1]] = false;
2687 if (tags[evtri[2]]) {
2688 meshDS->SetNodeOnFace(nodes[evtri[2]], tag);
2689 tags[evtri[2]] = false;
2691 tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]],
2692 nodes[evtri[0]], nodes[evtri[1]], nodes[evtri[2]]);
2695 tri = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
2697 meshDS->SetMeshElementOnShape(tri, tag);
2700 /* enumerate quadrangles */
2701 for(int it=1;it<=nq;it++) {
2702 SMDS_MeshFace* quad;
2703 mesh_get_quadrangle_vertices(msh, it, vtx);
2704 mesh_get_quadrangle_extra_vertices(msh, it, &type, evquad);
2705 mesh_get_quadrangle_tag(msh, it, &tag);
2707 meshDS->SetNodeOnFace(nodes[vtx[0]], tag);
2708 tags[vtx[0]] = false;
2711 meshDS->SetNodeOnFace(nodes[vtx[1]], tag);
2712 tags[vtx[1]] = false;
2715 meshDS->SetNodeOnFace(nodes[vtx[2]], tag);
2716 tags[vtx[2]] = false;
2719 meshDS->SetNodeOnFace(nodes[vtx[3]], tag);
2720 tags[vtx[3]] = false;
2722 if (type == MESHGEMS_MESH_ELEMENT_TYPE_QUAD9) {
2723 // QUADRATIC QUADRANGLE
2724 std::cout << "This is a quadratic quadrangle" << std::endl;
2725 if (tags[evquad[0]]) {
2726 meshDS->SetNodeOnFace(nodes[evquad[0]], tag);
2727 tags[evquad[0]] = false;
2729 if (tags[evquad[1]]) {
2730 meshDS->SetNodeOnFace(nodes[evquad[1]], tag);
2731 tags[evquad[1]] = false;
2733 if (tags[evquad[2]]) {
2734 meshDS->SetNodeOnFace(nodes[evquad[2]], tag);
2735 tags[evquad[2]] = false;
2737 if (tags[evquad[3]]) {
2738 meshDS->SetNodeOnFace(nodes[evquad[3]], tag);
2739 tags[evquad[3]] = false;
2741 if (tags[evquad[4]]) {
2742 meshDS->SetNodeOnFace(nodes[evquad[4]], tag);
2743 tags[evquad[4]] = false;
2745 quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]],
2746 nodes[evquad[0]], nodes[evquad[1]], nodes[evquad[2]], nodes[evquad[3]],
2750 quad = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
2752 meshDS->SetMeshElementOnShape(quad, tag);
2755 /* release the mesh object, the rest is released by cleaner */
2756 cadsurf_regain_mesh(css, msh);
2759 // Remove free nodes that can appear e.g. if "remove tiny edges"(IPAL53235)
2760 for(int iv=1;iv<=nv;iv++)
2761 if ( nodes[iv] && nodes[iv]->NbInverseElements() == 0 )
2762 meshDS->RemoveFreeNode( nodes[iv], 0, /*fromGroups=*/false );
2765 if ( needMerge ) // sew mesh computed by MG-CADSurf with pre-existing mesh
2767 SMESH_MeshEditor editor( &aMesh );
2768 SMESH_MeshEditor::TListOfListOfNodes nodeGroupsToMerge;
2769 TIDSortedElemSet segementsOnEdge;
2770 TSubMeshSet::iterator smIt;
2771 SMESHDS_SubMesh* smDS;
2773 // merge nodes on EDGE's with ones computed by MG-CADSurf
2774 for ( smIt = mergeSubmeshes.begin(); smIt != mergeSubmeshes.end(); ++smIt )
2776 if (! (smDS = *smIt) ) continue;
2777 getNodeGroupsToMerge( smDS, meshDS->IndexToShape((*smIt)->GetID()), nodeGroupsToMerge );
2779 SMDS_ElemIteratorPtr segIt = smDS->GetElements();
2780 while ( segIt->more() )
2781 segementsOnEdge.insert( segIt->next() );
2784 editor.MergeNodes( nodeGroupsToMerge );
2787 SMESH_MeshEditor::TListOfListOfElementsID equalSegments;
2788 editor.FindEqualElements( segementsOnEdge, equalSegments );
2789 editor.MergeElements( equalSegments );
2791 // remove excess segments created on the boundary of viscous layers
2792 const SMDS_TypeOfPosition onFace = SMDS_TOP_FACE;
2793 for ( int i = 1; i <= emap.Extent(); ++i )
2795 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( emap( i )))
2797 SMDS_ElemIteratorPtr segIt = smDS->GetElements();
2798 while ( segIt->more() )
2800 const SMDS_MeshElement* seg = segIt->next();
2801 if ( seg->GetNode(0)->GetPosition()->GetTypeOfPosition() == onFace ||
2802 seg->GetNode(1)->GetPosition()->GetTypeOfPosition() == onFace )
2803 meshDS->RemoveFreeElement( seg, smDS );
2810 // SetIsAlwaysComputed( true ) to sub-meshes of EDGEs w/o mesh
2811 for (int i = 1; i <= emap.Extent(); i++)
2812 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
2813 sm->SetIsAlwaysComputed( true );
2814 for (int i = 1; i <= pmap.Extent(); i++)
2815 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( pmap( i )))
2816 if ( !sm->IsMeshComputed() )
2817 sm->SetIsAlwaysComputed( true );
2819 // Set error to FACE's w/o elements
2820 SMESH_ComputeErrorName err = COMPERR_ALGO_FAILED;
2821 if ( _comment.empty() && status == STATUS_OK )
2823 err = COMPERR_WARNING;
2824 _comment = "No mesh elements assigned to a face";
2826 bool badFaceFound = false;
2827 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
2829 TopoDS_Face f = TopoDS::Face(face_iter.Current());
2830 SMESH_subMesh* sm = aMesh.GetSubMesh( f );
2831 if ( !sm->GetSubMeshDS() || sm->GetSubMeshDS()->NbElements() == 0 )
2833 int faceTag = sm->GetId();
2834 if ( faceTag != BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis ))
2836 // triangles are assigned to the first face of hyper-patch
2837 sm->SetIsAlwaysComputed( true );
2841 sm->GetComputeError().reset( new SMESH_ComputeError( err, _comment, this ));
2842 badFaceFound = true;
2846 if ( err == COMPERR_WARNING )
2850 if ( status != STATUS_OK && !badFaceFound ) {
2854 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
2856 if ( oldFEFlags > 0 )
2857 feenableexcept( oldFEFlags );
2858 feclearexcept( FE_ALL_EXCEPT );
2862 std::cout << "FacesWithSizeMap" << std::endl;
2863 FacesWithSizeMap.Statistics(std::cout);
2864 std::cout << "EdgesWithSizeMap" << std::endl;
2865 EdgesWithSizeMap.Statistics(std::cout);
2866 std::cout << "VerticesWithSizeMap" << std::endl;
2867 VerticesWithSizeMap.Statistics(std::cout);
2868 std::cout << "FacesWithEnforcedVertices" << std::endl;
2869 FacesWithEnforcedVertices.Statistics(std::cout);
2872 return ( status == STATUS_OK && !quadraticSubMeshAndViscousLayer );
2875 //================================================================================
2877 * \brief Compute a mesh basing on discrete CAD description
2879 //================================================================================
2881 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
2883 if ( aMesh.NbFaces() == 0 )
2884 return error( COMPERR_BAD_INPUT_MESH, "2D elements are missing" );
2886 context_t *ctx = context_new();
2887 if (!ctx) return error("Pb in context_new()");
2889 BLSURF_Cleaner cleaner( ctx );
2891 message_cb_user_data mcud;
2892 mcud._error = & this->SMESH_Algo::_comment;
2893 mcud._progress = & this->SMESH_Algo::_progress;
2895 _hypothesis ? _hypothesis->GetVerbosity() : BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
2896 meshgems_status_t ret = context_set_message_callback(ctx, message_cb, &mcud);
2897 if (ret != STATUS_OK) return error("Pb. in context_set_message_callback() ");
2899 cadsurf_session_t * css = cadsurf_session_new(ctx);
2900 if(!css) return error( "Pb. in cadsurf_session_new() " );
2904 // Fill an input mesh
2906 mesh_t * msh = meshgems_mesh_new_in_memory( ctx );
2907 if ( !msh ) return error("Pb. in meshgems_mesh_new_in_memory()");
2909 // mark nodes used by 2D elements
2910 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
2911 SMDS_NodeIteratorPtr nodeIt = meshDS->nodesIterator();
2912 while ( nodeIt->more() )
2914 const SMDS_MeshNode* n = nodeIt->next();
2915 n->setIsMarked( n->NbInverseElements( SMDSAbs_Face ));
2917 meshgems_mesh_set_vertex_count( msh, meshDS->NbNodes() );
2919 // set node coordinates
2920 if ( meshDS->NbNodes() != meshDS->MaxNodeID() )
2922 meshDS->compactMesh();
2924 SMESH_TNodeXYZ nXYZ;
2925 nodeIt = meshDS->nodesIterator();
2927 for ( i = 1; nodeIt->more(); ++i )
2929 nXYZ.Set( nodeIt->next() );
2930 meshgems_mesh_set_vertex_coordinates( msh, i, nXYZ._xyz );
2933 // set nodes of faces
2934 meshgems_mesh_set_triangle_count ( msh, meshDS->GetMeshInfo().NbTriangles() );
2935 meshgems_mesh_set_quadrangle_count( msh, meshDS->GetMeshInfo().NbQuadrangles() );
2936 meshgems_integer nodeIDs[4];
2937 meshgems_integer iT = 1, iQ = 1;
2938 SMDS_FaceIteratorPtr faceIt = meshDS->facesIterator();
2939 while ( faceIt->more() )
2941 const SMDS_MeshElement* face = faceIt->next();
2942 meshgems_integer nbNodes = face->NbCornerNodes();
2943 if ( nbNodes > 4 || face->IsPoly() ) continue;
2945 for ( i = 0; i < nbNodes; ++i )
2946 nodeIDs[i] = face->GetNode( i )->GetID();
2948 meshgems_mesh_set_triangle_vertices ( msh, iT++, nodeIDs );
2950 meshgems_mesh_set_quadrangle_vertices( msh, iQ++, nodeIDs );
2953 ret = cadsurf_set_mesh(css, msh);
2954 if ( ret != STATUS_OK ) return error("Pb in cadsurf_set_mesh()");
2959 SetParameters(_hypothesis, css, aMesh.GetShapeToMesh() );
2961 ret = cadsurf_compute_mesh(css);
2962 if ( ret != STATUS_OK ) return false;
2965 cadsurf_get_mesh(css, &omsh);
2966 if ( !omsh ) return error( "Pb. in cadsurf_get_mesh()" );
2969 // Update SALOME mesh
2971 // remove quadrangles and triangles
2972 for ( faceIt = meshDS->facesIterator(); faceIt->more(); )
2974 const SMDS_MeshElement* face = faceIt->next();
2975 if ( !face->IsPoly() )
2976 meshDS->RemoveFreeElement( face, /*sm=*/0, /*fromGroups=*/true );
2978 // remove edges that bound the just removed faces
2979 for ( SMDS_EdgeIteratorPtr edgeIt = meshDS->edgesIterator(); edgeIt->more(); )
2981 const SMDS_MeshElement* edge = edgeIt->next();
2982 const SMDS_MeshNode* n0 = edge->GetNode(0);
2983 const SMDS_MeshNode* n1 = edge->GetNode(1);
2984 if ( n0->isMarked() &&
2986 n0->NbInverseElements( SMDSAbs_Volume ) == 0 &&
2987 n1->NbInverseElements( SMDSAbs_Volume ) == 0 )
2988 meshDS->RemoveFreeElement( edge, /*sm=*/0, /*fromGroups=*/true );
2990 // remove nodes that just became free
2991 for ( nodeIt = meshDS->nodesIterator(); nodeIt->more(); )
2993 const SMDS_MeshNode* n = nodeIt->next();
2994 if ( n->isMarked() && n->NbInverseElements() == 0 )
2995 meshDS->RemoveFreeNode( n, /*sm=*/0, /*fromGroups=*/true );
2999 meshgems_integer nbvtx = 0, nodeID;
3000 meshgems_mesh_get_vertex_count( omsh, &nbvtx );
3001 meshgems_real xyz[3];
3002 for ( i = 1; i <= nbvtx; ++i )
3004 meshgems_mesh_get_vertex_coordinates( omsh, i, xyz );
3005 SMDS_MeshNode* n = meshDS->AddNode( xyz[0], xyz[1], xyz[2] );
3006 nodeID = n->GetID();
3007 meshgems_mesh_set_vertex_tag( omsh, i, &nodeID ); // save mapping of IDs in MG and SALOME meshes
3011 meshgems_integer nbtri = 0;
3012 meshgems_mesh_get_triangle_count( omsh, &nbtri );
3013 const SMDS_MeshNode* nodes[3];
3014 for ( i = 1; i <= nbtri; ++i )
3016 meshgems_mesh_get_triangle_vertices( omsh, i, nodeIDs );
3017 for ( int j = 0; j < 3; ++j )
3019 meshgems_mesh_get_vertex_tag( omsh, nodeIDs[j], &nodeID );
3020 nodes[j] = meshDS->FindNode( nodeID );
3022 meshDS->AddFace( nodes[0], nodes[1], nodes[2] );
3025 cadsurf_regain_mesh(css, omsh);
3027 // as we don't assign the new triangles to a shape (the pseudo-shape),
3028 // we mark the shape as always computed to avoid the error messages
3029 // that no elements assigned to the shape
3030 aMesh.GetSubMesh( aHelper->GetSubShape() )->SetIsAlwaysComputed( true );
3035 //================================================================================
3037 * \brief Terminates computation
3039 //================================================================================
3041 void BLSURFPlugin_BLSURF::CancelCompute()
3043 _compute_canceled = true;
3046 //=============================================================================
3050 //=============================================================================
3052 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS,
3053 const SMDS_MeshNode* node,
3054 const TopoDS_Shape& ed)
3056 const TopoDS_Edge edge = TopoDS::Edge(ed);
3058 gp_Pnt pnt(node->X(), node->Y(), node->Z());
3060 Standard_Real p0 = 0.0;
3061 Standard_Real p1 = 1.0;
3062 TopLoc_Location loc;
3063 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
3064 if ( curve.IsNull() )
3066 // issue 22499. Node at a sphere apex
3067 meshDS->SetNodeOnEdge(node, edge, p0);
3071 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
3072 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
3075 if ( proj.NbPoints() > 0 )
3077 pa = (double)proj.LowerDistanceParameter();
3078 // Issue 0020656. Move node if it is too far from edge
3079 gp_Pnt curve_pnt = curve->Value( pa );
3080 double dist2 = pnt.SquareDistance( curve_pnt );
3081 double tol = BRep_Tool::Tolerance( edge );
3082 if ( 1e-14 < dist2 && dist2 <= 1000*tol ) // large enough and within tolerance
3084 curve_pnt.Transform( loc );
3085 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
3089 meshDS->SetNodeOnEdge(node, edge, pa);
3092 /* Curve definition function See cad_curv_t in file meshgems/cad.h for
3094 * NOTE : if when your CAD systems evaluates second
3095 * order derivatives it also computes first order derivatives and
3096 * function evaluation, you can optimize this example by making only
3097 * one CAD call and filling the necessary uv, dt, dtt arrays.
3099 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
3101 /* t is given. It contains the t (time) 1D parametric coordintaes
3102 of the point PreCAD/MG-CADSurf is querying on the curve */
3104 /* user_data identifies the edge PreCAD/MG-CADSurf is querying
3105 * (see cad_edge_new later in this example) */
3106 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
3109 /* MG-CADSurf is querying the function evaluation */
3112 uv[0]=P.X(); uv[1]=P.Y();
3116 /* query for the first order derivatives */
3119 dt[0]=V1.X(); dt[1]=V1.Y();
3123 /* query for the second order derivatives */
3126 dtt[0]=V2.X(); dtt[1]=V2.Y();
3132 /* Surface definition function.
3133 * See cad_surf_t in file meshgems/cad.h for more information.
3134 * NOTE : if when your CAD systems evaluates second order derivatives it also
3135 * computes first order derivatives and function evaluation, you can optimize
3136 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
3139 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
3140 real *duu, real *duv, real *dvv, void *user_data)
3142 /* uv[2] is given. It contains the u,v coordinates of the point
3143 * PreCAD/MG-CADSurf is querying on the surface */
3145 /* user_data identifies the face PreCAD/MG-CADSurf is querying (see
3146 * cad_face_new later in this example)*/
3147 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
3151 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
3152 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
3159 geometry->D1(uv[0],uv[1],P,D1U,D1V);
3160 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
3161 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
3164 if(duu && duv && dvv){
3168 gp_Vec D2U,D2V,D2UV;
3170 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
3171 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
3172 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
3173 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
3180 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
3182 TId2ClsAttractorVec::iterator f2attVec;
3183 if (FaceId2PythonSmp.count(face_id) != 0){
3184 assert(Py_IsInitialized());
3185 PyGILState_STATE gstate;
3186 gstate = PyGILState_Ensure();
3187 PyObject* pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
3189 if ( pyresult != NULL) {
3190 result = PyFloat_AsDouble(pyresult);
3191 Py_DECREF(pyresult);
3196 string err_description="";
3197 PyObject* new_stderr = newPyStdOut(err_description);
3198 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3199 Py_INCREF(old_stderr);
3200 PySys_SetObject((char*)"stderr", new_stderr);
3202 PySys_SetObject((char*)"stderr", old_stderr);
3203 Py_DECREF(new_stderr);
3204 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
3205 result = *((real*)user_data);
3208 PyGILState_Release(gstate);
3210 else if (( f2attVec = FaceIndex2ClassAttractor.find(face_id)) != FaceIndex2ClassAttractor.end() && !f2attVec->second.empty())
3214 std::vector< BLSURFPlugin_Attractor* > & attVec = f2attVec->second;
3215 for ( size_t i = 0; i < attVec.size(); ++i )
3217 //result += attVec[i]->GetSize(uv[0],uv[1]);
3218 result = Min( result, attVec[i]->GetSize(uv[0],uv[1]));
3220 //*size = result / attVec.size(); // mean of sizes defined by all attractors
3224 *size = *((real*)user_data);
3226 // std::cout << "Size_on_surface sur la face " << face_id << " donne une size de: " << *size << std::endl;
3230 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
3232 if (EdgeId2PythonSmp.count(edge_id) != 0){
3233 assert(Py_IsInitialized());
3234 PyGILState_STATE gstate;
3235 gstate = PyGILState_Ensure();
3236 PyObject* pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
3238 if ( pyresult != NULL) {
3239 result = PyFloat_AsDouble(pyresult);
3240 Py_DECREF(pyresult);
3245 string err_description="";
3246 PyObject* new_stderr = newPyStdOut(err_description);
3247 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3248 Py_INCREF(old_stderr);
3249 PySys_SetObject((char*)"stderr", new_stderr);
3251 PySys_SetObject((char*)"stderr", old_stderr);
3252 Py_DECREF(new_stderr);
3253 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
3254 result = *((real*)user_data);
3257 PyGILState_Release(gstate);
3260 *size = *((real*)user_data);
3265 status_t size_on_vertex(integer point_id, real *size, void *user_data)
3267 if (VertexId2PythonSmp.count(point_id) != 0){
3268 assert(Py_IsInitialized());
3269 PyGILState_STATE gstate;
3270 gstate = PyGILState_Ensure();
3271 PyObject* pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
3273 if ( pyresult != NULL) {
3274 result = PyFloat_AsDouble(pyresult);
3275 Py_DECREF(pyresult);
3280 string err_description="";
3281 PyObject* new_stderr = newPyStdOut(err_description);
3282 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3283 Py_INCREF(old_stderr);
3284 PySys_SetObject((char*)"stderr", new_stderr);
3286 PySys_SetObject((char*)"stderr", old_stderr);
3287 Py_DECREF(new_stderr);
3288 MESSAGE("Can't evaluate f()" << " error is " << err_description);
3289 result = *((real*)user_data);
3292 PyGILState_Release(gstate);
3295 *size = *((real*)user_data);
3301 * The following function will be called for PreCAD/MG-CADSurf message
3302 * printing. See context_set_message_callback (later in this
3303 * template) for how to set user_data.
3305 status_t message_cb(message_t *msg, void *user_data)
3307 integer errnumber = 0;
3309 message_get_number(msg, &errnumber);
3310 message_get_description(msg, &desc);
3312 message_cb_user_data * mcud = (message_cb_user_data*)user_data;
3313 // Get all the error message and some warning messages related to license and periodicity
3314 if ( errnumber < 0 ||
3315 err.find("license" ) != string::npos ||
3316 err.find("periodicity") != string::npos )
3318 // remove ^A from the tail
3319 int len = strlen( desc );
3320 while (len > 0 && desc[len-1] != '\n')
3322 mcud->_error->append( desc, len );
3325 if ( errnumber == 3009001 )
3326 * mcud->_progress = atof( desc + 11 ) / 100.;
3327 if ( mcud->_verbosity > 0 )
3328 std::cout << desc << std::endl;
3333 /* This is the interrupt callback. PreCAD/MG-CADSurf will call this
3334 * function regularily. See the file meshgems/interrupt.h
3336 status_t interrupt_cb(integer *interrupt_status, void *user_data)
3338 integer you_want_to_continue = 1;
3339 BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
3340 you_want_to_continue = !tmp->computeCanceled();
3342 if(you_want_to_continue)
3344 *interrupt_status = INTERRUPT_CONTINUE;
3347 else /* you want to stop MG-CADSurf */
3349 *interrupt_status = INTERRUPT_STOP;
3350 return STATUS_ERROR;
3354 //=============================================================================
3358 //=============================================================================
3359 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
3360 const TopoDS_Shape& aShape,
3361 MapShapeNbElems& aResMap)
3363 double diagonal = aMesh.GetShapeDiagonalSize();
3364 double bbSegmentation = _gen->GetBoundaryBoxSegmentation();
3365 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
3366 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize(diagonal, bbSegmentation);
3367 bool _phySizeRel = BLSURFPlugin_Hypothesis::GetDefaultPhySizeRel();
3368 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
3369 double _angleMesh = BLSURFPlugin_Hypothesis::GetDefaultAngleMesh();
3370 bool _quadAllowed = BLSURFPlugin_Hypothesis::GetDefaultQuadAllowed();
3372 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
3373 _phySizeRel = _hypothesis->IsPhySizeRel();
3374 if ( _hypothesis->GetPhySize() > 0)
3375 _phySize = _phySizeRel ? diagonal*_hypothesis->GetPhySize() : _hypothesis->GetPhySize();
3376 //_geometricMesh = (int) hyp->GetGeometricMesh();
3377 if (_hypothesis->GetAngleMesh() > 0)
3378 _angleMesh = _hypothesis->GetAngleMesh();
3379 _quadAllowed = _hypothesis->GetQuadAllowed();
3381 //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
3382 // GetDefaultPhySize() sometimes leads to computation failure
3383 _phySize = aMesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
3386 bool IsQuadratic = _quadraticMesh;
3391 TopTools_DataMapOfShapeInteger EdgesMap;
3392 double fullLen = 0.0;
3393 double fullNbSeg = 0;
3394 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
3395 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3396 if( EdgesMap.IsBound(E) )
3398 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
3399 double aLen = SMESH_Algo::EdgeLength(E);
3402 if(_physicalMesh==1) {
3403 nb1d = (int)( aLen/_phySize + 1 );
3408 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
3409 double fullAng = 0.0;
3410 double dp = (l-f)/200;
3415 for(int j=2; j<=200; j++) {
3418 fullAng += fabs(V1.Angle(V2));
3422 nb1d = (int)( fullAng/_angleMesh + 1 );
3425 std::vector<int> aVec(SMDSEntity_Last);
3426 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
3427 if( IsQuadratic > 0 ) {
3428 aVec[SMDSEntity_Node] = 2*nb1d - 1;
3429 aVec[SMDSEntity_Quad_Edge] = nb1d;
3432 aVec[SMDSEntity_Node] = nb1d - 1;
3433 aVec[SMDSEntity_Edge] = nb1d;
3435 aResMap.insert(std::make_pair(sm,aVec));
3436 EdgesMap.Bind(E,nb1d);
3438 double ELen = fullLen/fullNbSeg;
3442 // try to evaluate as in MEFISTO
3443 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
3444 TopoDS_Face F = TopoDS::Face( exp.Current() );
3445 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
3447 BRepGProp::SurfaceProperties(F,G);
3448 double anArea = G.Mass();
3450 std::vector<int> nb1dVec;
3451 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
3452 int nbSeg = EdgesMap.Find(exp1.Current());
3454 nb1dVec.push_back( nbSeg );
3457 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
3458 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
3461 if ( nb1dVec.size() == 4 ) // quadrangle geom face
3463 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
3465 nbNodes = (n1 + 1) * (n2 + 1);
3470 nbTria = nbQuad = nbTria / 3 + 1;
3473 std::vector<int> aVec(SMDSEntity_Last,0);
3475 int nb1d_in = (nbTria*3 - nb1d) / 2;
3476 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3477 aVec[SMDSEntity_Quad_Triangle] = nbTria;
3478 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
3481 aVec[SMDSEntity_Node] = nbNodes;
3482 aVec[SMDSEntity_Triangle] = nbTria;
3483 aVec[SMDSEntity_Quadrangle] = nbQuad;
3485 aResMap.insert(std::make_pair(sm,aVec));
3492 BRepGProp::VolumeProperties(aShape,G);
3493 double aVolume = G.Mass();
3494 double tetrVol = 0.1179*ELen*ELen*ELen;
3495 int nbVols = int(aVolume/tetrVol);
3496 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3497 std::vector<int> aVec(SMDSEntity_Last);
3498 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
3500 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3501 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3504 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3505 aVec[SMDSEntity_Tetra] = nbVols;
3507 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
3508 aResMap.insert(std::make_pair(sm,aVec));