Salome HOME
[bos #40653][CEA] New mesh import export formats with meshio.
[modules/smesh.git] / src / StdMeshers / StdMeshers_ViscousLayers2D.cxx
index 2de6e474c556a93f01fa94c5d496232c637a5a60..b89a42c2d2cb385ce1264b424a71d5c2b1cec7fe 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2016  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2024  CEA, EDF, OPEN CASCADE
 //
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
@@ -38,6 +38,7 @@
 #include "SMESH_Group.hxx"
 #include "SMESH_HypoFilter.hxx"
 #include "SMESH_Mesh.hxx"
+#include "SMESH_MeshEditor.hxx"
 #include "SMESH_MesherHelper.hxx"
 #include "SMESH_ProxyMesh.hxx"
 #include "SMESH_Quadtree.hxx"
@@ -82,6 +83,7 @@
 #include <gp_Ax1.hxx>
 #include <gp_Vec.hxx>
 #include <gp_XY.hxx>
+#include <smIdType.hxx>
 
 #include <list>
 #include <string>
@@ -110,15 +112,15 @@ namespace VISCOUS_2D
     // Proxy sub-mesh of an EDGE. It contains nodes in _uvPtStructVec.
     struct _EdgeSubMesh : public SMESH_ProxyMesh::SubMesh
     {
-      _EdgeSubMesh(int index=0): SubMesh(index) {}
+      _EdgeSubMesh(const SMDS_Mesh* mesh, int index=0): SubMesh(mesh,index) {}
       //virtual int NbElements() const { return _elements.size()+1; }
-      virtual int NbNodes() const { return Max( 0, _uvPtStructVec.size()-2 ); }
+      virtual smIdType NbNodes() const { return Max( 0, _uvPtStructVec.size()-2 ); }
       void SetUVPtStructVec(UVPtStructVec& vec) { _uvPtStructVec.swap( vec ); }
       UVPtStructVec& GetUVPtStructVec() { return _uvPtStructVec; }
     };
     _ProxyMeshOfFace(const SMESH_Mesh& mesh): SMESH_ProxyMesh(mesh) {}
     _EdgeSubMesh* GetEdgeSubMesh(int ID) { return (_EdgeSubMesh*) getProxySubMesh(ID); }
-    virtual SubMesh* newSubmesh(int index=0) const { return new _EdgeSubMesh(index); }
+    virtual SubMesh* newSubmesh(int index=0) const { return new _EdgeSubMesh( GetMeshDS(), index); }
   };
   //--------------------------------------------------------------------------------
   /*!
@@ -149,7 +151,7 @@ namespace VISCOUS_2D
     // Treat events
     void ProcessEvent(const int          event,
                       const int          eventType,
-                      SMESH_subMesh*     subMesh,
+                      SMESH_subMesh*     /*subMesh*/,
                       EventListenerData* data,
                       const SMESH_Hypothesis*  /*hyp*/)
     {
@@ -176,7 +178,7 @@ namespace VISCOUS_2D
    */
   struct _Segment
   {
-    const gp_XY* _uv[2];       // poiter to _LayerEdge::_uvIn
+    const gp_XY* _uv[2];       // pointer to _LayerEdge::_uvIn
     int          _indexInLine; // position in _PolyLine
 
     _Segment() {}
@@ -238,9 +240,7 @@ namespace VISCOUS_2D
 
     bool SetNewLength( const double length );
 
-#ifdef _DEBUG_
-    int           _ID;
-#endif
+    int           _ID; // debug
   };
   //--------------------------------------------------------------------------------
   /*!
@@ -501,15 +501,15 @@ namespace VISCOUS_2D
 //================================================================================
 // StdMeshers_ViscousLayers hypothesis
 //
-StdMeshers_ViscousLayers2D::StdMeshers_ViscousLayers2D(int hypId, int studyId, SMESH_Gen* gen)
-  :StdMeshers_ViscousLayers(hypId, studyId, gen)
+StdMeshers_ViscousLayers2D::StdMeshers_ViscousLayers2D(int hypId, SMESH_Gen* gen)
+  :StdMeshers_ViscousLayers(hypId, gen)
 {
   _name = StdMeshers_ViscousLayers2D::GetHypType();
   _param_algo_dim = -2; // auxiliary hyp used by 2D algos
 }
 // --------------------------------------------------------------------------------
-bool StdMeshers_ViscousLayers2D::SetParametersByMesh(const SMESH_Mesh*   theMesh,
-                                                     const TopoDS_Shape& theShape)
+bool StdMeshers_ViscousLayers2D::SetParametersByMesh(const SMESH_Mesh*   /*theMesh*/,
+                                                     const TopoDS_Shape& /*theShape*/)
 {
   // TODO ???
   return false;
@@ -585,7 +585,7 @@ StdMeshers_ViscousLayers2D::CheckHypothesis(SMESH_Mesh&
       VISCOUS_2D::_ViscousBuilder2D builder( theMesh, face, hyps, hypShapes );
       builder._faceSideVec =
         StdMeshers_FaceSide::GetFaceWires( face, theMesh, true, error,
-                                           SMESH_ProxyMesh::Ptr(),
+                                           NULL, SMESH_ProxyMesh::Ptr(),
                                            /*theCheckVertexNodes=*/false);
       if ( error->IsOK() && !builder.findEdgesWithLayers())
       {
@@ -600,7 +600,7 @@ StdMeshers_ViscousLayers2D::CheckHypothesis(SMESH_Mesh&
 // --------------------------------------------------------------------------------
 void StdMeshers_ViscousLayers2D::RestoreListeners() const
 {
-  StudyContextStruct* sc = _gen->GetStudyContext( _studyId );
+  StudyContextStruct* sc = _gen->GetStudyContext();
   std::map < int, SMESH_Mesh * >::iterator i_smesh = sc->mapMesh.begin();
   for ( ; i_smesh != sc->mapMesh.end(); ++i_smesh )
   {
@@ -672,9 +672,10 @@ bool _ViscousBuilder2D::error(const string& text )
       _error->myAlgo = smError->myAlgo;
     smError = _error;
   }
-#ifdef _DEBUG_
-  cout << "_ViscousBuilder2D::error " << text << endl;
-#endif
+
+  if (SALOME::VerbosityActivated())
+    cout << "_ViscousBuilder2D::error " << text << endl;
+
   return false;
 }
 
@@ -686,7 +687,7 @@ bool _ViscousBuilder2D::error(const string& text )
 
 SMESH_ProxyMesh::Ptr _ViscousBuilder2D::Compute()
 {
-  _faceSideVec = StdMeshers_FaceSide::GetFaceWires( _face, *_mesh, true, _error);
+  _faceSideVec = StdMeshers_FaceSide::GetFaceWires( _face, *_mesh, true, _error, &_helper );
 
   if ( !_error->IsOK() )
     return _proxyMesh;
@@ -1263,7 +1264,7 @@ void _ViscousBuilder2D::adjustCommonEdge( _PolyLine& LL, _PolyLine& LR )
           if ( isR )
             LR._lEdges.erase( LR._lEdges.begin()+1, eIt );
           else
-            LL._lEdges.erase( eIt, --LL._lEdges.end() );
+            LL._lEdges.erase( eIt+1, --LL._lEdges.end() );
           // eIt = isR ? L._lEdges.begin()+1 : L._lEdges.end()-2;
           // for ( size_t i = 1; i < iLE; ++i, eIt += dIt )
           //   eIt->_isBlocked = true;
@@ -1318,7 +1319,7 @@ void _ViscousBuilder2D::setLayerEdgeData( _LayerEdge&                 lEdge,
     faceProj->Perform( p );
     if ( !faceProj->IsDone() || faceProj->NbPoints() < 1 )
       return setLayerEdgeData( lEdge, u, pcurve, curve, p, reverse, NULL );
-    Quantity_Parameter U,V;
+    Standard_Real U,V;
     faceProj->LowerDistanceParameters(U,V);
     lEdge._normal2D.SetCoord( U - uv.X(), V - uv.Y() );
     lEdge._normal2D.Normalize();
@@ -1337,9 +1338,9 @@ void _ViscousBuilder2D::setLayerEdgeData( _LayerEdge&                 lEdge,
   lEdge._ray.SetDirection( lEdge._normal2D );
   lEdge._isBlocked = false;
   lEdge._length2D  = 0;
-#ifdef _DEBUG_
-  lEdge._ID        = _nbLE++;
-#endif
+
+  if (SALOME::VerbosityActivated())
+    lEdge._ID        = _nbLE++;
 }
 
 //================================================================================
@@ -1998,7 +1999,7 @@ bool _ViscousBuilder2D::shrink()
         throw SALOME_Exception(SMESH_Comment("ViscousBuilder2D: not SMDS_TOP_EDGE node position: ")
                                << oldNode->GetPosition()->GetTypeOfPosition()
                                << " of node " << oldNode->GetID());
-      SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( oldNode->GetPosition() );
+      SMDS_EdgePositionPtr pos = oldNode->GetPosition();
       pos->SetUParameter( nodeDataVec[iP].param );
 
       gp_Pnt newP = curve.Value( nodeDataVec[iP].param );
@@ -2187,7 +2188,7 @@ bool _ViscousBuilder2D::refine()
     vector< double > segLen( L._lEdges.size() );
     segLen[0] = 0.0;
 
-    // check if length modification is usefull: look for _LayerEdge's
+    // check if length modification is useful: look for _LayerEdge's
     // with length limited due to collisions
     bool lenLimited = false;
     for ( size_t iLE = 1; ( iLE < L._lEdges.size()-1 && !lenLimited ); ++iLE )
@@ -2401,6 +2402,17 @@ bool _ViscousBuilder2D::refine()
       outerNodes.swap( innerNodes );
     }
 
+    // Add faces to a group
+    SMDS_MeshGroup* group = StdMeshers_ViscousLayers::CreateGroup( hyp->GetGroupName(),
+                                                                   *_helper.GetMesh(),
+                                                                   SMDSAbs_Face );
+    if ( group )
+    {
+      TIDSortedElemSet::iterator fIt = L._newFaces.begin();
+      for ( ; fIt != L._newFaces.end(); ++fIt )
+        group->Add( *fIt );
+    }
+
     // faces between not shared _LayerEdge's (at concave VERTEX)
     for ( int isR = 0; isR < 2; ++isR )
     {
@@ -2412,15 +2424,22 @@ bool _ViscousBuilder2D::refine()
       if ( lNodes.empty() || rNodes.empty() || lNodes.size() != rNodes.size() )
         continue;
 
+      const SMDS_MeshElement* face = 0;
       for ( size_t i = 1; i < lNodes.size(); ++i )
-        _helper.AddFace( lNodes[ i+prev ], rNodes[ i+prev ],
-                         rNodes[ i+cur ],  lNodes[ i+cur ]);
+      {
+        face = _helper.AddFace( lNodes[ i+prev ], rNodes[ i+prev ],
+                                rNodes[ i+cur ],  lNodes[ i+cur ]);
+        if ( group )
+          group->Add( face );
+      }
 
       const UVPtStruct& ptOnVertex = points[ isR ? L._lastPntInd : L._firstPntInd ];
       if ( isReverse )
-        _helper.AddFace( ptOnVertex.node, lNodes[ 0 ], rNodes[ 0 ]);
+        face = _helper.AddFace( ptOnVertex.node, lNodes[ 0 ], rNodes[ 0 ]);
       else
-        _helper.AddFace( ptOnVertex.node, rNodes[ 0 ], lNodes[ 0 ]);
+        face = _helper.AddFace( ptOnVertex.node, rNodes[ 0 ], lNodes[ 0 ]);
+      if ( group )
+        group->Add( face );
     }
 
     // Fill the _ProxyMeshOfFace