1 // Copyright (C) 2007-2020 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>
41 #include <SMDS_EdgePosition.hxx>
42 #include <SMESHDS_Group.hxx>
43 #include <SMESH_Gen.hxx>
44 #include <SMESH_Group.hxx>
45 #include <SMESH_Mesh.hxx>
46 #include <SMESH_MeshEditor.hxx>
47 #include <SMESH_MesherHelper.hxx>
48 #include <StdMeshers_FaceSide.hxx>
49 #include <StdMeshers_ViscousLayers2D.hxx>
50 #include <SMESH_File.hxx>
52 #include <utilities.h>
60 // OPENCASCADE includes
61 #include <BRepBndLib.hxx>
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 {0, 0, 0, 0} /* 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 {0, 0, 0, 0, 0} /* 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 PyVarObject_HEAD_INIT(NULL, 0)
151 sizeof(PyStdOut), /*tp_basicsize*/
154 (destructor)PyStdOut_dealloc, /*tp_dealloc*/
161 0, /*tp_as_sequence*/
166 PyObject_GenericGetAttr, /*tp_getattro*/
167 /* softspace is writable: we must supply tp_setattro */
168 PyObject_GenericSetAttr, /* tp_setattro */
170 Py_TPFLAGS_DEFAULT, /*tp_flags*/
174 0, /*tp_richcompare*/
175 0, /*tp_weaklistoffset*/
178 PyStdOut_methods, /*tp_methods*/
179 PyStdOut_memberlist, /*tp_members*/
197 0, /*tp_version_tag*/
201 PyObject * newPyStdOut( std::string& out )
203 PyStdOut* self = PyObject_New(PyStdOut, &PyStdOut_Type);
208 return (PyObject*)self;
213 ////////////////////////END PYTHON///////////////////////////
215 //////////////////MY MAPS////////////////////////////////////////
217 TopTools_IndexedMapOfShape FacesWithSizeMap;
218 std::map<int,string> FaceId2SizeMap;
219 TopTools_IndexedMapOfShape EdgesWithSizeMap;
220 std::map<int,string> EdgeId2SizeMap;
221 TopTools_IndexedMapOfShape VerticesWithSizeMap;
222 std::map<int,string> VertexId2SizeMap;
224 std::map<int,PyObject*> FaceId2PythonSmp;
225 std::map<int,PyObject*> EdgeId2PythonSmp;
226 std::map<int,PyObject*> VertexId2PythonSmp;
228 typedef std::map<int, std::vector< BLSURFPlugin_Attractor* > > TId2ClsAttractorVec;
229 TId2ClsAttractorVec FaceId2ClassAttractor;
230 TId2ClsAttractorVec FaceIndex2ClassAttractor;
231 std::map<int,std::vector<double> > FaceId2AttractorCoords;
234 TopTools_IndexedMapOfShape FacesWithEnforcedVertices;
235 std::map< int, BLSURFPlugin_Hypothesis::TEnfVertexCoordsList > FaceId2EnforcedVertexCoords;
236 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexCoords > EnfVertexCoords2ProjVertex;
237 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList > EnfVertexCoords2EnfVertexList;
238 SMESH_MesherHelper* theHelper;
240 bool HasSizeMapOnFace=false;
241 bool HasSizeMapOnEdge=false;
242 bool HasSizeMapOnVertex=false;
243 //bool HasAttractorOnFace=false;
245 //=============================================================================
249 //=============================================================================
251 BLSURFPlugin_BLSURF::BLSURFPlugin_BLSURF(int hypId,
254 : SMESH_2D_Algo(hypId, gen)
256 _name = theHasGEOM ? "MG-CADSurf" : "MG-CADSurf_NOGEOM";//"BLSURF";
257 _shapeType = (1 << TopAbs_FACE); // 1 bit /shape type
258 _compatibleHypothesis.push_back(BLSURFPlugin_Hypothesis::GetHypType(theHasGEOM));
260 _compatibleHypothesis.push_back(StdMeshers_ViscousLayers2D::GetHypType());
261 _requireDiscreteBoundary = false;
262 _onlyUnaryInput = false;
264 _supportSubmeshes = true;
265 _requireShape = theHasGEOM;
267 /* Initialize the Python interpreter */
268 assert(Py_IsInitialized());
269 PyGILState_STATE gstate;
270 gstate = PyGILState_Ensure();
273 main_mod = PyImport_AddModule("__main__");
276 main_dict = PyModule_GetDict(main_mod);
278 PyRun_SimpleString("from math import *");
279 PyGILState_Release(gstate);
281 FacesWithSizeMap.Clear();
282 FaceId2SizeMap.clear();
283 EdgesWithSizeMap.Clear();
284 EdgeId2SizeMap.clear();
285 VerticesWithSizeMap.Clear();
286 VertexId2SizeMap.clear();
287 FaceId2PythonSmp.clear();
288 EdgeId2PythonSmp.clear();
289 VertexId2PythonSmp.clear();
290 FaceId2AttractorCoords.clear();
291 FaceId2ClassAttractor.clear();
292 FaceIndex2ClassAttractor.clear();
293 FacesWithEnforcedVertices.Clear();
294 FaceId2EnforcedVertexCoords.clear();
295 EnfVertexCoords2ProjVertex.clear();
296 EnfVertexCoords2EnfVertexList.clear();
298 _compute_canceled = false;
301 //=============================================================================
305 //=============================================================================
307 BLSURFPlugin_BLSURF::~BLSURFPlugin_BLSURF()
312 //=============================================================================
316 //=============================================================================
318 bool BLSURFPlugin_BLSURF::CheckHypothesis
320 const TopoDS_Shape& aShape,
321 SMESH_Hypothesis::Hypothesis_Status& aStatus)
324 _haveViscousLayers = false;
326 list<const SMESHDS_Hypothesis*>::const_iterator itl;
327 const SMESHDS_Hypothesis* theHyp;
329 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape,
330 /*ignoreAuxiliary=*/false);
331 aStatus = SMESH_Hypothesis::HYP_OK;
334 return true; // can work with no hypothesis
337 for ( itl = hyps.begin(); itl != hyps.end() && ( aStatus == HYP_OK ); ++itl )
340 string hypName = theHyp->GetName();
341 if ( hypName == BLSURFPlugin_Hypothesis::GetHypType(true) ||
342 hypName == BLSURFPlugin_Hypothesis::GetHypType(false) )
344 _hypothesis = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
346 if ( _hypothesis->GetPhysicalMesh() == BLSURFPlugin_Hypothesis::DefaultSize &&
347 _hypothesis->GetGeometricMesh() == BLSURFPlugin_Hypothesis::DefaultGeom )
348 // hphy_flag = 0 and hgeo_flag = 0 is not allowed (spec)
349 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
351 else if ( hypName == StdMeshers_ViscousLayers2D::GetHypType() )
353 if ( !_haveViscousLayers )
355 if ( error( StdMeshers_ViscousLayers2D::CheckHypothesis( aMesh, aShape, aStatus )))
356 _haveViscousLayers = true;
361 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
364 return aStatus == SMESH_Hypothesis::HYP_OK;
367 //=============================================================================
369 * Pass parameters to MG-CADSurf
371 //=============================================================================
373 inline std::string val_to_string(double d)
375 std::ostringstream o;
380 inline std::string val_to_string_rel(double d)
382 std::ostringstream o;
388 inline std::string val_to_string(int i)
390 std::ostringstream o;
395 inline std::string val_to_string_rel(int i)
397 std::ostringstream o;
403 double _smp_phy_size;
404 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data);
405 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data);
406 status_t size_on_vertex(integer vertex_id, real *size, void *user_data);
413 /////////////////////////////////////////////////////////
415 projectionPoint getProjectionPoint(TopoDS_Face& theFace, const gp_Pnt& thePoint)
417 projectionPoint myPoint;
419 if ( theFace.IsNull() )
421 TopoDS_Shape foundFace, myShape = theHelper->GetSubShape();
422 TopTools_MapOfShape checkedFaces;
423 std::map< double, std::pair< TopoDS_Face, gp_Pnt2d > > dist2face;
425 for ( TopExp_Explorer exp ( myShape, TopAbs_FACE ); exp.More(); exp.Next())
427 const TopoDS_Face& face = TopoDS::Face( exp.Current() );
428 if ( !checkedFaces.Add( face )) continue;
430 // check distance to face
431 Handle(ShapeAnalysis_Surface) surface = theHelper->GetSurface( face );
432 gp_Pnt2d uv = surface->ValueOfUV( thePoint, Precision::Confusion());
433 double distance = surface->Gap();
434 if ( distance > Precision::Confusion() )
436 // the face is far, store for future analysis
437 dist2face.insert( std::make_pair( distance, std::make_pair( face, uv )));
441 // check location on the face
442 BRepClass_FaceClassifier FC( face, uv, BRep_Tool::Tolerance( face ));
443 if ( FC.State() == TopAbs_IN )
445 if ( !foundFace.IsNull() )
446 return myPoint; // thePoint seems to be TopAbs_ON
448 myPoint.uv = uv.XY();
449 myPoint.xyz = surface->Value( uv ).XYZ();
452 if ( FC.State() == TopAbs_ON )
456 if ( foundFace.IsNull() )
458 // find the closest face
459 std::map< double, std::pair< TopoDS_Face, gp_Pnt2d > >::iterator d2f = dist2face.begin();
460 for ( ; d2f != dist2face.end(); ++d2f )
462 const TopoDS_Face& face = d2f->second.first;
463 const gp_Pnt2d & uv = d2f->second.second;
464 BRepClass_FaceClassifier FC( face, uv, Precision::Confusion());
465 if ( FC.State() == TopAbs_IN )
468 myPoint.uv = uv.XY();
469 myPoint.xyz = theHelper->GetSurface( face )->Value( uv ).XYZ();
474 // set the resultShape
475 // if ( foundFace.IsNull() )
476 // throw SMESH_ComputeError(COMPERR_BAD_PARMETERS,
477 // "getProjectionPoint: can't find a face by a vertex");
478 theFace = TopoDS::Face( foundFace );
482 Handle(Geom_Surface) surface = BRep_Tool::Surface( theFace );
483 GeomAPI_ProjectPointOnSurf projector( thePoint, surface );
484 if ( !projector.IsDone() || projector.NbPoints()==0 )
485 throw SMESH_ComputeError(COMPERR_BAD_PARMETERS,
486 "getProjectionPoint: can't project a vertex to a face");
489 projector.LowerDistanceParameters(u,v);
490 myPoint.uv = gp_XY(u,v);
491 gp_Pnt aPnt = projector.NearestPoint();
492 myPoint.xyz = gp_XYZ(aPnt.X(),aPnt.Y(),aPnt.Z());
494 BRepClass_FaceClassifier FC( theFace, myPoint.uv, Precision::Confusion());
495 if ( FC.State() != TopAbs_IN )
502 /////////////////////////////////////////////////////////
503 TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
505 GEOM::GEOM_Object_var aGeomObj;
506 TopoDS_Shape S = TopoDS_Shape();
507 SALOMEDS::SObject_var aSObj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->FindObjectID( entry.c_str() );
508 if (!aSObj->_is_nil()) {
509 CORBA::Object_var obj = aSObj->GetObject();
510 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
513 if ( !aGeomObj->_is_nil() )
514 S = SMESH_Gen_i::GetSMESHGen()->GeomObjectToShape( aGeomObj.in() );
518 void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
520 BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
522 // Find the face and get the (u,v) values of the enforced vertex on the face
523 projectionPoint myPoint = getProjectionPoint(faceShape,aPnt);
524 if ( faceShape.IsNull() )
527 enf_coords.push_back(aPnt.X());
528 enf_coords.push_back(aPnt.Y());
529 enf_coords.push_back(aPnt.Z());
531 coords.push_back(myPoint.uv.X());
532 coords.push_back(myPoint.uv.Y());
533 coords.push_back(myPoint.xyz.X());
534 coords.push_back(myPoint.xyz.Y());
535 coords.push_back(myPoint.xyz.Z());
537 s_coords.push_back(myPoint.xyz.X());
538 s_coords.push_back(myPoint.xyz.Y());
539 s_coords.push_back(myPoint.xyz.Z());
541 // Save pair projected vertex / enf vertex
542 EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
543 pair<BLSURFPlugin_Hypothesis::TEnfVertexList::iterator,bool> ret;
544 BLSURFPlugin_Hypothesis::TEnfVertexList::iterator it;
545 ret = EnfVertexCoords2EnfVertexList[s_coords].insert(enfVertex);
546 if (ret.second == false) {
548 (*it)->grpName = enfVertex->grpName;
552 if (! FacesWithEnforcedVertices.Contains(faceShape)) {
553 key = FacesWithEnforcedVertices.Add(faceShape);
556 key = FacesWithEnforcedVertices.FindIndex(faceShape);
559 // If a node is already created by an attractor, do not create enforced vertex
560 int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
561 bool sameAttractor = false;
562 if (attractorKey >= 0)
563 if (FaceId2AttractorCoords.count(attractorKey) > 0)
564 if (FaceId2AttractorCoords[attractorKey] == coords)
565 sameAttractor = true;
567 if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
569 FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
572 if (! sameAttractor) {
573 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList ens;
575 FaceId2EnforcedVertexCoords[key] = ens;
580 /////////////////////////////////////////////////////////
581 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
583 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
586 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
588 for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
589 enfVertex = *enfVertexListIt;
590 // Case of manual coords
591 if (enfVertex->coords.size() != 0) {
592 aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
593 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
596 // Case of geom vertex coords
597 if (enfVertex->geomEntry != "") {
598 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
599 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
600 if (GeomType == TopAbs_VERTEX)
602 enfVertex->vertex = TopoDS::Vertex( GeomShape );
603 aPnt = BRep_Tool::Pnt( enfVertex->vertex );
604 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
607 if (GeomType == TopAbs_COMPOUND)
609 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next())
610 if (it.Value().ShapeType() == TopAbs_VERTEX)
612 enfVertex->vertex = TopoDS::Vertex( it.Value() );
613 aPnt = BRep_Tool::Pnt( enfVertex->vertex );
614 _createEnforcedVertexOnFace( TopoDS::Face(faceShape), aPnt, enfVertex);
621 /////////////////////////////////////////////////////////
622 void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction, double defaultSize)
624 double xa=0., ya=0., za=0.; // Coordinates of attractor point
625 double a, b; // Attractor parameter
627 bool createNode=false; // To create a node on attractor projection
629 const char *sep = ";";
630 // atIt->second has the following pattern:
631 // ATTRACTOR(xa;ya;za;a;b;True|False;d)
633 // xa;ya;za : coordinates of attractor
634 // a : desired size on attractor
635 // b : distance of influence of attractor
636 // d : distance until which the size remains constant
638 // We search the parameters in the string
640 pos1 = AttractorFunction.find(sep);
641 if (pos1!=string::npos)
642 xa = atof(AttractorFunction.substr(10, pos1-10).c_str());
644 pos2 = AttractorFunction.find(sep, pos1+1);
645 if (pos2!=string::npos) {
646 ya = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
650 pos2 = AttractorFunction.find(sep, pos1+1);
651 if (pos2!=string::npos) {
652 za = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
656 pos2 = AttractorFunction.find(sep, pos1+1);
657 if (pos2!=string::npos) {
658 a = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
662 pos2 = AttractorFunction.find(sep, pos1+1);
663 if (pos2!=string::npos) {
664 b = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
668 pos2 = AttractorFunction.find(sep, pos1+1);
669 if (pos2!=string::npos) {
670 string createNodeStr = AttractorFunction.substr(pos1+1, pos2-pos1-1);
671 createNode = (AttractorFunction.substr(pos1+1, pos2-pos1-1) == "True");
675 pos2 = AttractorFunction.find(")");
676 if (pos2!=string::npos) {
677 d = atof(AttractorFunction.substr(pos1+1, pos2-pos1-1).c_str());
680 // Get the (u,v) values of the attractor on the face
681 projectionPoint myPoint = getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
682 gp_XY uvPoint = myPoint.uv;
683 gp_XYZ xyzPoint = myPoint.xyz;
684 Standard_Real u0 = uvPoint.X();
685 Standard_Real v0 = uvPoint.Y();
686 Standard_Real x0 = xyzPoint.X();
687 Standard_Real y0 = xyzPoint.Y();
688 Standard_Real z0 = xyzPoint.Z();
689 std::vector<double> coords;
690 coords.push_back(u0);
691 coords.push_back(v0);
692 coords.push_back(x0);
693 coords.push_back(y0);
694 coords.push_back(z0);
695 // We construct the python function
696 ostringstream attractorFunctionStream;
697 attractorFunctionStream << "def f(u,v): return ";
698 attractorFunctionStream << defaultSize << "-(" << defaultSize <<"-" << a << ")";
699 //attractorFunctionStream << "*exp(-((u-("<<u0<<"))*(u-("<<u0<<"))+(v-("<<v0<<"))*(v-("<<v0<<")))/(" << b << "*" << b <<"))";
700 // rnc: make possible to keep the size constant until
701 // a defined distance. Distance is expressed as the positiv part
702 // of r-d where r is the distance to (u0,v0)
703 attractorFunctionStream << "*exp(-(0.5*(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"+abs(sqrt((u-"<<u0<<")**2+(v-"<<v0<<")**2)-"<<d<<"))/(" << b << "))**2)";
706 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
707 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
710 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
712 FaceId2SizeMap[key] =attractorFunctionStream.str();
714 FaceId2AttractorCoords[key] = coords;
716 // // Test for new attractors
717 // gp_Pnt myP(xyzPoint);
718 // TopoDS_Vertex myV = BRepBuilderAPI_MakeVertex(myP);
719 // BLSURFPlugin_Attractor myAttractor(TopoDS::Face(GeomShape),myV,200);
720 // myAttractor.SetParameters(a, defaultSize, b, d);
721 // myAttractor.SetType(1);
722 // FaceId2ClassAttractor[key] = myAttractor;
723 // if(FaceId2ClassAttractor[key].GetFace().IsNull()){
727 // One sub-shape to get ids from
728 BLSURFPlugin_BLSURF::TListOfIDs _getSubShapeIDsInMainShape(const TopoDS_Shape& theMainShape,
729 const TopoDS_Shape& theSubShape,
730 TopAbs_ShapeEnum theShapeType)
732 BLSURFPlugin_BLSURF::TListOfIDs face_ids;
734 TopTools_MapOfShape subShapes;
735 TopTools_IndexedMapOfShape anIndices;
736 TopExp::MapShapes(theMainShape, theShapeType, anIndices);
738 for (TopExp_Explorer face_iter(theSubShape,theShapeType);face_iter.More();face_iter.Next())
740 if ( subShapes.Add( face_iter.Current() )) // issue 23416
742 int face_id = anIndices.FindIndex( face_iter.Current() );
744 throw SALOME_Exception( "Periodicity: sub_shape not found in main_shape");
745 face_ids.push_back( face_id );
751 BLSURFPlugin_BLSURF::TListOfIDs _getSubShapeIDsInMainShape(SMESH_Mesh* theMesh,
752 TopoDS_Shape theSubShape,
753 TopAbs_ShapeEnum theShapeType)
755 BLSURFPlugin_BLSURF::TListOfIDs face_ids;
757 for (TopExp_Explorer face_iter(theSubShape,theShapeType);face_iter.More();face_iter.Next())
759 int face_id = theMesh->GetMeshDS()->ShapeToIndex(face_iter.Current());
761 throw SALOME_Exception ( "Periodicity: sub_shape not found in main_shape");
762 face_ids.push_back(face_id);
768 void BLSURFPlugin_BLSURF::addCoordsFromVertices(const std::vector<std::string> &theVerticesEntries, std::vector<double> &theVerticesCoords)
770 for (std::vector<std::string>::const_iterator it = theVerticesEntries.begin(); it != theVerticesEntries.end(); it++)
772 BLSURFPlugin_Hypothesis::TEntry theVertexEntry = *it;
773 addCoordsFromVertex(theVertexEntry, theVerticesCoords);
778 void BLSURFPlugin_BLSURF::addCoordsFromVertex(BLSURFPlugin_Hypothesis::TEntry theVertexEntry, std::vector<double> &theVerticesCoords)
780 if (theVertexEntry!="")
782 TopoDS_Shape aShape = entryToShape(theVertexEntry);
784 gp_Pnt aPnt = BRep_Tool::Pnt( TopoDS::Vertex( aShape ) );
785 double theX, theY, theZ;
790 theVerticesCoords.push_back(theX);
791 theVerticesCoords.push_back(theY);
792 theVerticesCoords.push_back(theZ);
796 /////////////////////////////////////////////////////////
797 void BLSURFPlugin_BLSURF::createPreCadFacesPeriodicity(TopoDS_Shape theGeomShape, const BLSURFPlugin_Hypothesis::TPreCadPeriodicity &preCadPeriodicity)
799 TopoDS_Shape geomShape1 = entryToShape(preCadPeriodicity.shape1Entry);
800 TopoDS_Shape geomShape2 = entryToShape(preCadPeriodicity.shape2Entry);
802 TListOfIDs theFace1_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape1, TopAbs_FACE);
803 TListOfIDs theFace2_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape2, TopAbs_FACE);
805 TPreCadPeriodicityIDs preCadFacesPeriodicityIDs;
806 preCadFacesPeriodicityIDs.shape1IDs = theFace1_ids;
807 preCadFacesPeriodicityIDs.shape2IDs = theFace2_ids;
809 addCoordsFromVertices(preCadPeriodicity.theSourceVerticesEntries, preCadFacesPeriodicityIDs.theSourceVerticesCoords);
810 addCoordsFromVertices(preCadPeriodicity.theTargetVerticesEntries, preCadFacesPeriodicityIDs.theTargetVerticesCoords);
812 _preCadFacesIDsPeriodicityVector.push_back(preCadFacesPeriodicityIDs);
815 /////////////////////////////////////////////////////////
816 void BLSURFPlugin_BLSURF::createPreCadEdgesPeriodicity(TopoDS_Shape theGeomShape, const BLSURFPlugin_Hypothesis::TPreCadPeriodicity &preCadPeriodicity)
818 TopoDS_Shape geomShape1 = entryToShape(preCadPeriodicity.shape1Entry);
819 TopoDS_Shape geomShape2 = entryToShape(preCadPeriodicity.shape2Entry);
821 TListOfIDs theEdge1_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape1, TopAbs_EDGE);
822 TListOfIDs theEdge2_ids = _getSubShapeIDsInMainShape(theGeomShape, geomShape2, TopAbs_EDGE);
824 TPreCadPeriodicityIDs preCadEdgesPeriodicityIDs;
825 preCadEdgesPeriodicityIDs.shape1IDs = theEdge1_ids;
826 preCadEdgesPeriodicityIDs.shape2IDs = theEdge2_ids;
828 addCoordsFromVertices(preCadPeriodicity.theSourceVerticesEntries, preCadEdgesPeriodicityIDs.theSourceVerticesCoords);
829 addCoordsFromVertices(preCadPeriodicity.theTargetVerticesEntries, preCadEdgesPeriodicityIDs.theTargetVerticesCoords);
831 _preCadEdgesIDsPeriodicityVector.push_back(preCadEdgesPeriodicityIDs);
835 /////////////////////////////////////////////////////////
837 void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
838 cadsurf_session_t * css,
839 const TopoDS_Shape& theGeomShape)
842 // Clear map so that it is not stored in the algorithm with old enforced vertices in it
843 FacesWithSizeMap.Clear();
844 FaceId2SizeMap.clear();
845 EdgesWithSizeMap.Clear();
846 EdgeId2SizeMap.clear();
847 VerticesWithSizeMap.Clear();
848 VertexId2SizeMap.clear();
849 FaceId2PythonSmp.clear();
850 EdgeId2PythonSmp.clear();
851 VertexId2PythonSmp.clear();
852 FaceId2AttractorCoords.clear();
853 FaceId2ClassAttractor.clear();
854 FaceIndex2ClassAttractor.clear();
855 FacesWithEnforcedVertices.Clear();
856 FaceId2EnforcedVertexCoords.clear();
857 EnfVertexCoords2ProjVertex.clear();
858 EnfVertexCoords2EnfVertexList.clear();
860 double diagonal = SMESH_Mesh::GetShapeDiagonalSize( theGeomShape );
861 double bbSegmentation = _gen->GetBoundaryBoxSegmentation();
862 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
863 int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
864 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize(diagonal, bbSegmentation);
865 bool _phySizeRel = BLSURFPlugin_Hypothesis::GetDefaultPhySizeRel();
866 double _minSize = BLSURFPlugin_Hypothesis::GetDefaultMinSize(diagonal);
867 bool _minSizeRel = BLSURFPlugin_Hypothesis::GetDefaultMinSizeRel();
868 double _maxSize = BLSURFPlugin_Hypothesis::GetDefaultMaxSize(diagonal);
869 bool _maxSizeRel = BLSURFPlugin_Hypothesis::GetDefaultMaxSizeRel();
870 double _use_gradation = BLSURFPlugin_Hypothesis::GetDefaultUseGradation();
871 double _gradation = BLSURFPlugin_Hypothesis::GetDefaultGradation();
872 double _use_volume_gradation = BLSURFPlugin_Hypothesis::GetDefaultUseVolumeGradation();
873 double _volume_gradation = BLSURFPlugin_Hypothesis::GetDefaultVolumeGradation();
874 BLSURFPlugin_Hypothesis::ElementType _elementType = BLSURFPlugin_Hypothesis::GetDefaultElementType();
875 double _angleMesh = BLSURFPlugin_Hypothesis::GetDefaultAngleMesh();
876 double _chordalError = BLSURFPlugin_Hypothesis::GetDefaultChordalError(diagonal);
877 bool _anisotropic = BLSURFPlugin_Hypothesis::GetDefaultAnisotropic();
878 double _anisotropicRatio = BLSURFPlugin_Hypothesis::GetDefaultAnisotropicRatio();
879 bool _removeTinyEdges = BLSURFPlugin_Hypothesis::GetDefaultRemoveTinyEdges();
880 double _tinyEdgeLength = BLSURFPlugin_Hypothesis::GetDefaultTinyEdgeLength(diagonal);
881 bool _optimiseTinyEdges = BLSURFPlugin_Hypothesis::GetDefaultOptimiseTinyEdges();
882 double _tinyEdgeOptimisLength = BLSURFPlugin_Hypothesis::GetDefaultTinyEdgeOptimisationLength(diagonal);
883 bool _correctSurfaceIntersec= BLSURFPlugin_Hypothesis::GetDefaultCorrectSurfaceIntersection();
884 double _corrSurfaceIntersCost = BLSURFPlugin_Hypothesis::GetDefaultCorrectSurfaceIntersectionMaxCost();
885 bool _badElementRemoval = BLSURFPlugin_Hypothesis::GetDefaultBadElementRemoval();
886 double _badElementAspectRatio = BLSURFPlugin_Hypothesis::GetDefaultBadElementAspectRatio();
887 bool _optimizeMesh = BLSURFPlugin_Hypothesis::GetDefaultOptimizeMesh();
888 bool _quadraticMesh = BLSURFPlugin_Hypothesis::GetDefaultQuadraticMesh();
889 int _verb = BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
890 //int _topology = BLSURFPlugin_Hypothesis::GetDefaultTopology();
891 bool _useSurfaceProximity = BLSURFPlugin_Hypothesis::GetDefaultUseSurfaceProximity ();
892 int _nbSurfaceProximityLayers = BLSURFPlugin_Hypothesis::GetDefaultNbSurfaceProximityLayers();
893 double _surfaceProximityRatio = BLSURFPlugin_Hypothesis::GetDefaultSurfaceProximityRatio ();
894 bool _useVolumeProximity = BLSURFPlugin_Hypothesis::GetDefaultUseVolumeProximity ();
895 int _nbVolumeProximityLayers = BLSURFPlugin_Hypothesis::GetDefaultNbVolumeProximityLayers ();
896 double _volumeProximityRatio = BLSURFPlugin_Hypothesis::GetDefaultVolumeProximityRatio ();
899 //int _precadMergeEdges = BLSURFPlugin_Hypothesis::GetDefaultPreCADMergeEdges();
900 //int _precadRemoveDuplicateCADFaces = BLSURFPlugin_Hypothesis::GetDefaultPreCADRemoveDuplicateCADFaces();
901 //int _precadProcess3DTopology = BLSURFPlugin_Hypothesis::GetDefaultPreCADProcess3DTopology();
902 //int _precadDiscardInput = BLSURFPlugin_Hypothesis::GetDefaultPreCADDiscardInput();
904 const BLSURFPlugin_Hypothesis::TPreCadPeriodicityVector preCadFacesPeriodicityVector = BLSURFPlugin_Hypothesis::GetPreCadFacesPeriodicityVector(hyp);
907 _physicalMesh = (int) hyp->GetPhysicalMesh();
908 _geometricMesh = (int) hyp->GetGeometricMesh();
909 if (hyp->GetPhySize() > 0) {
910 _phySize = hyp->GetPhySize();
911 // if user size is not explicitly specified, "relative" flag is ignored
912 _phySizeRel = hyp->IsPhySizeRel();
914 if (hyp->GetMinSize() > 0) {
915 _minSize = hyp->GetMinSize();
916 // if min size is not explicitly specified, "relative" flag is ignored
917 _minSizeRel = hyp->IsMinSizeRel();
919 if (hyp->GetMaxSize() > 0) {
920 _maxSize = hyp->GetMaxSize();
921 // if max size is not explicitly specified, "relative" flag is ignored
922 _maxSizeRel = hyp->IsMaxSizeRel();
924 _use_gradation = hyp->GetUseGradation();
925 if (hyp->GetGradation() > 0 && _use_gradation)
926 _gradation = hyp->GetGradation();
927 _use_volume_gradation = hyp->GetUseVolumeGradation();
928 if (hyp->GetVolumeGradation() > 0 && _use_volume_gradation )
929 _volume_gradation = hyp->GetVolumeGradation();
930 _elementType = hyp->GetElementType();
931 if (hyp->GetAngleMesh() > 0)
932 _angleMesh = hyp->GetAngleMesh();
933 if (hyp->GetChordalError() > 0)
934 _chordalError = hyp->GetChordalError();
935 _anisotropic = hyp->GetAnisotropic();
936 if (hyp->GetAnisotropicRatio() >= 0)
937 _anisotropicRatio = hyp->GetAnisotropicRatio();
938 _removeTinyEdges = hyp->GetRemoveTinyEdges();
939 if (hyp->GetTinyEdgeLength() > 0)
940 _tinyEdgeLength = hyp->GetTinyEdgeLength();
941 _optimiseTinyEdges = hyp->GetOptimiseTinyEdges();
942 if (hyp->GetTinyEdgeOptimisationLength() > 0)
943 _tinyEdgeOptimisLength = hyp->GetTinyEdgeOptimisationLength();
944 _correctSurfaceIntersec = hyp->GetCorrectSurfaceIntersection();
945 if (hyp->GetCorrectSurfaceIntersectionMaxCost() > 0)
946 _corrSurfaceIntersCost = hyp->GetCorrectSurfaceIntersectionMaxCost();
947 _badElementRemoval = hyp->GetBadElementRemoval();
948 if (hyp->GetBadElementAspectRatio() >= 0)
949 _badElementAspectRatio = hyp->GetBadElementAspectRatio();
950 _optimizeMesh = hyp->GetOptimizeMesh();
951 _quadraticMesh = hyp->GetQuadraticMesh();
952 _verb = hyp->GetVerbosity();
953 //_topology = (int) hyp->GetTopology();
955 //_precadMergeEdges = hyp->GetPreCADMergeEdges();
956 //_precadRemoveDuplicateCADFaces = hyp->GetPreCADRemoveDuplicateCADFaces();
957 //_precadProcess3DTopology = hyp->GetPreCADProcess3DTopology();
958 //_precadDiscardInput = hyp->GetPreCADDiscardInput();
959 _useSurfaceProximity = hyp->GetUseSurfaceProximity ();
960 _nbSurfaceProximityLayers = hyp->GetNbSurfaceProximityLayers();
961 _surfaceProximityRatio = hyp->GetSurfaceProximityRatio ();
962 _useVolumeProximity = hyp->GetUseVolumeProximity ();
963 _nbVolumeProximityLayers = hyp->GetNbVolumeProximityLayers ();
964 _volumeProximityRatio = hyp->GetVolumeProximityRatio ();
967 const BLSURFPlugin_Hypothesis::TOptionValues& opts = hyp->GetOptionValues();
968 BLSURFPlugin_Hypothesis::TOptionValues::const_iterator opIt;
969 for ( opIt = opts.begin(); opIt != opts.end(); ++opIt ){
970 MESSAGE("OptionValue: " << opIt->first.c_str() << ", value: " << opIt->second.c_str());
971 if ( !opIt->second.empty() ) {
972 // With MeshGems 2.4-5, there are issues with periodicity and multithread
973 // => As a temporary workaround, we enforce to use only one thread if periodicity is used.
974 if (opIt->first == "max_number_of_threads" && opIt->second != "1" && ! preCadFacesPeriodicityVector.empty()){
975 std::cout << "INFO: Disabling multithread to avoid periodicity issues" << std::endl;
976 set_param(css, opIt->first.c_str(), "1");
979 set_param(css, opIt->first.c_str(), opIt->second.c_str());
983 const BLSURFPlugin_Hypothesis::TOptionValues& custom_opts = hyp->GetCustomOptionValues();
984 for ( opIt = custom_opts.begin(); opIt != custom_opts.end(); ++opIt )
985 if ( !opIt->second.empty() ) {
986 set_param(css, opIt->first.c_str(), opIt->second.c_str());
989 const BLSURFPlugin_Hypothesis::TOptionValues& preCADopts = hyp->GetPreCADOptionValues();
990 for ( opIt = preCADopts.begin(); opIt != preCADopts.end(); ++opIt )
991 if ( !opIt->second.empty() ) {
992 set_param(css, opIt->first.c_str(), opIt->second.c_str());
995 if ( hyp->GetHyperPatches().size() < hyp->GetHyperPatchEntries().size() )
997 std::map< std::string, TopoDS_Shape > entryToShape;
998 FillEntryToShape( hyp, entryToShape );
999 const_cast<BLSURFPlugin_Hypothesis*>( hyp )->SetHyperPatchIDsByEntry( theGeomShape,
1004 if ( BLSURFPlugin_Hypothesis::HasPreCADOptions( hyp ))
1006 cadsurf_set_param(css, "use_precad", "yes" ); // for old versions
1008 // PreProcessor (formerly PreCAD) -- commented params are preCADoptions (since 0023307)
1009 //set_param(css, "merge_edges", _precadMergeEdges ? "yes" : "no");
1010 //set_param(css, "remove_duplicate_cad_faces", _precadRemoveDuplicateCADFaces ? "yes" : "no");
1011 //set_param(css, "process_3d_topology", _precadProcess3DTopology ? "1" : "0");
1012 //set_param(css, "discard_input_topology", _precadDiscardInput ? "1" : "0");
1013 //set_param(css, "max_number_of_points_per_patch", "1000000");
1015 bool useGradation = false;
1016 switch (_physicalMesh)
1018 case BLSURFPlugin_Hypothesis::PhysicalGlobalSize:
1019 set_param(css, "physical_size_mode", "global");
1020 set_param(css, "global_physical_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1021 //useGradation = true;
1023 case BLSURFPlugin_Hypothesis::PhysicalLocalSize:
1024 set_param(css, "physical_size_mode", "local");
1025 set_param(css, "global_physical_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1026 //useGradation = true;
1029 set_param(css, "physical_size_mode", "none");
1032 switch (_geometricMesh)
1034 case BLSURFPlugin_Hypothesis::GeometricalGlobalSize:
1035 set_param(css, "geometric_size_mode", "global");
1036 set_param(css, "geometric_approximation", val_to_string(_angleMesh).c_str());
1037 set_param(css, "chordal_error", val_to_string(_chordalError).c_str());
1038 //useGradation = true;
1040 case BLSURFPlugin_Hypothesis::GeometricalLocalSize:
1041 set_param(css, "geometric_size_mode", "local");
1042 set_param(css, "geometric_approximation", val_to_string(_angleMesh).c_str());
1043 set_param(css, "chordal_error", val_to_string(_chordalError).c_str());
1044 //useGradation = true;
1047 set_param(css, "geometric_size_mode", "none");
1050 if ( hyp && hyp->GetPhySize() > 0 ) {
1051 // user size is explicitly specified via hypothesis parameters
1052 // min and max sizes should be compared with explicitly specified user size
1053 // - compute absolute min size
1054 double mins = _minSizeRel ? _minSize * diagonal : _minSize;
1055 // - min size should not be greater than user size
1056 if ( _phySize < mins )
1057 set_param(css, "min_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1059 set_param(css, "min_size", _minSizeRel ? val_to_string_rel(_minSize).c_str() : val_to_string(_minSize).c_str());
1060 // - compute absolute max size
1061 double maxs = _maxSizeRel ? _maxSize * diagonal : _maxSize;
1062 // - max size should not be less than user size
1063 if ( _phySize > maxs )
1064 set_param(css, "max_size", _phySizeRel ? val_to_string_rel(_phySize).c_str() : val_to_string(_phySize).c_str());
1066 set_param(css, "max_size", _maxSizeRel ? val_to_string_rel(_maxSize).c_str() : val_to_string(_maxSize).c_str());
1069 // user size is not explicitly specified
1070 // - if minsize is not explicitly specified, we pass default value computed automatically, in this case "relative" flag is ignored
1071 set_param(css, "min_size", _minSizeRel ? val_to_string_rel(_minSize).c_str() : val_to_string(_minSize).c_str());
1072 // - if maxsize is not explicitly specified, we pass default value computed automatically, in this case "relative" flag is ignored
1073 set_param(css, "max_size", _maxSizeRel ? val_to_string_rel(_maxSize).c_str() : val_to_string(_maxSize).c_str());
1075 useGradation = true; // bos #18758
1076 // anisotropic and quadrangle mesh requires disabling gradation
1077 if ( _anisotropic && _elementType != BLSURFPlugin_Hypothesis::Triangles )
1078 useGradation = false; // limitation of V1.3
1079 if ( useGradation && _use_gradation )
1080 set_param(css, "gradation", val_to_string(_gradation).c_str());
1081 if ( useGradation && _use_volume_gradation )
1082 set_param(css, "volume_gradation", val_to_string(_volume_gradation).c_str());
1084 // New since MeshGems 2.5: add full_quad
1085 const char * element_generation = "";
1086 switch ( _elementType )
1088 case BLSURFPlugin_Hypothesis::Triangles:
1089 element_generation = "triangle";
1091 case BLSURFPlugin_Hypothesis::QuadrangleDominant:
1092 element_generation = "quad_dominant";
1094 case BLSURFPlugin_Hypothesis::Quadrangles:
1095 element_generation = "full_quad";
1099 set_param(css, "element_generation", element_generation);
1102 set_param(css, "metric", _anisotropic ? "anisotropic" : "isotropic");
1104 set_param(css, "anisotropic_ratio", val_to_string(_anisotropicRatio).c_str());
1105 set_param(css, "remove_tiny_edges", _removeTinyEdges ? "1" : "0");
1106 if ( _removeTinyEdges )
1107 set_param(css, "tiny_edge_length", val_to_string(_tinyEdgeLength).c_str());
1108 set_param(css, "optimise_tiny_edges", _optimiseTinyEdges ? "1" : "0");
1109 if ( _optimiseTinyEdges )
1110 set_param(css, "tiny_edge_optimisation_length", val_to_string(_tinyEdgeOptimisLength).c_str());
1111 set_param(css, "correct_surface_intersections", _correctSurfaceIntersec ? "1" : "0");
1112 if ( _correctSurfaceIntersec )
1113 set_param(css, "surface_intersections_processing_max_cost", val_to_string(_corrSurfaceIntersCost ).c_str());
1114 set_param(css, "force_bad_surface_element_removal", _badElementRemoval ? "1" : "0");
1115 if ( _badElementRemoval )
1116 set_param(css, "bad_surface_element_aspect_ratio", val_to_string(_badElementAspectRatio).c_str());
1117 set_param(css, "optimisation", _optimizeMesh ? "yes" : "no");
1118 set_param(css, "element_order", _quadraticMesh ? "quadratic" : "linear");
1119 set_param(css, "verbose", val_to_string(_verb).c_str());
1121 set_param(css, "use_surface_proximity", _useSurfaceProximity ? "yes" : "no" );
1122 if ( _useSurfaceProximity )
1124 set_param(css, "surface_proximity_layers", SMESH_Comment( _nbSurfaceProximityLayers ));
1125 set_param(css, "surface_proximity_ratio", SMESH_Comment( _surfaceProximityRatio ));
1127 set_param(css, "use_volume_proximity", _useVolumeProximity ? "yes" : "no" );
1128 if ( _useVolumeProximity )
1130 set_param(css, "volume_proximity_layers", SMESH_Comment( _nbVolumeProximityLayers ));
1131 set_param(css, "volume_proximity_ratio", SMESH_Comment( _volumeProximityRatio ));
1134 _smp_phy_size = _phySizeRel ? _phySize*diagonal : _phySize;
1136 std::cout << "_smp_phy_size = " << _smp_phy_size << std::endl;
1138 if (_physicalMesh == BLSURFPlugin_Hypothesis::PhysicalLocalSize)
1140 TopoDS_Shape GeomShape;
1141 TopoDS_Shape AttShape;
1142 TopAbs_ShapeEnum GeomType;
1144 // Standard Size Maps
1146 const BLSURFPlugin_Hypothesis::TSizeMap sizeMaps = BLSURFPlugin_Hypothesis::GetSizeMapEntries(hyp);
1147 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator smIt = sizeMaps.begin();
1148 for ( ; smIt != sizeMaps.end(); ++smIt ) {
1149 if ( !smIt->second.empty() ) {
1150 GeomShape = entryToShape(smIt->first);
1151 GeomType = GeomShape.ShapeType();
1154 if (GeomType == TopAbs_COMPOUND) {
1155 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1157 if (it.Value().ShapeType() == TopAbs_FACE){
1158 HasSizeMapOnFace = true;
1159 if (! FacesWithSizeMap.Contains(TopoDS::Face(it.Value()))) {
1160 key = FacesWithSizeMap.Add(TopoDS::Face(it.Value()));
1163 key = FacesWithSizeMap.FindIndex(TopoDS::Face(it.Value()));
1165 FaceId2SizeMap[key] = smIt->second;
1168 if (it.Value().ShapeType() == TopAbs_EDGE){
1169 HasSizeMapOnEdge = true;
1170 HasSizeMapOnFace = true;
1171 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(it.Value()))) {
1172 key = EdgesWithSizeMap.Add(TopoDS::Edge(it.Value()));
1175 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(it.Value()));
1177 EdgeId2SizeMap[key] = smIt->second;
1179 // Group of vertices
1180 if (it.Value().ShapeType() == TopAbs_VERTEX){
1181 HasSizeMapOnVertex = true;
1182 HasSizeMapOnEdge = true;
1183 HasSizeMapOnFace = true;
1184 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(it.Value()))) {
1185 key = VerticesWithSizeMap.Add(TopoDS::Vertex(it.Value()));
1188 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(it.Value()));
1190 VertexId2SizeMap[key] = smIt->second;
1195 if (GeomType == TopAbs_FACE){
1196 HasSizeMapOnFace = true;
1197 if (! FacesWithSizeMap.Contains(TopoDS::Face(GeomShape))) {
1198 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape));
1201 key = FacesWithSizeMap.FindIndex(TopoDS::Face(GeomShape));
1203 FaceId2SizeMap[key] = smIt->second;
1206 if (GeomType == TopAbs_EDGE){
1207 HasSizeMapOnEdge = true;
1208 HasSizeMapOnFace = true;
1209 if (! EdgesWithSizeMap.Contains(TopoDS::Edge(GeomShape))) {
1210 key = EdgesWithSizeMap.Add(TopoDS::Edge(GeomShape));
1213 key = EdgesWithSizeMap.FindIndex(TopoDS::Edge(GeomShape));
1215 EdgeId2SizeMap[key] = smIt->second;
1218 if (GeomType == TopAbs_VERTEX){
1219 HasSizeMapOnVertex = true;
1220 HasSizeMapOnEdge = true;
1221 HasSizeMapOnFace = true;
1222 if (! VerticesWithSizeMap.Contains(TopoDS::Vertex(GeomShape))) {
1223 key = VerticesWithSizeMap.Add(TopoDS::Vertex(GeomShape));
1226 key = VerticesWithSizeMap.FindIndex(TopoDS::Vertex(GeomShape));
1228 VertexId2SizeMap[key] = smIt->second;
1236 // TODO appeler le constructeur des attracteurs directement ici
1237 // if ( !_phySizeRel ) {
1238 const BLSURFPlugin_Hypothesis::TSizeMap attractors = BLSURFPlugin_Hypothesis::GetAttractorEntries(hyp);
1239 BLSURFPlugin_Hypothesis::TSizeMap::const_iterator atIt = attractors.begin();
1240 for ( ; atIt != attractors.end(); ++atIt ) {
1241 if ( !atIt->second.empty() ) {
1242 GeomShape = entryToShape(atIt->first);
1243 GeomType = GeomShape.ShapeType();
1245 if (GeomType == TopAbs_COMPOUND){
1246 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1247 if (it.Value().ShapeType() == TopAbs_FACE){
1248 HasSizeMapOnFace = true;
1249 createAttractorOnFace(it.Value(), atIt->second, _phySizeRel ? _phySize*diagonal : _phySize);
1254 if (GeomType == TopAbs_FACE){
1255 HasSizeMapOnFace = true;
1256 createAttractorOnFace(GeomShape, atIt->second, _phySizeRel ? _phySize*diagonal : _phySize);
1259 if (GeomType == TopAbs_EDGE){
1260 HasSizeMapOnEdge = true;
1261 HasSizeMapOnFace = true;
1262 EdgeId2SizeMap[TopoDS::Edge(GeomShape).HashCode(IntegerLast())] = atIt->second;
1264 if (GeomType == TopAbs_VERTEX){
1265 HasSizeMapOnVertex = true;
1266 HasSizeMapOnEdge = true;
1267 HasSizeMapOnFace = true;
1268 VertexId2SizeMap[TopoDS::Vertex(GeomShape).HashCode(IntegerLast())] = atIt->second;
1276 // temporary commented out for testing
1278 // - Fill in the BLSURFPlugin_Hypothesis::TAttractorMap map in the hypothesis
1279 // - Uncomment and complete this part to construct the attractors from the attractor shape and the passed parameters on each face of the map
1280 // - To do this use the public methodss: SetParameters(several double parameters) and SetType(int type)
1282 // - Construct the attractors with an empty dist. map in the hypothesis
1283 // - 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()
1284 // -> define a bool _mapbuilt in the class that is set to false by default and set to true when calling _buildmap() OK
1286 theNbAttractors = 0;
1287 const BLSURFPlugin_Hypothesis::TAttractorMap class_attractors = BLSURFPlugin_Hypothesis::GetClassAttractorEntries(hyp);
1289 BLSURFPlugin_Hypothesis::TAttractorMap::const_iterator AtIt = class_attractors.begin();
1290 for ( ; AtIt != class_attractors.end(); ++AtIt ) {
1291 if ( !AtIt->second->Empty() ) {
1292 GeomShape = entryToShape(AtIt->first);
1293 if ( !SMESH_MesherHelper::IsSubShape( GeomShape, theGeomShape ))
1295 AttShape = AtIt->second->GetAttractorShape();
1296 GeomType = GeomShape.ShapeType();
1298 // if (GeomType == TopAbs_COMPOUND){
1299 // for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1300 // if (it.Value().ShapeType() == TopAbs_FACE){
1301 // HasAttractorOnFace = true;
1302 // myAttractor = BLSURFPluginAttractor(GeomShape, AttShape);
1307 if (GeomType == TopAbs_FACE
1308 && (AttShape.ShapeType() == TopAbs_VERTEX
1309 || AttShape.ShapeType() == TopAbs_EDGE
1310 || AttShape.ShapeType() == TopAbs_WIRE
1311 || AttShape.ShapeType() == TopAbs_COMPOUND) ){
1312 HasSizeMapOnFace = true;
1314 key = FacesWithSizeMap.Add(TopoDS::Face(GeomShape) );
1316 FaceId2ClassAttractor[key].push_back( AtIt->second );
1320 MESSAGE("Wrong shape type !!")
1328 // Enforced Vertices
1330 const BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap entryEnfVertexListMap = BLSURFPlugin_Hypothesis::GetAllEnforcedVerticesByFace(hyp);
1331 BLSURFPlugin_Hypothesis::TFaceEntryEnfVertexListMap::const_iterator enfIt = entryEnfVertexListMap.begin();
1332 for ( ; enfIt != entryEnfVertexListMap.end(); ++enfIt ) {
1333 if ( !enfIt->second.empty() ) {
1334 GeomShape = entryToShape(enfIt->first);
1335 if ( GeomShape.IsNull() )
1337 createEnforcedVertexOnFace( GeomShape, enfIt->second );
1340 else if ( GeomShape.ShapeType() == TopAbs_COMPOUND)
1342 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1343 if (it.Value().ShapeType() == TopAbs_FACE){
1344 HasSizeMapOnFace = true;
1345 createEnforcedVertexOnFace(it.Value(), enfIt->second);
1349 else if ( GeomShape.ShapeType() == TopAbs_FACE)
1351 HasSizeMapOnFace = true;
1352 createEnforcedVertexOnFace(GeomShape, enfIt->second);
1357 // Internal vertices
1358 bool useInternalVertexAllFaces = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFaces(hyp);
1359 if (useInternalVertexAllFaces) {
1360 std::string grpName = BLSURFPlugin_Hypothesis::GetInternalEnforcedVertexAllFacesGroup(hyp);
1362 TopExp_Explorer exp (theGeomShape, TopAbs_FACE);
1363 for (; exp.More(); exp.Next()){
1364 TopExp_Explorer exp_face (exp.Current(), TopAbs_VERTEX, TopAbs_EDGE);
1365 for (; exp_face.More(); exp_face.Next())
1367 // Get coords of vertex
1368 // Check if current coords is already in enfVertexList
1369 // If coords not in enfVertexList, add new enfVertex
1370 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(exp_face.Current()));
1371 BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex = new BLSURFPlugin_Hypothesis::TEnfVertex();
1372 enfVertex->coords.push_back(aPnt.X());
1373 enfVertex->coords.push_back(aPnt.Y());
1374 enfVertex->coords.push_back(aPnt.Z());
1375 enfVertex->name = "";
1376 enfVertex->faceEntries.clear();
1377 enfVertex->geomEntry = "";
1378 enfVertex->grpName = grpName;
1379 enfVertex->vertex = TopoDS::Vertex( exp_face.Current() );
1380 _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()), aPnt, enfVertex);
1381 HasSizeMapOnFace = true;
1386 cadsurf_set_sizemap_iso_cad_face(css, size_on_surface, &_smp_phy_size);
1388 if (HasSizeMapOnEdge){
1389 cadsurf_set_sizemap_iso_cad_edge(css, size_on_edge, &_smp_phy_size);
1391 if (HasSizeMapOnVertex){
1392 cadsurf_set_sizemap_iso_cad_point(css, size_on_vertex, &_smp_phy_size);
1399 _preCadFacesIDsPeriodicityVector.clear();
1400 _preCadEdgesIDsPeriodicityVector.clear();
1402 for (std::size_t i = 0; i<preCadFacesPeriodicityVector.size(); i++){
1403 createPreCadFacesPeriodicity(theGeomShape, preCadFacesPeriodicityVector[i]);
1406 const BLSURFPlugin_Hypothesis::TPreCadPeriodicityVector preCadEdgesPeriodicityVector = BLSURFPlugin_Hypothesis::GetPreCadEdgesPeriodicityVector(hyp);
1408 for (std::size_t i = 0; i<preCadEdgesPeriodicityVector.size(); i++){
1409 createPreCadEdgesPeriodicity(theGeomShape, preCadEdgesPeriodicityVector[i]);
1413 //================================================================================
1415 * \brief Throws an exception if a parameter name is wrong
1417 //================================================================================
1419 void BLSURFPlugin_BLSURF::set_param(cadsurf_session_t *css,
1420 const char * option_name,
1421 const char * option_value)
1423 status_t status = cadsurf_set_param(css, option_name, option_value );
1425 if ( _hypothesis && _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
1426 cout << option_name << " = " << option_value << endl;
1428 if ( status != MESHGEMS_STATUS_OK )
1430 if ( status == MESHGEMS_STATUS_UNKNOWN_PARAMETER ) {
1431 throw SALOME_Exception
1432 ( SMESH_Comment("Invalid name of CADSURF parameter: ") << option_name );
1434 else if ( status == MESHGEMS_STATUS_NOLICENSE )
1435 throw SALOME_Exception
1436 ( "No valid license available" );
1438 throw SALOME_Exception
1439 ( SMESH_Comment("Either wrong name or unacceptable value of CADSURF parameter '")
1440 << option_name << "': " << option_value);
1446 // --------------------------------------------------------------------------
1448 * \brief Class correctly terminating usage of MG-CADSurf library at destruction
1450 struct BLSURF_Cleaner
1453 cadsurf_session_t* _css;
1457 BLSURF_Cleaner(context_t * ctx,
1458 cadsurf_session_t* css=0,
1469 Clean( /*exceptContext=*/false );
1471 void Clean(const bool exceptContext)
1475 cadsurf_session_delete(_css); _css = 0;
1477 // #if BLSURF_VERSION_LONG >= "3.1.1"
1478 // // if(geo_sizemap_e)
1479 // // distene_sizemap_delete(geo_sizemap_e);
1480 // // if(geo_sizemap_f)
1481 // // distene_sizemap_delete(geo_sizemap_f);
1482 // if(iso_sizemap_p)
1483 // distene_sizemap_delete(iso_sizemap_p);
1484 // if(iso_sizemap_e)
1485 // distene_sizemap_delete(iso_sizemap_e);
1486 // if(iso_sizemap_f)
1487 // distene_sizemap_delete(iso_sizemap_f);
1489 // // if(clean_geo_sizemap_e)
1490 // // distene_sizemap_delete(clean_geo_sizemap_e);
1491 // // if(clean_geo_sizemap_f)
1492 // // distene_sizemap_delete(clean_geo_sizemap_f);
1493 // if(clean_iso_sizemap_p)
1494 // distene_sizemap_delete(clean_iso_sizemap_p);
1495 // if(clean_iso_sizemap_e)
1496 // distene_sizemap_delete(clean_iso_sizemap_e);
1497 // if(clean_iso_sizemap_f)
1498 // distene_sizemap_delete(clean_iso_sizemap_f);
1501 cad_delete(_cad); _cad = 0;
1502 dcad_delete(_dcad); _dcad = 0;
1503 if ( !exceptContext )
1505 context_delete(_ctx); _ctx = 0;
1511 // --------------------------------------------------------------------------
1512 // comparator to sort nodes and sub-meshes
1513 struct ShapeTypeCompare
1515 // sort nodes by position in the following order:
1516 // SMDS_TOP_FACE=2, SMDS_TOP_EDGE=1, SMDS_TOP_VERTEX=0, SMDS_TOP_3DSPACE=3
1517 bool operator()( const SMDS_MeshNode* n1, const SMDS_MeshNode* n2 ) const
1519 // NEW ORDER: nodes earlier added to sub-mesh are considered "less"
1520 //return n1->getIdInShape() < n2->getIdInShape();
1521 return n1->GetID() < n2->GetID(); // earlier created nodes have less IDs
1523 // sort sub-meshes in order: EDGE, VERTEX
1524 bool operator()( const SMESHDS_SubMesh* s1, const SMESHDS_SubMesh* s2 ) const
1526 int isVertex1 = ( s1 && s1->NbElements() == 0 );
1527 int isVertex2 = ( s2 && s2->NbElements() == 0 );
1528 if ( isVertex1 == isVertex2 )
1530 return isVertex1 < isVertex2;
1534 //================================================================================
1536 * \brief Fills groups of nodes to be merged
1538 //================================================================================
1540 void getNodeGroupsToMerge( const SMESHDS_SubMesh* smDS,
1541 const TopoDS_Shape& shape,
1542 SMESH_MeshEditor::TListOfListOfNodes& nodeGroupsToMerge)
1544 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1545 switch ( shape.ShapeType() )
1547 case TopAbs_VERTEX: {
1548 std::list< const SMDS_MeshNode* > nodes;
1549 while ( nIt->more() )
1550 nodes.push_back( nIt->next() );
1551 if ( nodes.size() > 1 )
1552 nodeGroupsToMerge.push_back( nodes );
1556 std::multimap< double, const SMDS_MeshNode* > u2node;
1557 while ( nIt->more() )
1559 const SMDS_MeshNode* n = nIt->next();
1560 if ( SMDS_EdgePositionPtr ePos = n->GetPosition() )
1561 u2node.insert( make_pair( ePos->GetUParameter(), n ));
1563 if ( u2node.size() < 2 ) return;
1565 //double tol = (( u2node.rbegin()->first - u2node.begin()->first ) / 20.) / u2node.size();
1567 BRep_Tool::Range( TopoDS::Edge( shape ), f,l );
1568 double tol = (( l - f ) / 10.) / u2node.size(); // 10. - adjusted for #17262
1570 std::multimap< double, const SMDS_MeshNode* >::iterator un2, un1;
1571 for ( un2 = u2node.begin(), un1 = un2++; un2 != u2node.end(); un1 = un2++ )
1573 if (( un2->first - un1->first ) <= tol )
1575 std::list< const SMDS_MeshNode* > nodes;
1576 nodes.push_back( un1->second );
1577 while (( un2->first - un1->first ) <= tol )
1579 nodes.push_back( un2->second );
1580 if ( ++un2 == u2node.end()) {
1585 // make nodes created on the boundary of viscous layer replace nodes created
1586 // by MG-CADSurf as their SMDS_Position is more correct
1587 nodes.sort( ShapeTypeCompare() );
1588 nodeGroupsToMerge.push_back( nodes );
1595 // SMESH_MeshEditor::TListOfListOfNodes::const_iterator nll = nodeGroupsToMerge.begin();
1596 // for ( ; nll != nodeGroupsToMerge.end(); ++nll )
1598 // cout << "Merge ";
1599 // const std::list< const SMDS_MeshNode* >& nl = *nll;
1600 // std::list< const SMDS_MeshNode* >::const_iterator nIt = nl.begin();
1601 // for ( ; nIt != nl.end(); ++nIt )
1602 // cout << (*nIt) << " ";
1608 //================================================================================
1610 * \brief A temporary mesh used to compute mesh on a proxy FACE
1612 //================================================================================
1614 struct TmpMesh: public SMESH_Mesh
1616 typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
1617 TN2NMap _tmp2origNN;
1618 TopoDS_Face _proxyFace;
1622 _myMeshDS = new SMESHDS_Mesh( _id, true );
1624 //--------------------------------------------------------------------------------
1626 * \brief Creates a FACE bound by viscous layers and mesh each its EDGE with 1 segment
1628 //--------------------------------------------------------------------------------
1630 const TopoDS_Face& makeProxyFace( SMESH_ProxyMesh::Ptr& viscousMesh,
1631 const TopoDS_Face& origFace)
1633 SMESH_Mesh* origMesh = viscousMesh->GetMesh();
1635 SMESH_MesherHelper helper( *origMesh );
1636 helper.SetSubShape( origFace );
1637 const bool hasSeam = helper.HasRealSeam();
1639 // get data of nodes on inner boundary of viscous layers
1641 TSideVector wireVec = StdMeshers_FaceSide::GetFaceWires(origFace, *origMesh,
1642 /*skipMediumNodes = */true,
1643 err, &helper, viscousMesh );
1644 if ( err && err->IsKO() )
1645 throw *err.get(); // it should be caught at SMESH_subMesh
1647 // proxy nodes and corresponding tmp VERTEXes
1648 std::vector<const SMDS_MeshNode*> origNodes;
1649 std::vector<TopoDS_Vertex> tmpVertex;
1651 // create a proxy FACE
1652 TopoDS_Face origFaceCopy = TopoDS::Face( origFace.EmptyCopied() );
1653 BRepBuilderAPI_MakeFace newFace( origFaceCopy );
1654 bool hasPCurves = false;
1655 for ( size_t iW = 0; iW != wireVec.size(); ++iW )
1657 StdMeshers_FaceSidePtr& wireData = wireVec[iW];
1658 const UVPtStructVec& wirePoints = wireData->GetUVPtStruct();
1659 if ( wirePoints.size() < 3 )
1662 BRepBuilderAPI_MakePolygon polygon;
1663 const size_t i0 = tmpVertex.size();
1664 for ( size_t iN = 0; iN < wirePoints.size(); ++iN )
1666 polygon.Add( SMESH_TNodeXYZ( wirePoints[ iN ].node ));
1667 origNodes.push_back( wirePoints[ iN ].node );
1668 tmpVertex.push_back( polygon.LastVertex() );
1670 // check presence of a pcurve
1671 checkPCurve( polygon, origFaceCopy, hasPCurves, &wirePoints[ iN-1 ] );
1673 tmpVertex[ i0 ] = polygon.FirstVertex(); // polygon.LastVertex()==NULL for 1 vertex in wire
1675 if ( !polygon.IsDone() )
1676 throw SALOME_Exception("BLSURFPlugin_BLSURF: BRepBuilderAPI_MakePolygon failed");
1677 TopoDS_Wire wire = polygon;
1679 wire = updateSeam( wire, origNodes );
1680 newFace.Add( wire );
1682 _proxyFace = newFace;
1684 // set a new shape to mesh
1685 TopoDS_Compound auxCompoundToMesh;
1686 BRep_Builder shapeBuilder;
1687 shapeBuilder.MakeCompound( auxCompoundToMesh );
1688 shapeBuilder.Add( auxCompoundToMesh, _proxyFace );
1689 shapeBuilder.Add( auxCompoundToMesh, origMesh->GetShapeToMesh() );
1691 ShapeToMesh( auxCompoundToMesh );
1694 // Make input mesh for MG-CADSurf: segments on EDGE's of newFace
1696 // make nodes and fill in _tmp2origNN
1698 SMESHDS_Mesh* tmpMeshDS = GetMeshDS();
1699 for ( size_t i = 0; i < origNodes.size(); ++i )
1701 GetSubMesh( tmpVertex[i] )->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1702 if ( const SMDS_MeshNode* tmpN = SMESH_Algo::VertexNode( tmpVertex[i], tmpMeshDS ))
1703 _tmp2origNN.insert( _tmp2origNN.end(), make_pair( tmpN, origNodes[i] ));
1704 // else -- it can be a seam vertex replaced by updateSeam()
1705 // throw SALOME_Exception("BLSURFPlugin_BLSURF: a proxy vertex not meshed");
1709 TopoDS_Vertex v1, v2;
1710 for ( TopExp_Explorer edge( _proxyFace, TopAbs_EDGE ); edge.More(); edge.Next() )
1712 const TopoDS_Edge& E = TopoDS::Edge( edge.Current() );
1713 TopExp::Vertices( E, v1, v2 );
1714 const SMDS_MeshNode* n1 = SMESH_Algo::VertexNode( v1, tmpMeshDS );
1715 const SMDS_MeshNode* n2 = SMESH_Algo::VertexNode( v2, tmpMeshDS );
1717 if ( SMDS_MeshElement* seg = tmpMeshDS->AddEdge( n1, n2 ))
1718 tmpMeshDS->SetMeshElementOnShape( seg, E );
1724 //--------------------------------------------------------------------------------
1726 * \brief Add pcurve to the last edge of a wire
1728 //--------------------------------------------------------------------------------
1730 void checkPCurve( BRepBuilderAPI_MakePolygon& wire,
1731 const TopoDS_Face& face,
1733 const uvPtStruct * wirePoints )
1737 TopoDS_Edge edge = wire.Edge();
1738 if ( edge.IsNull() ) return;
1740 if ( BRep_Tool::CurveOnSurface(edge, face, f, l))
1745 gp_XY p1 = wirePoints[ 0 ].UV(), p2 = wirePoints[ 1 ].UV();
1746 Handle(Geom2d_Line) pcurve = new Geom2d_Line( p1, gp_Dir2d( p2 - p1 ));
1747 BRep_Builder().UpdateEdge( edge, Handle(Geom_Curve)(), Precision::Confusion() );
1748 BRep_Builder().UpdateEdge( edge, pcurve, face, Precision::Confusion() );
1749 BRep_Builder().Range( edge, 0, ( p2 - p1 ).Modulus() );
1750 // cout << "n1 = mesh.AddNode( " << p1.X()*10 << ", " << p1.Y() << ", 0 )" << endl
1751 // << "n2 = mesh.AddNode( " << p2.X()*10 << ", " << p2.Y() << ", 0 )" << endl
1752 // << "mesh.AddEdge( [ n1, n2 ] )" << endl;
1755 //--------------------------------------------------------------------------------
1757 * \brief Replace coincident EDGEs with reversed copies.
1759 //--------------------------------------------------------------------------------
1761 TopoDS_Wire updateSeam( const TopoDS_Wire& wire,
1762 const std::vector<const SMDS_MeshNode*>& nodesOfVertices )
1764 BRepBuilderAPI_MakeWire newWire;
1766 typedef NCollection_DataMap<SMESH_TLink, TopoDS_Edge, SMESH_TLink > TSeg2EdgeMap;
1767 TSeg2EdgeMap seg2EdgeMap;
1769 TopoDS_Iterator edgeIt( wire );
1770 for ( int iSeg = 1; edgeIt.More(); edgeIt.Next(), ++iSeg )
1772 SMESH_TLink link( nodesOfVertices[ iSeg-1 ], nodesOfVertices[ iSeg ]);
1773 TopoDS_Edge edge( TopoDS::Edge( edgeIt.Value() ));
1775 TopoDS_Edge* edgeInMap = seg2EdgeMap.Bound( link, edge );
1776 bool isSeam = ( *edgeInMap != edge );
1779 edgeInMap->Reverse();
1782 newWire.Add( edge );
1787 //--------------------------------------------------------------------------------
1789 * \brief Fill in the origMesh with faces computed by MG-CADSurf in this tmp mesh
1791 //--------------------------------------------------------------------------------
1793 void FillInOrigMesh( SMESH_Mesh& origMesh,
1794 const TopoDS_Face& origFace )
1796 SMESH_MesherHelper helper( origMesh );
1797 helper.SetSubShape( origFace );
1798 helper.SetElementsOnShape( true );
1800 SMESH_MesherHelper tmpHelper( *this );
1801 tmpHelper.SetSubShape( _proxyFace );
1803 // iterate over tmp faces and copy them in origMesh
1804 const SMDS_MeshNode* nodes[27];
1805 const SMDS_MeshNode* nullNode = 0;
1807 SMDS_FaceIteratorPtr fIt = GetMeshDS()->facesIterator();
1808 while ( fIt->more() )
1810 const SMDS_MeshElement* f = fIt->next();
1811 SMDS_ElemIteratorPtr nIt = f->nodesIterator();
1813 for ( ; nIt->more(); ++nbN )
1815 const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nIt->next() );
1816 TN2NMap::iterator n2nIt =
1817 _tmp2origNN.insert( _tmp2origNN.end(), make_pair( n, nullNode ));
1818 if ( !n2nIt->second ) {
1820 gp_XY uv = tmpHelper.GetNodeUV( _proxyFace, n );
1821 n2nIt->second = helper.AddNode( xyz[0], xyz[1], xyz[2], uv.X(), uv.Y() );
1823 nodes[ nbN ] = n2nIt->second;
1826 case 3: helper.AddFace( nodes[0], nodes[1], nodes[2] ); break;
1827 // case 6: helper.AddFace( nodes[0], nodes[1], nodes[2],
1828 // nodes[3], nodes[4], nodes[5]); break;
1829 case 4: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3] ); break;
1830 // case 9: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3],
1831 // nodes[4], nodes[5], nodes[6], nodes[7], nodes[8]); break;
1832 // case 8: helper.AddFace( nodes[0], nodes[1], nodes[2], nodes[3],
1833 // nodes[4], nodes[5], nodes[6], nodes[7]); break;
1840 * \brief A structure holding an error description and a verbisity level
1842 struct message_cb_user_data
1844 std::string * _error;
1851 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data);
1852 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
1853 real *duu, real *duv, real *dvv, void *user_data);
1854 status_t message_cb(message_t *msg, void *user_data);
1855 status_t interrupt_cb(integer *interrupt_status, void *user_data);
1857 //=============================================================================
1861 //=============================================================================
1863 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
1865 // Fix problem with locales
1866 Kernel_Utils::Localizer aLocalizer;
1868 this->SMESH_Algo::_progress = 1e-3; // prevent progress advancment while computing attractors
1870 bool viscousLayersMade =
1871 ( aShape.ShapeType() == TopAbs_FACE &&
1872 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( aShape ), aMesh ));
1874 if ( !viscousLayersMade )
1875 if ( !compute( aMesh, aShape, /*allowSubMeshClearing=*/true ))
1878 if ( _haveViscousLayers || viscousLayersMade )
1880 // Compute viscous layers
1882 TopTools_MapOfShape map;
1883 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1885 const TopoDS_Face& F = TopoDS::Face(face_iter.Current());
1886 if ( !map.Add( F )) continue;
1887 SMESH_ProxyMesh::Ptr viscousMesh = StdMeshers_ViscousLayers2D::Compute( aMesh, F );
1889 return false; // error in StdMeshers_ViscousLayers2D::Compute()
1891 // Compute MG-CADSurf mesh on viscous layers
1893 if ( viscousMesh->NbProxySubMeshes() > 0 )
1896 const TopoDS_Face& proxyFace = tmpMesh.makeProxyFace( viscousMesh, F );
1897 if ( !compute( tmpMesh, proxyFace, /*allowSubMeshClearing=*/false ))
1899 tmpMesh.FillInOrigMesh( aMesh, F );
1903 // Re-compute MG-CADSurf mesh on the rest faces if the mesh was cleared
1905 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
1907 const TopoDS_Face& F = TopoDS::Face(face_iter.Current());
1908 SMESH_subMesh* fSM = aMesh.GetSubMesh( F );
1909 if ( fSM->IsMeshComputed() ) continue;
1911 if ( !compute( aMesh, aShape, /*allowSubMeshClearing=*/true ))
1919 //=============================================================================
1923 //=============================================================================
1925 bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh& aMesh,
1926 const TopoDS_Shape& aShape,
1927 bool allowSubMeshClearing)
1929 /* create a distene context (generic object) */
1930 status_t status = STATUS_ERROR;
1932 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1933 SMESH_MesherHelper helper( aMesh ), helperWithShape( aMesh );
1934 myHelper = theHelper = & helperWithShape;
1935 // do not call helper.IsQuadraticSubMesh() because sub-meshes
1936 // may be cleaned and helper.myTLinkNodeMap gets invalid in such a case
1937 bool haveQuadraticSubMesh = helperWithShape.IsQuadraticSubMesh( aShape );
1938 bool quadraticSubMeshAndViscousLayer = false;
1939 bool needMerge = false;
1940 typedef set< SMESHDS_SubMesh*, ShapeTypeCompare > TSubMeshSet;
1941 TSubMeshSet edgeSubmeshes;
1942 TSubMeshSet& mergeSubmeshes = edgeSubmeshes;
1943 double existingPhySize = 0;
1945 TopTools_IndexedMapOfShape pmap, emap, fmap;
1947 TopTools_IndexedDataMapOfShapeListOfShape e2ffmap;
1948 TopExp::MapShapesAndAncestors( aShape, TopAbs_EDGE, TopAbs_FACE, e2ffmap );
1950 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
1952 feclearexcept( FE_ALL_EXCEPT );
1953 int oldFEFlags = fedisableexcept( FE_ALL_EXCEPT );
1956 context_t *ctx = context_new();
1958 /* Set the message callback in the working context */
1959 message_cb_user_data mcud;
1960 mcud._error = & this->SMESH_Algo::_comment;
1961 mcud._progress = & this->SMESH_Algo::_progress;
1963 _hypothesis ? _hypothesis->GetVerbosity() : BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
1964 context_set_message_callback(ctx, message_cb, &mcud);
1966 /* set the interruption callback */
1967 _compute_canceled = false;
1968 context_set_interrupt_callback(ctx, interrupt_cb, this);
1970 /* create the CAD object we will work on. It is associated to the context ctx. */
1971 cad_t *c = cad_new(ctx);
1972 dcad_t *dcad = dcad_new(c);
1974 // To enable multithreading
1975 cad_set_thread_safety(c, 1);
1977 /* Now fill the CAD object with data from your CAD
1978 * environement. This is the most complex part of a successfull
1984 cadsurf_session_t *css = cadsurf_session_new(ctx);
1986 // an object that correctly deletes all cadsurf objects at destruction
1987 BLSURF_Cleaner cleaner( ctx,css,c,dcad );
1989 SetParameters(_hypothesis, css, aShape);
1991 haveQuadraticSubMesh = haveQuadraticSubMesh || (_hypothesis != NULL && _hypothesis->GetQuadraticMesh());
1992 helper.SetIsQuadratic( haveQuadraticSubMesh );
1994 // To remove as soon as quadratic mesh is allowed - BEGIN
1995 // GDD: Viscous layer is not allowed with quadratic mesh
1996 if (_haveViscousLayers && haveQuadraticSubMesh ) {
1997 quadraticSubMeshAndViscousLayer = true;
1998 _haveViscousLayers = !haveQuadraticSubMesh;
1999 _comment += "Warning: Viscous layer is not possible with a quadratic mesh, it is ignored.";
2000 error(COMPERR_WARNING, _comment);
2002 // To remove as soon as quadratic mesh is allowed - END
2004 // needed to prevent the opencascade memory managmement from freeing things
2005 vector<Handle(Geom2d_Curve)> curves;
2006 vector<Handle(Geom_Surface)> surfaces;
2010 FaceId2PythonSmp.clear();
2011 EdgeId2PythonSmp.clear();
2012 VertexId2PythonSmp.clear();
2014 /****************************************************************************************
2016 *****************************************************************************************/
2018 string bad_end = "return";
2020 TopTools_IndexedMapOfShape _map;
2021 TopExp::MapShapes(aShape,TopAbs_VERTEX,_map);
2022 int ienf = _map.Extent();
2024 assert(Py_IsInitialized());
2025 PyGILState_STATE gstate;
2027 string theSizeMapStr;
2029 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
2031 TopoDS_Face f = TopoDS::Face(face_iter.Current());
2033 SMESH_subMesh* fSM = aMesh.GetSubMesh( f );
2034 if ( !fSM->IsEmpty() ) continue; // skip already meshed FACE with viscous layers
2036 // make INTERNAL face oriented FORWARD (issue 0020993)
2037 if (f.Orientation() != TopAbs_FORWARD && f.Orientation() != TopAbs_REVERSED )
2038 f.Orientation(TopAbs_FORWARD);
2040 iface = fmap.Add(f);
2042 surfaces.push_back(BRep_Tool::Surface(f));
2044 /* create an object representing the face for cadsurf */
2045 /* where face_id is an integer identifying the face.
2046 * surf_function is the function that defines the surface
2047 * (For this face, it will be called by cadsurf with your_face_object_ptr
2048 * as last parameter.
2050 cad_face_t *fce = cad_face_new(c, iface, surf_fun, surfaces.back().get());
2052 /* by default a face has no tag (color).
2053 The following call sets it to the same value as the Geom module ID : */
2054 int faceTag = meshDS->ShapeToIndex(f);
2055 faceTag = BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis );
2056 cad_face_set_tag(fce, faceTag);
2058 /* Set face orientation (optional if you want a well oriented output mesh)*/
2059 if(f.Orientation() != TopAbs_FORWARD)
2060 cad_face_set_orientation(fce, CAD_ORIENTATION_REVERSED);
2062 cad_face_set_orientation(fce, CAD_ORIENTATION_FORWARD);
2064 if (HasSizeMapOnFace /*&& !use_precad*/) //22903: use_precad seems not to interfere
2066 // -----------------
2068 // -----------------
2069 faceKey = FacesWithSizeMap.FindIndex(f);
2072 if (FaceId2SizeMap.find(faceKey)!=FaceId2SizeMap.end())
2074 theSizeMapStr = FaceId2SizeMap[faceKey];
2075 // check if function ends with "return"
2076 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2078 // Expr To Python function, verification is performed at validation in GUI
2079 gstate = PyGILState_Ensure();
2080 PyObject * obj = NULL;
2081 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2083 PyObject * func = NULL;
2084 func = PyObject_GetAttrString(main_mod, "f");
2085 FaceId2PythonSmp[iface]=func;
2086 FaceId2SizeMap.erase(faceKey);
2087 PyGILState_Release(gstate);
2090 // Specific size map = Attractor
2091 std::map<int,std::vector<double> >::iterator attractor_iter = FaceId2AttractorCoords.begin();
2093 for (; attractor_iter != FaceId2AttractorCoords.end(); ++attractor_iter) {
2094 if (attractor_iter->first == faceKey)
2096 double xyzCoords[3] = {attractor_iter->second[2],
2097 attractor_iter->second[3],
2098 attractor_iter->second[4]};
2100 gp_Pnt P(xyzCoords[0],xyzCoords[1],xyzCoords[2]);
2101 BRepClass_FaceClassifier scl(f,P,1e-7);
2102 scl.Perform(f, P, 1e-7);
2103 TopAbs_State result = scl.State();
2104 if ( result == TopAbs_OUT )
2105 MESSAGE("Point is out of face: node is not created");
2106 if ( result == TopAbs_UNKNOWN )
2107 MESSAGE("Point position on face is unknown: node is not created");
2108 if ( result == TopAbs_ON )
2109 MESSAGE("Point is on border of face: node is not created");
2110 if ( result == TopAbs_IN )
2112 // Point is inside face and not on border
2113 double uvCoords[2] = {attractor_iter->second[0],attractor_iter->second[1]};
2115 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
2116 cad_point_set_tag(point_p, ienf);
2118 FaceId2AttractorCoords.erase(faceKey);
2122 // -----------------
2124 // -----------------
2125 TId2ClsAttractorVec::iterator clAttractor_iter = FaceId2ClassAttractor.find(faceKey);
2126 if (clAttractor_iter != FaceId2ClassAttractor.end()){
2127 std::vector< BLSURFPlugin_Attractor* > & attVec = clAttractor_iter->second;
2128 for ( size_t i = 0; i < attVec.size(); ++i )
2129 if ( !attVec[i]->IsMapBuilt() ) {
2130 std::cout<<"Compute " << theNbAttractors-- << "-th attractor" <<std::endl;
2131 attVec[i]->BuildMap();
2133 FaceIndex2ClassAttractor[iface].swap( attVec );
2134 FaceId2ClassAttractor.erase(clAttractor_iter);
2136 } // if (HasSizeMapOnFace && !use_precad)
2138 // ------------------
2139 // Enforced Vertices
2140 // ------------------
2141 faceKey = FacesWithEnforcedVertices.FindIndex(f);
2142 std::map<int,BLSURFPlugin_Hypothesis::TEnfVertexCoordsList >::const_iterator evmIt = FaceId2EnforcedVertexCoords.find(faceKey);
2143 if (evmIt != FaceId2EnforcedVertexCoords.end())
2145 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList evl = evmIt->second;
2146 BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
2147 for (; evlIt != evl.end(); ++evlIt)
2149 double uvCoords[2] = { evlIt->at(0), evlIt->at(1) };
2151 cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
2153 BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
2154 xyzCoords.push_back(evlIt->at(2));
2155 xyzCoords.push_back(evlIt->at(3));
2156 xyzCoords.push_back(evlIt->at(4));
2157 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(xyzCoords);
2158 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end() &&
2159 !enfCoordsIt->second.empty() )
2161 // to merge nodes of an INTERNAL vertex belonging to several faces
2162 TopoDS_Vertex v = (*enfCoordsIt->second.begin() )->vertex;
2163 if ( v.IsNull() ) v = (*enfCoordsIt->second.rbegin())->vertex;
2164 if ( !v.IsNull() && meshDS->ShapeToIndex( v ) > 0 )
2166 tag = pmap.Add( v );
2167 SMESH_subMesh* vSM = aMesh.GetSubMesh( v );
2168 vSM->ComputeStateEngine( SMESH_subMesh::COMPUTE );
2169 mergeSubmeshes.insert( vSM->GetSubMeshDS() );
2170 // //if ( tag != pmap.Extent() )
2171 // needMerge = true;
2174 if ( tag == 0 ) tag = ienf;
2175 cad_point_set_tag(point_p, tag);
2177 FaceId2EnforcedVertexCoords.erase(faceKey);
2181 /****************************************************************************************
2183 now create the edges associated to this face
2184 *****************************************************************************************/
2186 for (TopExp_Explorer edge_iter(f,TopAbs_EDGE);edge_iter.More();edge_iter.Next())
2188 TopoDS_Edge e = TopoDS::Edge(edge_iter.Current());
2189 int ic = emap.FindIndex(e);
2194 curves.push_back(BRep_Tool::CurveOnSurface(e, f, tmin, tmax));
2196 if (HasSizeMapOnEdge){
2197 edgeKey = EdgesWithSizeMap.FindIndex(e);
2198 if (EdgeId2SizeMap.find(edgeKey)!=EdgeId2SizeMap.end())
2200 theSizeMapStr = EdgeId2SizeMap[edgeKey];
2201 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2203 // Expr To Python function, verification is performed at validation in GUI
2204 gstate = PyGILState_Ensure();
2205 PyObject * obj = NULL;
2206 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2208 PyObject * func = NULL;
2209 func = PyObject_GetAttrString(main_mod, "f");
2210 EdgeId2PythonSmp[ic]=func;
2211 EdgeId2SizeMap.erase(edgeKey);
2212 PyGILState_Release(gstate);
2215 /* data of nodes existing on the edge */
2216 StdMeshers_FaceSidePtr nodeData;
2217 SMESH_subMesh* sm = aMesh.GetSubMesh( e );
2218 if ( !sm->IsEmpty() )
2220 // SMESH_subMeshIteratorPtr subsmIt = sm->getDependsOnIterator( /*includeSelf=*/true,
2221 // /*complexFirst=*/false);
2222 // while ( subsmIt->more() )
2223 // edgeSubmeshes.insert( subsmIt->next()->GetSubMeshDS() );
2224 edgeSubmeshes.insert( sm->GetSubMeshDS() );
2226 nodeData.reset( new StdMeshers_FaceSide( f, e, &aMesh, /*isForwrd = */true,
2227 /*ignoreMedium=*/haveQuadraticSubMesh));
2228 if ( nodeData->MissVertexNode() )
2229 return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
2231 const std::vector<UVPtStruct>& nodeDataVec = nodeData->GetUVPtStruct();
2232 if ( !nodeDataVec.empty() )
2234 if ( Abs( nodeDataVec[0].param - tmin ) > Abs( nodeDataVec.back().param - tmin ))
2236 nodeData->Reverse();
2237 nodeData->GetUVPtStruct(); // nodeData recomputes nodeDataVec
2239 // tmin and tmax can change in case of viscous layer on an adjacent edge
2240 tmin = nodeDataVec.front().param;
2241 tmax = nodeDataVec.back().param;
2243 existingPhySize += nodeData->Length() / ( nodeDataVec.size() - 1 );
2247 cout << "---------------- Invalid nodeData" << endl;
2252 /* attach the edge to the current cadsurf face */
2253 cad_edge_t *edg = cad_edge_new(fce, ic, tmin, tmax, curv_fun, curves.back().get());
2255 /* by default an edge has no tag (color).
2256 The following call sets it to the same value as the edge_id : */
2257 // IMP23368. Do not set tag to an EDGE shared by FACEs of a hyper-patch
2258 bool isInHyperPatch = false;
2260 std::set< int > faceTags, faceIDs;
2261 TopTools_ListIteratorOfListOfShape fIt( e2ffmap.FindFromKey( e ));
2262 for ( ; fIt.More(); fIt.Next() )
2264 int faceTag = meshDS->ShapeToIndex( fIt.Value() );
2265 if ( !faceIDs.insert( faceTag ).second )
2266 continue; // a face encounters twice for a seam edge
2267 int hpTag = BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis );
2268 if ( !faceTags.insert( hpTag ).second )
2270 isInHyperPatch = true;
2275 if ( !isInHyperPatch )
2276 cad_edge_set_tag(edg, ic);
2278 /* by default, an edge does not necessalry appear in the resulting mesh,
2279 unless the following property is set :
2281 cad_edge_set_property(edg, EDGE_PROPERTY_SOFT_REQUIRED);
2283 /* by default an edge is a boundary edge */
2284 if (e.Orientation() == TopAbs_INTERNAL)
2285 cad_edge_set_property(edg, EDGE_PROPERTY_INTERNAL);
2287 // pass existing nodes of sub-meshes to MG-CADSurf
2290 const std::vector<UVPtStruct>& nodeDataVec = nodeData->GetUVPtStruct();
2291 const int nbNodes = nodeDataVec.size();
2293 dcad_edge_discretization_t *dedge;
2294 dcad_get_edge_discretization(dcad, edg, &dedge);
2295 dcad_edge_discretization_set_vertex_count( dedge, nbNodes );
2297 // cout << endl << " EDGE " << ic << endl;
2298 // cout << "tmin = "<<tmin << ", tmax = "<< tmax << endl;
2299 for ( int iN = 0; iN < nbNodes; ++iN )
2301 const UVPtStruct& nData = nodeDataVec[ iN ];
2302 double t = nData.param;
2303 real uv[2] = { nData.u, nData.v };
2304 SMESH_TNodeXYZ nXYZ( nData.node );
2305 // cout << "\tt = " << t
2306 // << "\t uv = ( " << uv[0] << ","<< uv[1] << " ) "
2307 // << "\t u = " << nData.param
2308 // << "\t ID = " << nData.node->GetID() << endl;
2309 dcad_edge_discretization_set_vertex_coordinates( dedge, iN+1, t, uv, nXYZ.ChangeData() );
2311 TopoDS_Shape v = helper.GetSubShapeByNode( nodeDataVec[0].node, meshDS );
2312 if ( !v.IsNull() && v.ShapeType() == TopAbs_VERTEX )
2313 dcad_edge_discretization_set_vertex_tag( dedge, 1, pmap.Add( v ));
2315 v = helper.GetSubShapeByNode( nodeDataVec.back().node, meshDS );
2316 if ( !v.IsNull() && v.ShapeType() == TopAbs_VERTEX )
2317 dcad_edge_discretization_set_vertex_tag( dedge, nbNodes, pmap.Add( v ));
2319 dcad_edge_discretization_set_property(dedge, DISTENE_DCAD_PROPERTY_REQUIRED);
2322 /****************************************************************************************
2324 *****************************************************************************************/
2328 gp_Pnt2d e0 = curves.back()->Value(tmin);
2329 gp_Pnt ee0 = surfaces.back()->Value(e0.X(), e0.Y());
2330 Standard_Real d1=0,d2=0;
2333 for (TopExp_Explorer ex_edge(e ,TopAbs_VERTEX); ex_edge.More(); ex_edge.Next()) {
2334 TopoDS_Vertex v = TopoDS::Vertex(ex_edge.Current());
2338 d1 = ee0.SquareDistance(BRep_Tool::Pnt(v));
2341 d2 = ee0.SquareDistance(BRep_Tool::Pnt(v));
2343 *ip = pmap.FindIndex(v);
2346 // SMESH_subMesh* sm = aMesh.GetSubMesh(v);
2347 // if ( sm->IsMeshComputed() )
2348 // edgeSubmeshes.insert( sm->GetSubMeshDS() );
2351 // std::string aFileName = "fmap_vertex_";
2352 // aFileName.append(val_to_string(*ip));
2353 // aFileName.append(".brep");
2354 // BRepTools::Write(v,aFileName.c_str());
2356 if (HasSizeMapOnVertex){
2357 vertexKey = VerticesWithSizeMap.FindIndex(v);
2358 if (VertexId2SizeMap.find(vertexKey)!=VertexId2SizeMap.end()){
2359 theSizeMapStr = VertexId2SizeMap[vertexKey];
2360 if (theSizeMapStr.find(bad_end) == (theSizeMapStr.size()-bad_end.size()-1))
2362 // Expr To Python function, verification is performed at validation in GUI
2363 gstate = PyGILState_Ensure();
2364 PyObject * obj = NULL;
2365 obj= PyRun_String(theSizeMapStr.c_str(), Py_file_input, main_dict, NULL);
2367 PyObject * func = NULL;
2368 func = PyObject_GetAttrString(main_mod, "f");
2369 VertexId2PythonSmp[*ip]=func;
2370 VertexId2SizeMap.erase(vertexKey); // do not erase if using a vector
2371 PyGILState_Release(gstate);
2376 // should not happen
2377 MESSAGE("An edge does not have 2 extremities.");
2380 // This defines the curves extremity connectivity
2381 cad_edge_set_extremities(edg, ip1, ip2);
2382 /* set the tag (color) to the same value as the extremity id : */
2383 cad_edge_set_extremities_tag(edg, ip1, ip2);
2386 cad_edge_set_extremities(edg, ip2, ip1);
2387 cad_edge_set_extremities_tag(edg, ip2, ip1);
2393 ///////////////////////
2395 ///////////////////////
2397 if (! _preCadFacesIDsPeriodicityVector.empty())
2399 for (std::size_t i=0; i < _preCadFacesIDsPeriodicityVector.size(); i++){
2400 std::vector<int> theFace1_ids = _preCadFacesIDsPeriodicityVector[i].shape1IDs;
2401 std::vector<int> theFace2_ids = _preCadFacesIDsPeriodicityVector[i].shape2IDs;
2402 int* theFace1_ids_c = &theFace1_ids[0];
2403 int* theFace2_ids_c = &theFace2_ids[0];
2404 std::ostringstream o;
2405 o << "_preCadFacesIDsPeriodicityVector[" << i << "] = [";
2406 for (std::size_t j=0; j < theFace1_ids.size(); j++)
2407 o << theFace1_ids[j] << ", ";
2409 for (std::size_t j=0; j < theFace2_ids.size(); j++)
2410 o << theFace2_ids[j] << ", ";
2412 // if ( _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
2413 // cout << o.str() << endl;
2414 if (_preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords.empty())
2416 // If no source points, call periodicity without transformation function
2417 meshgems_cad_periodicity_transformation_t periodicity_transformation = NULL;
2418 status = cad_add_face_multiple_periodicity_with_transformation_function(c, theFace1_ids_c, theFace1_ids.size(),
2419 theFace2_ids_c, theFace2_ids.size(), periodicity_transformation, NULL);
2420 if(status != STATUS_OK)
2421 cout << "cad_add_face_multiple_periodicity_with_transformation_function failed with error code " << status << "\n";
2425 // get the transformation vertices
2426 double* theSourceVerticesCoords_c = &_preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords[0];
2427 double* theTargetVerticesCoords_c = &_preCadFacesIDsPeriodicityVector[i].theTargetVerticesCoords[0];
2428 int nbSourceVertices = _preCadFacesIDsPeriodicityVector[i].theSourceVerticesCoords.size()/3;
2429 int nbTargetVertices = _preCadFacesIDsPeriodicityVector[i].theTargetVerticesCoords.size()/3;
2431 status = cad_add_face_multiple_periodicity_with_transformation_function_by_points(c, theFace1_ids_c, theFace1_ids.size(),
2432 theFace2_ids_c, theFace2_ids.size(), theSourceVerticesCoords_c, nbSourceVertices, theTargetVerticesCoords_c, nbTargetVertices);
2433 if(status != STATUS_OK)
2434 cout << "cad_add_face_multiple_periodicity_with_transformation_function_by_points failed with error code " << status << "\n";
2439 if (! _preCadEdgesIDsPeriodicityVector.empty())
2441 for (std::size_t i=0; i < _preCadEdgesIDsPeriodicityVector.size(); i++){
2442 std::vector<int> theEdge1_ids = _preCadEdgesIDsPeriodicityVector[i].shape1IDs;
2443 std::vector<int> theEdge2_ids = _preCadEdgesIDsPeriodicityVector[i].shape2IDs;
2444 // Use the address of the first element of the vector to initialize the array
2445 int* theEdge1_ids_c = &theEdge1_ids[0];
2446 int* theEdge2_ids_c = &theEdge2_ids[0];
2448 std::ostringstream o;
2449 o << "_preCadEdgesIDsPeriodicityVector[" << i << "] = [";
2450 for (std::size_t j=0; j < theEdge1_ids.size(); j++)
2451 o << theEdge1_ids[j] << ", ";
2453 for (std::size_t j=0; j < theEdge2_ids.size(); j++)
2454 o << theEdge2_ids[j] << ", ";
2456 // if ( _hypothesis->GetVerbosity() > _hypothesis->GetDefaultVerbosity() )
2457 // cout << o.str() << endl;
2459 if (_preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords.empty())
2461 // If no source points, call periodicity without transformation function
2462 meshgems_cad_periodicity_transformation_t periodicity_transformation = NULL;
2463 status = cad_add_edge_multiple_periodicity_with_transformation_function(c, theEdge1_ids_c, theEdge1_ids.size(),
2464 theEdge2_ids_c, theEdge2_ids.size(), periodicity_transformation, NULL);
2465 if(status != STATUS_OK)
2466 cout << "cad_add_edge_multiple_periodicity_with_transformation_function failed with error code " << status << "\n";
2470 // get the transformation vertices
2471 double* theSourceVerticesCoords_c = &_preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords[0];
2472 double* theTargetVerticesCoords_c = &_preCadEdgesIDsPeriodicityVector[i].theTargetVerticesCoords[0];
2473 int nbSourceVertices = _preCadEdgesIDsPeriodicityVector[i].theSourceVerticesCoords.size()/3;
2474 int nbTargetVertices = _preCadEdgesIDsPeriodicityVector[i].theTargetVerticesCoords.size()/3;
2476 status = cad_add_edge_multiple_periodicity_with_transformation_function_by_points(c, theEdge1_ids_c, theEdge1_ids.size(),
2477 theEdge2_ids_c, theEdge2_ids.size(), theSourceVerticesCoords_c, nbSourceVertices, theTargetVerticesCoords_c, nbTargetVertices);
2478 if(status != STATUS_OK)
2479 cout << "cad_add_edge_multiple_periodicity_with_transformation_function_by_points failed with error code " << status << "\n";
2484 if ( !_hypothesis && !edgeSubmeshes.empty() && existingPhySize != 0 )
2486 // prevent failure due to the default PhySize incompatible with size of existing 1D mesh
2487 // and with face size
2488 // double minFaceSize = existingPhySize / edgeSubmeshes.size();
2489 // for ( int iF = 1; iF <= fmap.Extent(); ++iF )
2492 // BRepBndLib::Add( fmap( iF ), box );
2493 // gp_XYZ delta = box.CornerMax().XYZ() - box.CornerMin().XYZ();
2494 // std::sort( delta.ChangeData(), delta.ChangeData() + 3 );
2495 // minFaceSize = Min( minFaceSize, delta.Coord(2) );
2497 // set_param(css, "global_physical_size", val_to_string( minFaceSize * 0.5 ).c_str());
2498 // set_param(css, "max_size", val_to_string( minFaceSize * 5 ).c_str());
2501 // TODO: be able to use a mesh in input.
2502 // See imsh usage in Products/templates/mg-cadsurf_template_common.cpp
2503 // => cadsurf_set_mesh
2505 // Use the original dcad
2506 cadsurf_set_dcad(css, dcad);
2508 // Use the original cad
2509 cadsurf_set_cad(css, c);
2511 std::cout << std::endl;
2512 std::cout << "Beginning of Surface Mesh generation" << std::endl;
2513 std::cout << std::endl;
2518 status = cadsurf_compute_mesh(css);
2521 catch ( std::exception& exc ) {
2522 _comment += exc.what();
2524 catch (Standard_Failure& ex) {
2525 _comment += ex.DynamicType()->Name();
2526 if ( ex.GetMessageString() && strlen( ex.GetMessageString() )) {
2528 _comment += ex.GetMessageString();
2532 if ( _comment.empty() )
2533 _comment = "Exception in cadsurf_compute_mesh()";
2536 std::cout << std::endl;
2537 std::cout << "End of Surface Mesh generation" << std::endl;
2538 std::cout << std::endl;
2541 cadsurf_get_mesh(css, &msh);
2542 if ( !msh || STATUS_IS_ERROR( status ))
2544 /* release the mesh object */
2545 cadsurf_regain_mesh(css, msh);
2546 return error(_comment);
2549 // Clear mesh from already meshed edges if possible else
2550 // remember that merge is needed
2551 TSubMeshSet::iterator smIt = edgeSubmeshes.begin();
2552 for ( ; smIt != edgeSubmeshes.end(); ++smIt ) // loop on already meshed EDGEs
2554 SMESHDS_SubMesh* smDS = *smIt;
2555 if ( !smDS ) continue;
2556 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2559 const SMDS_MeshNode* n = nIt->next();
2560 if ( n->NbInverseElements( SMDSAbs_Face ) > 0 )
2562 needMerge = true; // to correctly sew with viscous mesh
2563 // add existing medium nodes to helper
2564 if ( aMesh.NbEdges( ORDER_QUADRATIC ) > 0 )
2566 SMDS_ElemIteratorPtr edgeIt = smDS->GetElements();
2567 while ( edgeIt->more() )
2568 helper.AddTLinks( static_cast<const SMDS_MeshEdge*>(edgeIt->next()));
2573 if ( allowSubMeshClearing )
2575 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
2576 while ( eIt->more() ) meshDS->RemoveFreeElement( eIt->next(), 0 );
2577 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2578 while ( nIt->more() ) meshDS->RemoveFreeNode( nIt->next(), 0 );
2587 std::string GMFFileName = BLSURFPlugin_Hypothesis::GetDefaultGMFFile();
2589 GMFFileName = _hypothesis->GetGMFFile();
2590 if (GMFFileName != "") {
2591 bool asciiFound = (GMFFileName.find(".mesh", GMFFileName.length()-5) != std::string::npos);
2592 bool binaryFound = (GMFFileName.find(".meshb",GMFFileName.length()-6) != std::string::npos);
2593 if (!asciiFound && !binaryFound)
2594 GMFFileName.append(".mesh");
2595 mesh_write_mesh(msh, GMFFileName.c_str());
2598 /* retrieve mesh data (see meshgems/mesh.h) */
2599 integer nv, ne, nt, nq, vtx[4], tag, nb_tag;
2600 integer *evedg, *evtri, *evquad, *tags_buff, type;
2603 mesh_get_vertex_count(msh, &nv);
2604 mesh_get_edge_count(msh, &ne);
2605 mesh_get_triangle_count(msh, &nt);
2606 mesh_get_quadrangle_count(msh, &nq);
2608 evedg = (integer *)mesh_calloc_generic_buffer(msh);
2609 evtri = (integer *)mesh_calloc_generic_buffer(msh);
2610 evquad = (integer *)mesh_calloc_generic_buffer(msh);
2611 tags_buff = (integer*)mesh_calloc_generic_buffer(msh);
2613 std::vector<const SMDS_MeshNode*> nodes(nv+1);
2614 std::vector<bool> tags(nv+1);
2616 /* enumerated vertices */
2617 for(int iv=1;iv<=nv;iv++) {
2618 mesh_get_vertex_coordinates(msh, iv, xyz);
2619 mesh_get_vertex_tag(msh, iv, &tag);
2620 // Issue 0020656. Use vertex coordinates
2622 if ( tag > 0 && tag <= pmap.Extent() ) {
2623 TopoDS_Vertex v = TopoDS::Vertex(pmap(tag));
2624 double tol = BRep_Tool::Tolerance( v );
2625 gp_Pnt p = BRep_Tool::Pnt( v );
2626 if ( p.IsEqual( gp_Pnt( xyz[0], xyz[1], xyz[2]), 1e3*tol))
2627 xyz[0] = p.X(), xyz[1] = p.Y(), xyz[2] = p.Z();
2629 tag = 0; // enforced or attracted vertex
2630 nodes[iv] = SMESH_Algo::VertexNode( v, meshDS );
2633 nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
2635 // Create group of enforced vertices if requested
2636 BLSURFPlugin_Hypothesis::TEnfVertexCoords projVertex;
2638 projVertex.push_back((double)xyz[0]);
2639 projVertex.push_back((double)xyz[1]);
2640 projVertex.push_back((double)xyz[2]);
2641 std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(projVertex);
2642 if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end())
2644 BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfListIt = enfCoordsIt->second.begin();
2645 BLSURFPlugin_Hypothesis::TEnfVertex *currentEnfVertex;
2646 for (; enfListIt != enfCoordsIt->second.end(); ++enfListIt) {
2647 currentEnfVertex = (*enfListIt);
2648 if (currentEnfVertex->grpName != "") {
2649 bool groupDone = false;
2650 SMESH_Mesh::GroupIteratorPtr grIt = aMesh.GetGroups();
2651 while (grIt->more()) {
2652 SMESH_Group * group = grIt->next();
2653 if ( !group ) continue;
2654 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
2655 if ( !groupDS ) continue;
2656 if ( groupDS->GetType()==SMDSAbs_Node && currentEnfVertex->grpName.compare(group->GetName())==0) {
2657 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
2658 aGroupDS->SMDSGroup().Add(nodes[iv]);
2659 // How can I inform the hypothesis ?
2660 // _hypothesis->AddEnfVertexNodeID(currentEnfVertex->grpName,nodes[iv]->GetID());
2667 SMESH_Group* aGroup = aMesh.AddGroup( SMDSAbs_Node, currentEnfVertex->grpName.c_str() );
2668 aGroup->SetName( currentEnfVertex->grpName.c_str() );
2669 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
2670 aGroupDS->SMDSGroup().Add(nodes[iv]);
2674 throw SALOME_Exception(LOCALIZED("An enforced vertex node was not added to a group"));
2677 MESSAGE("Group name is empty: '"<<currentEnfVertex->grpName<<"' => group is not created");
2681 // internal points are tagged to zero
2682 if(tag > 0 && tag <= pmap.Extent() ){
2683 meshDS->SetNodeOnVertex(nodes[iv], TopoDS::Vertex(pmap(tag)));
2690 /* enumerate edges */
2691 for(int it=1;it<=ne;it++) {
2693 mesh_get_edge_vertices(msh, it, vtx);
2694 mesh_get_edge_extra_vertices(msh, it, &type, evedg);
2695 mesh_get_edge_tag(msh, it, &tag);
2697 // If PreCAD performed some cleaning operations (remove tiny edges,
2698 // merge edges ...) an output tag can indeed represent several original tags.
2699 // Get the initial tags corresponding to the output tag and redefine the tag as
2700 // the last of the two initial tags (else the output tag is out of emap and hasn't any meaning)
2701 mesh_get_composite_tag_definition(msh, tag, &nb_tag, tags_buff);
2703 tag=tags_buff[nb_tag-1];
2704 if ( tag < 1 || tag > emap.Extent() )
2706 std::cerr << "MG-CADSurf BUG:::: Edge tag " << tag
2707 << " does not point to a CAD edge (nb edges " << emap.Extent() << ")" << std::endl;
2711 Set_NodeOnEdge(meshDS, nodes[vtx[0]], emap(tag));
2712 tags[vtx[0]] = false;
2715 Set_NodeOnEdge(meshDS, nodes[vtx[1]], emap(tag));
2716 tags[vtx[1]] = false;
2718 if (type == MESHGEMS_MESH_ELEMENT_TYPE_EDGE3) {
2720 if (tags[evedg[0]]) {
2721 Set_NodeOnEdge(meshDS, nodes[evedg[0]], emap(tag));
2722 tags[evedg[0]] = false;
2724 edg = meshDS->AddEdge(nodes[vtx[0]], nodes[vtx[1]], nodes[evedg[0]]);
2727 edg = helper.AddEdge(nodes[vtx[0]], nodes[vtx[1]]);
2729 meshDS->SetMeshElementOnShape(edg, TopoDS::Edge(emap(tag)));
2732 /* enumerate triangles */
2733 for(int it=1;it<=nt;it++) {
2735 mesh_get_triangle_vertices(msh, it, vtx);
2736 mesh_get_triangle_extra_vertices(msh, it, &type, evtri);
2737 mesh_get_triangle_tag(msh, it, &tag);
2739 meshDS->SetNodeOnFace(nodes[vtx[0]], tag);
2740 tags[vtx[0]] = false;
2743 meshDS->SetNodeOnFace(nodes[vtx[1]], tag);
2744 tags[vtx[1]] = false;
2747 meshDS->SetNodeOnFace(nodes[vtx[2]], tag);
2748 tags[vtx[2]] = false;
2750 if (type == MESHGEMS_MESH_ELEMENT_TYPE_TRIA6) {
2751 // QUADRATIC TRIANGLE
2752 if (tags[evtri[0]]) {
2753 meshDS->SetNodeOnFace(nodes[evtri[0]], tag);
2754 tags[evtri[0]] = false;
2756 if (tags[evtri[1]]) {
2757 meshDS->SetNodeOnFace(nodes[evtri[1]], tag);
2758 tags[evtri[1]] = false;
2760 if (tags[evtri[2]]) {
2761 meshDS->SetNodeOnFace(nodes[evtri[2]], tag);
2762 tags[evtri[2]] = false;
2764 tri = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]],
2765 nodes[evtri[0]], nodes[evtri[1]], nodes[evtri[2]]);
2768 if ( helper.GetIsQuadratic() )
2769 helper.SetSubShape( tag );
2770 tri = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]]);
2772 meshDS->SetMeshElementOnShape(tri, tag);
2775 /* enumerate quadrangles */
2776 for(int it=1;it<=nq;it++) {
2777 SMDS_MeshFace* quad;
2778 mesh_get_quadrangle_vertices(msh, it, vtx);
2779 mesh_get_quadrangle_extra_vertices(msh, it, &type, evquad);
2780 mesh_get_quadrangle_tag(msh, it, &tag);
2782 meshDS->SetNodeOnFace(nodes[vtx[0]], tag);
2783 tags[vtx[0]] = false;
2786 meshDS->SetNodeOnFace(nodes[vtx[1]], tag);
2787 tags[vtx[1]] = false;
2790 meshDS->SetNodeOnFace(nodes[vtx[2]], tag);
2791 tags[vtx[2]] = false;
2794 meshDS->SetNodeOnFace(nodes[vtx[3]], tag);
2795 tags[vtx[3]] = false;
2797 if (type == MESHGEMS_MESH_ELEMENT_TYPE_QUAD9) {
2798 // QUADRATIC QUADRANGLE
2799 if (tags[evquad[0]]) {
2800 meshDS->SetNodeOnFace(nodes[evquad[0]], tag);
2801 tags[evquad[0]] = false;
2803 if (tags[evquad[1]]) {
2804 meshDS->SetNodeOnFace(nodes[evquad[1]], tag);
2805 tags[evquad[1]] = false;
2807 if (tags[evquad[2]]) {
2808 meshDS->SetNodeOnFace(nodes[evquad[2]], tag);
2809 tags[evquad[2]] = false;
2811 if (tags[evquad[3]]) {
2812 meshDS->SetNodeOnFace(nodes[evquad[3]], tag);
2813 tags[evquad[3]] = false;
2815 if (tags[evquad[4]]) {
2816 meshDS->SetNodeOnFace(nodes[evquad[4]], tag);
2817 tags[evquad[4]] = false;
2819 quad = meshDS->AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]],
2820 nodes[evquad[0]], nodes[evquad[1]], nodes[evquad[2]], nodes[evquad[3]],
2824 quad = helper.AddFace(nodes[vtx[0]], nodes[vtx[1]], nodes[vtx[2]], nodes[vtx[3]]);
2826 meshDS->SetMeshElementOnShape(quad, tag);
2829 /* release the mesh object, the rest is released by cleaner */
2830 cadsurf_regain_mesh(css, msh);
2833 // Remove free nodes that can appear e.g. if "remove tiny edges"(IPAL53235)
2834 for(int iv=1;iv<=nv;iv++)
2835 if ( nodes[iv] && nodes[iv]->NbInverseElements() == 0 )
2836 meshDS->RemoveFreeNode( nodes[iv], 0, /*fromGroups=*/false );
2839 if ( needMerge ) // sew mesh computed by MG-CADSurf with pre-existing mesh
2841 SMESH_MeshEditor editor( &aMesh );
2842 SMESH_MeshEditor::TListOfListOfNodes nodeGroupsToMerge;
2843 TIDSortedElemSet segementsOnEdge;
2844 TSubMeshSet::iterator smIt;
2845 SMESHDS_SubMesh* smDS;
2847 // merge nodes on EDGE's with ones computed by MG-CADSurf
2848 for ( smIt = mergeSubmeshes.begin(); smIt != mergeSubmeshes.end(); ++smIt )
2850 if (! (smDS = *smIt) ) continue;
2851 getNodeGroupsToMerge( smDS, meshDS->IndexToShape((*smIt)->GetID()), nodeGroupsToMerge );
2853 SMDS_ElemIteratorPtr segIt = smDS->GetElements();
2854 while ( segIt->more() )
2855 segementsOnEdge.insert( segIt->next() );
2858 editor.MergeNodes( nodeGroupsToMerge );
2861 SMESH_MeshEditor::TListOfListOfElementsID equalSegments;
2862 editor.FindEqualElements( segementsOnEdge, equalSegments );
2863 editor.MergeElements( equalSegments );
2865 // remove excess segments created on the boundary of viscous layers
2866 const SMDS_TypeOfPosition onFace = SMDS_TOP_FACE;
2867 for ( int i = 1; i <= emap.Extent(); ++i )
2869 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( emap( i )))
2871 SMDS_ElemIteratorPtr segIt = smDS->GetElements();
2872 while ( segIt->more() )
2874 const SMDS_MeshElement* seg = segIt->next();
2875 if ( seg->GetNode(0)->GetPosition()->GetTypeOfPosition() == onFace ||
2876 seg->GetNode(1)->GetPosition()->GetTypeOfPosition() == onFace )
2877 meshDS->RemoveFreeElement( seg, smDS );
2884 // SetIsAlwaysComputed( true ) to sub-meshes of EDGEs w/o mesh
2885 for (int i = 1; i <= emap.Extent(); i++)
2886 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( emap( i )))
2887 sm->SetIsAlwaysComputed( true );
2888 for (int i = 1; i <= pmap.Extent(); i++)
2889 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( pmap( i )))
2890 if ( !sm->IsMeshComputed() )
2891 sm->SetIsAlwaysComputed( true );
2893 // Set error to FACE's w/o elements
2894 SMESH_ComputeErrorName err = COMPERR_ALGO_FAILED;
2895 if ( _comment.empty() && status == STATUS_OK )
2897 err = COMPERR_WARNING;
2898 _comment = "No mesh elements assigned to a face";
2900 bool badFaceFound = false;
2901 for (TopExp_Explorer face_iter(aShape,TopAbs_FACE);face_iter.More();face_iter.Next())
2903 TopoDS_Face f = TopoDS::Face(face_iter.Current());
2904 SMESH_subMesh* sm = aMesh.GetSubMesh( f );
2905 if ( !sm->GetSubMeshDS() || sm->GetSubMeshDS()->NbElements() == 0 )
2907 int faceTag = sm->GetId();
2908 if ( faceTag != BLSURFPlugin_Hypothesis::GetHyperPatchTag( faceTag, _hypothesis ))
2910 // triangles are assigned to the first face of hyper-patch
2911 sm->SetIsAlwaysComputed( true );
2915 sm->GetComputeError().reset( new SMESH_ComputeError( err, _comment, this ));
2916 badFaceFound = true;
2920 if ( err == COMPERR_WARNING )
2924 if ( status != STATUS_OK && !badFaceFound ) {
2928 // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
2930 if ( oldFEFlags > 0 )
2931 feenableexcept( oldFEFlags );
2932 feclearexcept( FE_ALL_EXCEPT );
2936 std::cout << "FacesWithSizeMap" << std::endl;
2937 FacesWithSizeMap.Statistics(std::cout);
2938 std::cout << "EdgesWithSizeMap" << std::endl;
2939 EdgesWithSizeMap.Statistics(std::cout);
2940 std::cout << "VerticesWithSizeMap" << std::endl;
2941 VerticesWithSizeMap.Statistics(std::cout);
2942 std::cout << "FacesWithEnforcedVertices" << std::endl;
2943 FacesWithEnforcedVertices.Statistics(std::cout);
2946 return ( status == STATUS_OK && !quadraticSubMeshAndViscousLayer );
2949 //================================================================================
2951 * \brief Compute a mesh basing on discrete CAD description
2953 //================================================================================
2955 bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
2957 if ( aMesh.NbFaces() == 0 )
2958 return error( COMPERR_BAD_INPUT_MESH, "2D elements are missing" );
2960 context_t *ctx = context_new();
2961 if (!ctx) return error("Pb in context_new()");
2963 BLSURF_Cleaner cleaner( ctx );
2965 message_cb_user_data mcud;
2966 mcud._error = & this->SMESH_Algo::_comment;
2967 mcud._progress = & this->SMESH_Algo::_progress;
2969 _hypothesis ? _hypothesis->GetVerbosity() : BLSURFPlugin_Hypothesis::GetDefaultVerbosity();
2970 meshgems_status_t ret = context_set_message_callback(ctx, message_cb, &mcud);
2971 if (ret != STATUS_OK) return error("Pb. in context_set_message_callback() ");
2973 cadsurf_session_t * css = cadsurf_session_new(ctx);
2974 if(!css) return error( "Pb. in cadsurf_session_new() " );
2978 // Fill an input mesh
2980 mesh_t * msh = meshgems_mesh_new_in_memory( ctx );
2981 if ( !msh ) return error("Pb. in meshgems_mesh_new_in_memory()");
2983 // mark nodes used by 2D elements
2984 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
2985 SMDS_NodeIteratorPtr nodeIt = meshDS->nodesIterator();
2986 while ( nodeIt->more() )
2988 const SMDS_MeshNode* n = nodeIt->next();
2989 n->setIsMarked( n->NbInverseElements( SMDSAbs_Face ));
2991 meshgems_mesh_set_vertex_count( msh, meshDS->NbNodes() );
2993 // set node coordinates
2994 if ( meshDS->NbNodes() != meshDS->MaxNodeID() )
2996 meshDS->CompactMesh();
2999 nodeIt = meshDS->nodesIterator();
3001 for ( i = 1; nodeIt->more(); ++i )
3003 nXYZ.Set( nodeIt->next() );
3004 meshgems_mesh_set_vertex_coordinates( msh, i, nXYZ.ChangeData() );
3007 // set nodes of faces
3008 meshgems_mesh_set_triangle_count ( msh, meshDS->GetMeshInfo().NbTriangles() );
3009 meshgems_mesh_set_quadrangle_count( msh, meshDS->GetMeshInfo().NbQuadrangles() );
3010 meshgems_integer nodeIDs[4];
3011 meshgems_integer iT = 1, iQ = 1;
3012 SMDS_FaceIteratorPtr faceIt = meshDS->facesIterator();
3013 while ( faceIt->more() )
3015 const SMDS_MeshElement* face = faceIt->next();
3016 meshgems_integer nbNodes = face->NbCornerNodes();
3017 if ( nbNodes > 4 || face->IsPoly() ) continue;
3019 for ( i = 0; i < nbNodes; ++i )
3020 nodeIDs[i] = face->GetNode( i )->GetID();
3022 meshgems_mesh_set_triangle_vertices ( msh, iT++, nodeIDs );
3024 meshgems_mesh_set_quadrangle_vertices( msh, iQ++, nodeIDs );
3027 ret = cadsurf_set_mesh(css, msh);
3028 if ( ret != STATUS_OK ) return error("Pb in cadsurf_set_mesh()");
3033 SetParameters(_hypothesis, css, aMesh.GetShapeToMesh() );
3035 ret = cadsurf_compute_mesh(css);
3036 if ( ret != STATUS_OK ) return false;
3039 cadsurf_get_mesh(css, &omsh);
3040 if ( !omsh ) return error( "Pb. in cadsurf_get_mesh()" );
3043 // Update SALOME mesh
3045 // remove quadrangles and triangles
3046 for ( faceIt = meshDS->facesIterator(); faceIt->more(); )
3048 const SMDS_MeshElement* face = faceIt->next();
3049 if ( !face->IsPoly() )
3050 meshDS->RemoveFreeElement( face, /*sm=*/0, /*fromGroups=*/true );
3052 // remove edges that bound the just removed faces
3053 for ( SMDS_EdgeIteratorPtr edgeIt = meshDS->edgesIterator(); edgeIt->more(); )
3055 const SMDS_MeshElement* edge = edgeIt->next();
3056 const SMDS_MeshNode* n0 = edge->GetNode(0);
3057 const SMDS_MeshNode* n1 = edge->GetNode(1);
3058 if ( n0->isMarked() &&
3060 n0->NbInverseElements( SMDSAbs_Volume ) == 0 &&
3061 n1->NbInverseElements( SMDSAbs_Volume ) == 0 )
3062 meshDS->RemoveFreeElement( edge, /*sm=*/0, /*fromGroups=*/true );
3064 // remove nodes that just became free
3065 for ( nodeIt = meshDS->nodesIterator(); nodeIt->more(); )
3067 const SMDS_MeshNode* n = nodeIt->next();
3068 if ( n->isMarked() && n->NbInverseElements() == 0 )
3069 meshDS->RemoveFreeNode( n, /*sm=*/0, /*fromGroups=*/true );
3073 meshgems_integer nbvtx = 0, nodeID;
3074 meshgems_mesh_get_vertex_count( omsh, &nbvtx );
3075 meshgems_real xyz[3];
3076 for ( i = 1; i <= nbvtx; ++i )
3078 meshgems_mesh_get_vertex_coordinates( omsh, i, xyz );
3079 SMDS_MeshNode* n = meshDS->AddNode( xyz[0], xyz[1], xyz[2] );
3080 nodeID = n->GetID();
3081 meshgems_mesh_set_vertex_tag( omsh, i, &nodeID ); // save mapping of IDs in MG and SALOME meshes
3085 meshgems_integer nbtri = 0;
3086 meshgems_mesh_get_triangle_count( omsh, &nbtri );
3087 const SMDS_MeshNode* nodes[4];
3088 for ( i = 1; i <= nbtri; ++i )
3090 meshgems_mesh_get_triangle_vertices( omsh, i, nodeIDs );
3091 for ( int j = 0; j < 3; ++j )
3093 meshgems_mesh_get_vertex_tag( omsh, nodeIDs[j], &nodeID );
3094 nodes[j] = meshDS->FindNode( nodeID );
3096 meshDS->AddFace( nodes[0], nodes[1], nodes[2] );
3100 meshgems_integer nbquad = 0;
3101 meshgems_mesh_get_quadrangle_count( omsh, &nbquad );
3102 for ( i = 1; i <= nbquad; ++i )
3104 meshgems_mesh_get_quadrangle_vertices( omsh, i, nodeIDs );
3105 for ( int j = 0; j < 4; ++j )
3107 meshgems_mesh_get_vertex_tag( omsh, nodeIDs[j], &nodeID );
3108 nodes[j] = meshDS->FindNode( nodeID );
3110 meshDS->AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
3115 std::string GMFFileName = _hypothesis->GetGMFFile();
3116 if ( !GMFFileName.empty() )
3118 bool asciiFound = (GMFFileName.find(".mesh", GMFFileName.size()-5) != std::string::npos);
3119 bool binaryFound = (GMFFileName.find(".meshb",GMFFileName.size()-6) != std::string::npos);
3120 if ( !asciiFound && !binaryFound )
3121 GMFFileName.append(".mesh");
3122 mesh_write_mesh(msh, GMFFileName.c_str());
3126 cadsurf_regain_mesh(css, omsh);
3128 // as we don't assign the new triangles to a shape (the pseudo-shape),
3129 // we mark the shape as always computed to avoid the error messages
3130 // that no elements assigned to the shape
3131 aMesh.GetSubMesh( aHelper->GetSubShape() )->SetIsAlwaysComputed( true );
3136 //================================================================================
3138 * \brief Terminates computation
3140 //================================================================================
3142 void BLSURFPlugin_BLSURF::CancelCompute()
3144 _compute_canceled = true;
3147 //=============================================================================
3151 //=============================================================================
3153 void BLSURFPlugin_BLSURF::Set_NodeOnEdge(SMESHDS_Mesh* meshDS,
3154 const SMDS_MeshNode* node,
3155 const TopoDS_Shape& ed)
3157 const TopoDS_Edge edge = TopoDS::Edge(ed);
3159 gp_Pnt pnt(node->X(), node->Y(), node->Z());
3161 Standard_Real p0 = 0.0;
3162 Standard_Real p1 = 1.0;
3163 TopLoc_Location loc;
3164 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, loc, p0, p1);
3165 if ( curve.IsNull() )
3167 // issue 22499. Node at a sphere apex
3168 meshDS->SetNodeOnEdge(node, edge, p0);
3172 if ( !loc.IsIdentity() ) pnt.Transform( loc.Transformation().Inverted() );
3173 GeomAPI_ProjectPointOnCurve proj(pnt, curve, p0, p1);
3176 if ( proj.NbPoints() > 0 )
3178 pa = (double)proj.LowerDistanceParameter();
3179 // Issue 0020656. Move node if it is too far from edge
3180 gp_Pnt curve_pnt = curve->Value( pa );
3181 double dist2 = pnt.SquareDistance( curve_pnt );
3182 double tol = BRep_Tool::Tolerance( edge );
3183 if ( 1e-14 < dist2 && dist2 <= 1000*tol ) // large enough and within tolerance
3185 curve_pnt.Transform( loc );
3186 meshDS->MoveNode( node, curve_pnt.X(), curve_pnt.Y(), curve_pnt.Z() );
3190 meshDS->SetNodeOnEdge(node, edge, pa);
3193 /* Curve definition function See cad_curv_t in file meshgems/cad.h for
3195 * NOTE : if when your CAD systems evaluates second
3196 * order derivatives it also computes first order derivatives and
3197 * function evaluation, you can optimize this example by making only
3198 * one CAD call and filling the necessary uv, dt, dtt arrays.
3200 status_t curv_fun(real t, real *uv, real *dt, real *dtt, void *user_data)
3202 /* t is given. It contains the t (time) 1D parametric coordintaes
3203 of the point PreCAD/MG-CADSurf is querying on the curve */
3205 /* user_data identifies the edge PreCAD/MG-CADSurf is querying
3206 * (see cad_edge_new later in this example) */
3207 const Geom2d_Curve*pargeo = (const Geom2d_Curve*) user_data;
3210 /* MG-CADSurf is querying the function evaluation */
3213 uv[0]=P.X(); uv[1]=P.Y();
3217 /* query for the first order derivatives */
3220 dt[0]=V1.X(); dt[1]=V1.Y();
3224 /* query for the second order derivatives */
3227 dtt[0]=V2.X(); dtt[1]=V2.Y();
3233 /* Surface definition function.
3234 * See cad_surf_t in file meshgems/cad.h for more information.
3235 * NOTE : if when your CAD systems evaluates second order derivatives it also
3236 * computes first order derivatives and function evaluation, you can optimize
3237 * this example by making only one CAD call and filling the necessary xyz, du, dv, etc..
3240 status_t surf_fun(real *uv, real *xyz, real*du, real *dv,
3241 real *duu, real *duv, real *dvv, void *user_data)
3243 /* uv[2] is given. It contains the u,v coordinates of the point
3244 * PreCAD/MG-CADSurf is querying on the surface */
3246 /* user_data identifies the face PreCAD/MG-CADSurf is querying (see
3247 * cad_face_new later in this example)*/
3248 const Geom_Surface* geometry = (const Geom_Surface*) user_data;
3252 P=geometry->Value(uv[0],uv[1]); // S.D0(U,V,P);
3253 xyz[0]=P.X(); xyz[1]=P.Y(); xyz[2]=P.Z();
3260 geometry->D1(uv[0],uv[1],P,D1U,D1V);
3261 du[0]=D1U.X(); du[1]=D1U.Y(); du[2]=D1U.Z();
3262 dv[0]=D1V.X(); dv[1]=D1V.Y(); dv[2]=D1V.Z();
3265 if(duu && duv && dvv){
3269 gp_Vec D2U,D2V,D2UV;
3271 geometry->D2(uv[0],uv[1],P,D1U,D1V,D2U,D2V,D2UV);
3272 duu[0]=D2U.X(); duu[1]=D2U.Y(); duu[2]=D2U.Z();
3273 duv[0]=D2UV.X(); duv[1]=D2UV.Y(); duv[2]=D2UV.Z();
3274 dvv[0]=D2V.X(); dvv[1]=D2V.Y(); dvv[2]=D2V.Z();
3281 status_t size_on_surface(integer face_id, real *uv, real *size, void *user_data)
3283 TId2ClsAttractorVec::iterator f2attVec;
3284 if (FaceId2PythonSmp.count(face_id) != 0) {
3285 assert(Py_IsInitialized());
3286 PyGILState_STATE gstate;
3287 gstate = PyGILState_Ensure();
3288 PyObject* pyresult = PyObject_CallFunction(FaceId2PythonSmp[face_id],(char*)"(f,f)",uv[0],uv[1]);
3290 if ( pyresult != NULL) {
3291 result = PyFloat_AsDouble(pyresult);
3292 Py_DECREF(pyresult);
3297 string err_description="";
3298 PyObject* new_stderr = newPyStdOut(err_description);
3299 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3300 Py_INCREF(old_stderr);
3301 PySys_SetObject((char*)"stderr", new_stderr);
3303 PySys_SetObject((char*)"stderr", old_stderr);
3304 Py_DECREF(new_stderr);
3305 MESSAGE("Can't evaluate f(" << uv[0] << "," << uv[1] << ")" << " error is " << err_description);
3306 result = *((real*)user_data);
3309 PyGILState_Release(gstate);
3311 else if (( f2attVec = FaceIndex2ClassAttractor.find(face_id)) != FaceIndex2ClassAttractor.end() && !f2attVec->second.empty())
3315 std::vector< BLSURFPlugin_Attractor* > & attVec = f2attVec->second;
3316 for ( size_t i = 0; i < attVec.size(); ++i )
3318 //result += attVec[i]->GetSize(uv[0],uv[1]);
3319 result = Min( result, attVec[i]->GetSize(uv[0],uv[1]));
3321 //*size = result / attVec.size(); // mean of sizes defined by all attractors
3325 *size = *((real*)user_data);
3327 // std::cout << "Size_on_surface sur la face " << face_id << " donne une size de: " << *size << std::endl;
3331 status_t size_on_edge(integer edge_id, real t, real *size, void *user_data)
3333 if (EdgeId2PythonSmp.count(edge_id) != 0){
3334 assert(Py_IsInitialized());
3335 PyGILState_STATE gstate;
3336 gstate = PyGILState_Ensure();
3337 PyObject* pyresult = PyObject_CallFunction(EdgeId2PythonSmp[edge_id],(char*)"(f)",t);
3339 if ( pyresult != NULL) {
3340 result = PyFloat_AsDouble(pyresult);
3341 Py_DECREF(pyresult);
3346 string err_description="";
3347 PyObject* new_stderr = newPyStdOut(err_description);
3348 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3349 Py_INCREF(old_stderr);
3350 PySys_SetObject((char*)"stderr", new_stderr);
3352 PySys_SetObject((char*)"stderr", old_stderr);
3353 Py_DECREF(new_stderr);
3354 MESSAGE("Can't evaluate f(" << t << ")" << " error is " << err_description);
3355 result = *((real*)user_data);
3358 PyGILState_Release(gstate);
3361 *size = *((real*)user_data);
3366 status_t size_on_vertex(integer point_id, real *size, void *user_data)
3368 if (VertexId2PythonSmp.count(point_id) != 0){
3369 assert(Py_IsInitialized());
3370 PyGILState_STATE gstate;
3371 gstate = PyGILState_Ensure();
3372 PyObject* pyresult = PyObject_CallFunction(VertexId2PythonSmp[point_id],(char*)"");
3374 if ( pyresult != NULL) {
3375 result = PyFloat_AsDouble(pyresult);
3376 Py_DECREF(pyresult);
3381 string err_description="";
3382 PyObject* new_stderr = newPyStdOut(err_description);
3383 PyObject* old_stderr = PySys_GetObject((char*)"stderr");
3384 Py_INCREF(old_stderr);
3385 PySys_SetObject((char*)"stderr", new_stderr);
3387 PySys_SetObject((char*)"stderr", old_stderr);
3388 Py_DECREF(new_stderr);
3389 MESSAGE("Can't evaluate f()" << " error is " << err_description);
3390 result = *((real*)user_data);
3393 PyGILState_Release(gstate);
3396 *size = *((real*)user_data);
3402 * The following function will be called for PreCAD/MG-CADSurf message
3403 * printing. See context_set_message_callback (later in this
3404 * template) for how to set user_data.
3406 status_t message_cb(message_t *msg, void *user_data)
3408 integer errnumber = 0;
3410 message_get_number(msg, &errnumber);
3411 message_get_description(msg, &desc);
3413 message_cb_user_data * mcud = (message_cb_user_data*)user_data;
3414 // Get all the error message and some warning messages related to license and periodicity
3415 if ( errnumber < 0 ||
3416 err.find("license" ) != string::npos ||
3417 err.find("periodicity") != string::npos )
3419 // remove ^A from the tail
3420 int len = strlen( desc );
3421 while (len > 0 && desc[len-1] != '\n')
3423 mcud->_error->append( desc, len );
3426 if ( errnumber == 3009001 )
3427 * mcud->_progress = atof( desc + 11 ) / 100.;
3428 if ( mcud->_verbosity > 0 )
3429 std::cout << desc << std::endl;
3434 /* This is the interrupt callback. PreCAD/MG-CADSurf will call this
3435 * function regularily. See the file meshgems/interrupt.h
3437 status_t interrupt_cb(integer *interrupt_status, void *user_data)
3439 integer you_want_to_continue = 1;
3440 BLSURFPlugin_BLSURF* tmp = (BLSURFPlugin_BLSURF*)user_data;
3441 you_want_to_continue = !tmp->computeCanceled();
3443 if(you_want_to_continue)
3445 *interrupt_status = INTERRUPT_CONTINUE;
3448 else /* you want to stop MG-CADSurf */
3450 *interrupt_status = INTERRUPT_STOP;
3451 return STATUS_ERROR;
3455 //=============================================================================
3459 //=============================================================================
3460 bool BLSURFPlugin_BLSURF::Evaluate(SMESH_Mesh& aMesh,
3461 const TopoDS_Shape& aShape,
3462 MapShapeNbElems& aResMap)
3464 double diagonal = aMesh.GetShapeDiagonalSize();
3465 double bbSegmentation = _gen->GetBoundaryBoxSegmentation();
3466 int _physicalMesh = BLSURFPlugin_Hypothesis::GetDefaultPhysicalMesh();
3467 double _phySize = BLSURFPlugin_Hypothesis::GetDefaultPhySize(diagonal, bbSegmentation);
3468 bool _phySizeRel = BLSURFPlugin_Hypothesis::GetDefaultPhySizeRel();
3469 //int _geometricMesh = BLSURFPlugin_Hypothesis::GetDefaultGeometricMesh();
3470 double _angleMesh = BLSURFPlugin_Hypothesis::GetDefaultAngleMesh();
3471 BLSURFPlugin_Hypothesis::ElementType _elementType = BLSURFPlugin_Hypothesis::GetDefaultElementType();
3473 _physicalMesh = (int) _hypothesis->GetPhysicalMesh();
3474 _phySizeRel = _hypothesis->IsPhySizeRel();
3475 if ( _hypothesis->GetPhySize() > 0)
3476 _phySize = _phySizeRel ? diagonal*_hypothesis->GetPhySize() : _hypothesis->GetPhySize();
3477 //_geometricMesh = (int) hyp->GetGeometricMesh();
3478 if (_hypothesis->GetAngleMesh() > 0)
3479 _angleMesh = _hypothesis->GetAngleMesh();
3480 _elementType = _hypothesis->GetElementType();
3482 //0020968: EDF1545 SMESH: Problem in the creation of a mesh group on geometry
3483 // GetDefaultPhySize() sometimes leads to computation failure
3484 _phySize = aMesh.GetShapeDiagonalSize() / _gen->GetBoundaryBoxSegmentation();
3487 bool IsQuadratic = _quadraticMesh;
3492 TopTools_DataMapOfShapeInteger EdgesMap;
3493 double fullLen = 0.0;
3494 double fullNbSeg = 0;
3495 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
3496 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3497 if( EdgesMap.IsBound(E) )
3499 SMESH_subMesh *sm = aMesh.GetSubMesh(E);
3500 double aLen = SMESH_Algo::EdgeLength(E);
3503 if(_physicalMesh==1) {
3504 nb1d = (int)( aLen/_phySize + 1 );
3509 Handle(Geom_Curve) C = BRep_Tool::Curve(E,f,l);
3510 double fullAng = 0.0;
3511 double dp = (l-f)/200;
3516 for(int j=2; j<=200; j++) {
3519 fullAng += fabs(V1.Angle(V2));
3523 nb1d = (int)( fullAng/_angleMesh + 1 );
3526 std::vector<int> aVec(SMDSEntity_Last);
3527 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
3528 if( IsQuadratic > 0 ) {
3529 aVec[SMDSEntity_Node] = 2*nb1d - 1;
3530 aVec[SMDSEntity_Quad_Edge] = nb1d;
3533 aVec[SMDSEntity_Node] = nb1d - 1;
3534 aVec[SMDSEntity_Edge] = nb1d;
3536 aResMap.insert(std::make_pair(sm,aVec));
3537 EdgesMap.Bind(E,nb1d);
3539 double ELen = fullLen/fullNbSeg;
3543 // try to evaluate as in MEFISTO
3544 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
3545 TopoDS_Face F = TopoDS::Face( exp.Current() );
3546 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
3548 BRepGProp::SurfaceProperties(F,G);
3549 double anArea = G.Mass();
3551 std::vector<int> nb1dVec;
3552 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next()) {
3553 int nbSeg = EdgesMap.Find(exp1.Current());
3555 nb1dVec.push_back( nbSeg );
3558 int nbTria = (int) ( anArea/( ELen*ELen*sqrt(3.) / 4 ) );
3559 int nbNodes = (int) ( ( nbTria*3 - (nb1d-1)*2 ) / 6 + 1 );
3560 if ( _elementType != BLSURFPlugin_Hypothesis::Quadrangles )
3562 if ( nb1dVec.size() == 4 ) // quadrangle geom face
3564 int n1 = nb1dVec[0], n2 = nb1dVec[ nb1dVec[1] == nb1dVec[0] ? 2 : 1 ];
3566 nbNodes = (n1 + 1) * (n2 + 1);
3571 nbTria = nbQuad = nbTria / 3 + 1;
3574 std::vector<int> aVec(SMDSEntity_Last,0);
3576 int nb1d_in = (nbTria*3 - nb1d) / 2;
3577 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3578 aVec[SMDSEntity_Quad_Triangle] = nbTria;
3579 aVec[SMDSEntity_Quad_Quadrangle] = nbQuad;
3582 aVec[SMDSEntity_Node] = nbNodes;
3583 aVec[SMDSEntity_Triangle] = nbTria;
3584 aVec[SMDSEntity_Quadrangle] = nbQuad;
3586 aResMap.insert(std::make_pair(sm,aVec));
3593 BRepGProp::VolumeProperties(aShape,G);
3594 double aVolume = G.Mass();
3595 double tetrVol = 0.1179*ELen*ELen*ELen;
3596 int nbVols = int(aVolume/tetrVol);
3597 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3598 std::vector<int> aVec(SMDSEntity_Last);
3599 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
3601 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3602 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3605 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3606 aVec[SMDSEntity_Tetra] = nbVols;
3608 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
3609 aResMap.insert(std::make_pair(sm,aVec));
3614 //================================================================================
3616 * \brief Find TopoDS_Shape for each hyper-patch study entry in a hypothesis
3618 //================================================================================
3620 void BLSURFPlugin_BLSURF::FillEntryToShape( const BLSURFPlugin_Hypothesis* hyp,
3621 std::map< std::string, TopoDS_Shape > & entryToShape )
3623 SMESH_Gen_i* smeshGen = SMESH_Gen_i::GetSMESHGen();
3624 for ( const ::BLSURFPlugin_Hypothesis::THyperPatchEntries& entries : hyp->GetHyperPatchEntries() )
3625 for ( const std::string& entry : entries )
3627 GEOM::GEOM_Object_var go = smeshGen->GetGeomObjectByEntry( entry );
3628 TopoDS_Shape shape = smeshGen->GeomObjectToShape( go );
3629 if ( !shape.IsNull() )
3630 entryToShape.insert({ entry, shape });