]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
[bos #42217][EDF 28921] Horseshoe with bodyfitting.
authorKonstantin Leontev <Konstantin.LEONTEV@opencascade.com>
Mon, 5 Aug 2024 12:44:07 +0000 (13:44 +0100)
committerKonstantin Leontev <Konstantin.LEONTEV@opencascade.com>
Mon, 5 Aug 2024 12:44:07 +0000 (13:44 +0100)
Nodes linked to a null nodes without any possible splits now being cleared befor splits initialization. This prevents building non-planar faces in cases when the source geometry surface coincident with an intermediate grid plane.

src/StdMeshers/StdMeshers_Cartesian_3D_Hexahedron.cxx
src/StdMeshers/StdMeshers_Cartesian_3D_Hexahedron.hxx

index 677dd676101b85151d96838c37cd4e433b98e30f..66409b01a3c12988aec951eb7c91cc9b3fd55c95 100644 (file)
@@ -393,6 +393,28 @@ void Hexahedron::_Node::Add( const E_IntersectPoint* ip )
     _node = _intPoint->_node = node;
 }
 
+void Hexahedron::_Node::clear()
+{
+  _node = nullptr;
+  _boundaryCornerNode = nullptr;
+  _intPoint = nullptr;
+  _usedInFace = nullptr;
+  _isInternalFlags = IS_NOT_INTERNAL;
+}
+
+void Hexahedron::_Link::clear()
+{
+  for (std::size_t i = 0; i < nodesNum; ++i)
+  {
+    if (_nodes[i])
+      _nodes[i]->clear();
+  }
+
+  _fIntPoints.clear();
+  _fIntNodes.clear();
+  _splits.clear();
+}
+
 //================================================================================
 /*!
   * \brief Return IDs of SOLIDs interfering with this Hexahedron
@@ -606,6 +628,8 @@ void Hexahedron::init( size_t i, size_t j, size_t k, const Solid* solid )
     // this method can be called in parallel, so use own helper
     SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
 
+    clearNodesLinkedToNull(solid, helper);
+
     // Create sub-links (_Link::_splits) by splitting links with _Link::_fIntPoints
     // ---------------------------------------------------------------
     _Link split;
@@ -798,6 +822,104 @@ void Hexahedron::init( size_t i, size_t j, size_t k, const Solid* solid )
 
 } // init( _i, _j, _k )
 
+//================================================================================
+/*!
+  * \brief Clear nodes linked to null by not splitted links.
+  * Exclude from computation all nodes those linked to null without any splits in links.
+  * This prevents building non-planar faces in cases when the source geometry surface
+  * coincident with an intermediate grid plane.
+  */
+void Hexahedron::clearNodesLinkedToNull(const Solid* solid, SMESH_MesherHelper& helper)
+{
+  for (std::size_t i = 0; i < HEX_LINKS_NUM; ++i)
+  {
+    _Link& link = _hexLinks[i];
+
+    if (isSplittedLink(solid, helper, link))
+      continue;
+
+    // Links where both nodes are valid should have itself as a split.
+    // Check if we have at least one null node here for a chance if we have
+    // some unexpected edge case.
+    _Node* n1 = link._nodes[0];
+    _Node* n2 = link._nodes[1];
+
+    if (!n1->_node || !n2->_node)
+    {
+      link.clear();
+    }
+  }
+}
+
+//================================================================================
+/*!
+  * \brief Returns true if the link will be splitted on the Hexaedron initialization
+  */
+bool Hexahedron::isSplittedLink(const Solid* solid, SMESH_MesherHelper& helper, const Hexahedron::_Link& linkIn) const
+{
+  _Link split;
+
+  _Link link = linkIn;
+  link._fIntNodes.clear();
+  link._fIntNodes.reserve( link._fIntPoints.size() );
+
+  std::vector<_Node> intNodes;
+  intNodes.reserve(link._fIntPoints.size());
+
+  for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
+    if ( solid->ContainsAny( link._fIntPoints[i]->_faceIDs ))
+    {
+      intNodes.emplace_back(nullptr, link._fIntPoints[i]);
+      link._fIntNodes.push_back( & intNodes.back() );
+    }
+
+  link._splits.clear();
+  split._nodes[ 0 ] = link._nodes[0];
+  bool isOut = ( ! link._nodes[0]->Node() );
+  bool checkTransition = false;
+  for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
+  {
+    const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
+    if ( !isGridNode ) // intersection non-coincident with a grid node
+    {
+      if ( split._nodes[ 0 ]->Node() && !isOut )
+      {
+        return true;
+      }
+      split._nodes[ 0 ] = link._fIntNodes[i];
+      checkTransition = true;
+    }
+    else // FACE intersection coincident with a grid node (at link ends)
+    {
+      checkTransition = ( i == 0 && link._nodes[0]->Node() );
+    }
+    if ( checkTransition )
+    {
+      const vector< TGeomID >& faceIDs = link._fIntNodes[i]->_intPoint->_faceIDs;
+      if ( _grid->IsInternal( faceIDs.back() ))
+        isOut = false;
+      else if ( faceIDs.size() > 1 || _eIntPoints.size() > 0 )
+        isOut = isOutPoint( link, i, helper, solid );
+      else
+      {
+        bool okTransi = _grid->IsCorrectTransition( faceIDs[0], solid );
+        switch ( link._fIntNodes[i]->FaceIntPnt()->_transition ) {
+        case Trans_OUT: isOut = okTransi;  break;
+        case Trans_IN : isOut = !okTransi; break;
+        default:
+          isOut = isOutPoint( link, i, helper, solid );
+        }
+      }
+    }
+  }
+  if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
+  {
+    return true;
+  }
+
+  return false;
+}
+
 //================================================================================
 /*!
   * \brief Compute mesh volumes resulted from intersection of the Hexahedron
index 5676ac8563748af34ed3a1973301080560477a90..ade329d8de36a793267a9cb950ca5a45124399d0 100644 (file)
@@ -156,17 +156,22 @@ namespace Cartesian3D
       }
 
       void Add( const StdMeshers::Cartesian3D::E_IntersectPoint* ip );
+      void clear();
     };
 
     // --------------------------------------------------------------------------------
     struct _Link // link connecting two _Node's
     {
-      _Node* _nodes[2];
-      _Face* _faces[2]; // polygons sharing a link
+      static const std::size_t nodesNum = 2;
+
+      _Node* _nodes[nodesNum];
+      _Face* _faces[nodesNum]; // polygons sharing a link
       std::vector< const StdMeshers::Cartesian3D::F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
       std::vector< _Node* > _fIntNodes;   // _Node's at _fIntPoints
       std::vector< _Link >  _splits;
       _Link(): _nodes{ 0, 0 }, _faces{ 0, 0 } {}
+
+      void clear();
     };
 
     // --------------------------------------------------------------------------------
@@ -388,9 +393,12 @@ namespace Cartesian3D
     };
 
     // topology of a hexahedron
-    _Node _hexNodes [8];
-    _Link _hexLinks [12];
-    _Face _hexQuads [6];
+    static const std::size_t HEX_NODES_NUM = 8;
+    static const std::size_t HEX_LINKS_NUM = 12;
+    static const std::size_t HEX_QUADS_NUM = 6;
+    _Node _hexNodes [HEX_NODES_NUM];
+    _Link _hexLinks [HEX_LINKS_NUM];
+    _Face _hexQuads [HEX_QUADS_NUM];
 
     // faces resulted from hexahedron intersection
     std::vector< _Face > _polygons;
@@ -426,6 +434,8 @@ namespace Cartesian3D
     Hexahedron(const Hexahedron& other, size_t i, size_t j, size_t k, int cellID );
     void init( size_t i, size_t j, size_t k, const Solid* solid=0 );
     void init( size_t i );
+    void clearNodesLinkedToNull(const Solid* solid, SMESH_MesherHelper& helper);
+    bool isSplittedLink(const Solid* solid, SMESH_MesherHelper& helper, const Hexahedron::_Link& linkIn) const;
     void setIJK( size_t i );
     /*Auxiliary methods to extract operations from monolitic compute method*/
     void defineHexahedralFaces( const Solid* solid, const IsInternalFlag intFlag );