Salome HOME
[TetraIntersect] Avoid multiple divisions in calculatePolygonBarycenter
[tools/medcoupling.git] / src / INTERP_KERNEL / CurveIntersector.txx
index 44df669bbc5e55a0c0a2228d7d426cad35d252f5..2864e7f8f8e1ab639c355c4bc4b27066fec4dbbf 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2014  CEA/DEN, EDF R&D
+// Copyright (C) 2007-2024  CEA, EDF
 //
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 // Author : Anthony Geay (CEA/DEN)
-#ifndef __CURVEINTERSECTOR_TXX__
-#define __CURVEINTERSECTOR_TXX__
+
+#pragma once
 
 #include "CurveIntersector.hxx"
 #include "InterpolationUtils.hxx"
+#include "PointLocatorAlgos.txx"
 
 #include <limits>
 
@@ -77,7 +78,7 @@ namespace INTERP_KERNEL
     int ibox=0;
     for(long icell=0; icell<nbelems; icell++)
       {
-        int nb_nodes_per_elem = conn_index[icell+1]-conn_index[icell];
+        ConnType nb_nodes_per_elem = conn_index[icell+1]-conn_index[icell];
         //initializing bounding box limits
         for(int idim=0; idim<SPACEDIM; idim++)
           {
@@ -85,7 +86,7 @@ namespace INTERP_KERNEL
             bbox[2*SPACEDIM*ibox+2*idim+1] = -std::numeric_limits<double>::max();
           }
         //updating the bounding box with each node of the element
-        for (int j=0; j<nb_nodes_per_elem; j++)
+        for (ConnType j=0; j<nb_nodes_per_elem; j++)
           {
             const double* coord_node = coords + 
               SPACEDIM*OTT<ConnType,numPol>
@@ -103,12 +104,9 @@ namespace INTERP_KERNEL
       }
   }
 
-  //================================================================================
   /*!
-    Computes the bouding box of a given element. iP in numPol mode.
+    Computes the bounding box of a given element. iP in numPol mode.
   */
-  //================================================================================
-
   template<class MyMeshType, class MyMatrix>
   void CurveIntersector<MyMeshType,MyMatrix>::getElemBB (double*           bb,
                                                          const MyMeshType& mesh,
@@ -138,31 +136,21 @@ namespace INTERP_KERNEL
           }
       }
   }
-
-  //================================================================================
-  /*! Readjusts a set of bounding boxes so that they are extended
-    in all dimensions for avoiding missing interesting intersections
-
-    \param bbox vector containing the bounding boxes
-  */
-  //================================================================================
-
+  
+  /*!
+   * \param [in] startOfSeg - input coming from intersectSegments or intersectSegmentsInternal
+   * \param [in] endOfSeg - input coming from intersectSegments or intersectSegmentsInternal. NO Assume about sort 
+   * \param [in] pt - position of point that the method computes the bary coords for.
+   */
   template<class MyMeshType, class MyMatrix>
-  void CurveIntersector<MyMeshType,MyMatrix>::adjustBoundingBoxes (std::vector<double>& bbox,
-                                                                   double adjustmentEpsAbs)
+  void CurveIntersector<MyMeshType,MyMatrix>::ComputeBaryCoordsOf(double startOfSeg, double endOfSeg, double pt, double& startPos, double& endPos)
   {
-    long size = bbox.size()/(2*SPACEDIM);
-    for (int i=0; i<size; i++)
-      {
-        for(int idim=0; idim<SPACEDIM; idim++)
-          {
-            bbox[i*2*SPACEDIM+2*idim  ] -= adjustmentEpsAbs;
-            bbox[i*2*SPACEDIM+2*idim+1] += adjustmentEpsAbs;
-          }
-      }
+    double deno(endOfSeg-startOfSeg);
+    startPos = (endOfSeg-pt)/deno;
+    startPos = std::max(startPos,0.); startPos = std::min(startPos,1.);
+    endPos=1.-startPos; 
   }
 
-  //================================================================================
   /*!
    * @param icellT id in target mesh in format of MyMeshType.
    * @param coordsT output val that stores coordinates of the target cell
@@ -170,14 +158,10 @@ namespace INTERP_KERNEL
    * @return true if segment is quadratic and in this case coordinates of medium node
    * are placed in the middle of coordsT
    */
-  //================================================================================
-
   template<class MyMeshType, class MyMatrix>
-  bool CurveIntersector<MyMeshType,MyMatrix>::getRealTargetCoordinates
-  (ConnType icellT, std::vector<double>& coordsT)
+  bool CurveIntersector<MyMeshType,MyMatrix>::getRealTargetCoordinates(ConnType icellT, std::vector<double>& coordsT) const
   {
-    int nbNodesT = _connIndexT[OTT<ConnType,numPol>::ind2C(icellT)+1] -
-      _connIndexT[OTT<ConnType,numPol>::ind2C(icellT)];
+    ConnType nbNodesT(_connIndexT[OTT<ConnType,numPol>::ind2C(icellT)+1] - _connIndexT[OTT<ConnType,numPol>::ind2C(icellT)]);
     coordsT.resize(SPACEDIM*nbNodesT);
     for (ConnType iT=0; iT<nbNodesT; iT++)
       {
@@ -195,22 +179,27 @@ namespace INTERP_KERNEL
       }
     return false;
   }
+  
+  template<class MyMeshType, class MyMatrix>
+  typename MyMeshType::MyConnType CurveIntersector<MyMeshType,MyMatrix>::getNodeIdOfTargetCellAt(ConnType icellT, ConnType nodeIdInCellT) const
+  {
+    ConnType nbNodesT(_connIndexT[OTT<ConnType,numPol>::ind2C(icellT)+1] - _connIndexT[OTT<ConnType,numPol>::ind2C(icellT)]);
+    if(nodeIdInCellT>=0 && nodeIdInCellT<nbNodesT)
+      return OTT<ConnType,numPol>::coo2C(_connectT[OTT<ConnType,numPol>::conn2C(_connIndexT[OTT<ConnType,numPol>::ind2C(icellT)]+nodeIdInCellT)]);
+    else
+      throw Exception("getNodeIdOfTargetCellAt : error in nodeId in cell");
+  }
 
-  //================================================================================
   /*!
    * @param icellS id in source mesh in format of MyMeshType.
    * @param coordsS output val that stores coordinates of the source cell automatically resized to the right length.
    * @return true if segment is quadratic and in this case coordinates of medium node
    * are placed in the middle of coordsS
    */
-  //================================================================================
-
   template<class MyMeshType, class MyMatrix>
-  bool CurveIntersector<MyMeshType,MyMatrix>::getRealSourceCoordinates
-  (ConnType icellS, std::vector<double>& coordsS)
+  bool CurveIntersector<MyMeshType,MyMatrix>::getRealSourceCoordinates(ConnType icellS, std::vector<double>& coordsS) const
   {
-    int nbNodesS = _connIndexS[OTT<ConnType,numPol>::ind2C(icellS)+1] -
-      _connIndexS[OTT<ConnType,numPol>::ind2C(icellS)];
+    ConnType nbNodesS = _connIndexS[OTT<ConnType,numPol>::ind2C(icellS)+1] - _connIndexS[OTT<ConnType,numPol>::ind2C(icellS)];
     coordsS.resize(SPACEDIM*nbNodesS);
     for(ConnType iS=0; iS<nbNodesS; iS++)
       {
@@ -229,24 +218,31 @@ namespace INTERP_KERNEL
     return false;
   }
 
-  //================================================================================
+  template<class MyMeshType, class MyMatrix>
+  typename MyMeshType::MyConnType CurveIntersector<MyMeshType,MyMatrix>::getNodeIdOfSourceCellAt(ConnType icellS, ConnType nodeIdInCellS) const
+  {
+    ConnType nbNodesS(_connIndexS[OTT<ConnType,numPol>::ind2C(icellS)+1] - _connIndexS[OTT<ConnType,numPol>::ind2C(icellS)]);
+    if(nodeIdInCellS>=0 && nodeIdInCellS<nbNodesS)
+      return OTT<ConnType,numPol>::coo2C(_connectS[OTT<ConnType,numPol>::conn2C(_connIndexS[OTT<ConnType,numPol>::ind2C(icellS)]+nodeIdInCellS)]);
+    else
+      throw Exception("getNodeIdOfSourceCellAt : error in nodeId in cell");
+  }
+
   /*!
    * \brief Return dual segments of given segment
    *  \param icell - given segment in C mode
    *  \param mesh - mesh
    *  \param segments - dual segments
    */
-  //================================================================================
-
   template<class MyMeshType, class MyMatrix>
   void CurveIntersector<MyMeshType,MyMatrix>::getDualSegments(ConnType                   icell,
                                                               const MyMeshType&          mesh,
                                                               std::vector<TDualSegment>& segments)
   {
     // get coordinates of cell nodes
-    int nbNodes;
-    std::vector<double> ncoords;
-    std::vector<int>    nodeIds;
+    ConnType nbNodes;
+    std::vector<double>   ncoords;
+    std::vector<ConnType> nodeIds;
     {
       const ConnType *connect   = mesh.getConnectivityPtr();
       const ConnType *connIndex = mesh.getConnectivityIndexPtr();
@@ -257,7 +253,7 @@ namespace INTERP_KERNEL
       ncoords.resize(SPACEDIM*nbNodes);
       nodeIds.resize(nbNodes);
 
-      for(int i=0; i<nbNodes; i++)
+      for(ConnType i=0; i<nbNodes; i++)
         for(int idim=0; idim<SPACEDIM; idim++)
           {
             nodeIds[i] = connect[OTT<ConnType,numPol>::conn2C(connIndex[OTT<ConnType,numPol>::ind2C(icell)]+i)];
@@ -274,7 +270,7 @@ namespace INTERP_KERNEL
     // fill segments
     segments.clear();
     segments.reserve( 2*nbNodes );
-    for(int i=0; i<nbNodes-1; i++)
+    for(ConnType i=0; i<nbNodes-1; i++)
       {
         segments.push_back(TDualSegment());
         TDualSegment& seg1 = segments.back();
@@ -299,19 +295,45 @@ namespace INTERP_KERNEL
           }
       }
   }
-
-  //================================================================================
-  /*!
-   * \brief Return length of intersection of two segments
-   */
-  //================================================================================
+  
+  template<class MyMeshType, class MyMatrix>
+  bool CurveIntersector<MyMeshType,MyMatrix>::projectionThis(const double *coordsT, const double *coordsS, double& xs0, double& xs1, double& xt) const
+  {    
+    enum { X=0, Y };
+    switch(SPACEDIM)
+      {
+      case 1:
+        {
+          xt  = coordsT[0];
+          xs0 = coordsS[0]; xs1 = coordsS[1];
+          return true;
+        }
+      case 2:
+        {
+          const double *s0(coordsS),*s1(coordsS + 2);
+          double s01[2] = { s1[X]-s0[X], s1[Y]-s0[Y] }; // src segment direction
+          double sSize = sqrt( s01[X]*s01[X] + s01[Y]*s01[Y] ); // src segment size
+          if( sSize < this->_precision )
+            return false;
+          s01[X] /= sSize; s01[Y] /= sSize; // normalize s01
+          double t[2] = { coordsT[X]-s0[X], coordsT[Y]-s0[Y] };
+          xs0 = 0. ; xs1 = sSize; xt = s01[X]*t[X] + s01[Y]*t[Y];
+          double proj_t_on_s[2] = { s0[X]+xt*s01[X], s0[Y]+xt*s01[Y] };
+          double dist_t_s_vect[2] = { coordsT[X]-proj_t_on_s[X], coordsT[Y]-proj_t_on_s[Y] };
+          double dist_t_s = sqrt( dist_t_s_vect[X]*dist_t_s_vect[X]+dist_t_s_vect[Y]*dist_t_s_vect[Y] );
+          return dist_t_s < this->_precision;
+        }
+      default:
+        throw Exception("CurveIntersector::projectionThis : space dimension of mesh must be 1 or 2");
+      }
+  }
 
   template<class MyMeshType, class MyMatrix>
-  double CurveIntersector<MyMeshType,MyMatrix>::intersectSegments(double *Coords_T,
-                                                                  double *Coords_S)
+  bool CurveIntersector<MyMeshType,MyMatrix>::projectionThis(const double *coordsT, const double *coordsS,
+                                                             double& xs0, double& xs1, double& xt0, double& xt1) const
   {
-    double xt0 = Coords_T[0], xt1 = Coords_T[1];
-    double xs0 = Coords_S[0], xs1 = Coords_S[1];
+    xt0 = coordsT[0]; xt1 = coordsT[1];
+    xs0 = coordsS[0]; xs1 = coordsS[1];
     if ( SPACEDIM == 2 )
       {
         // Pass 2D->1D
@@ -320,15 +342,16 @@ namespace INTERP_KERNEL
 
         // check if two segments overlap in 2D within tolerance
 
-        double* t0 = Coords_T;
-        double* t1 = Coords_T + 2;
+        const double* t0 = coordsT;
+        const double* t1 = coordsT + 2;
         double t01[2] = { t1[X]-t0[X], t1[Y]-t0[Y] }; // tgt segment direction
         double tSize = sqrt( t01[X]*t01[X] + t01[Y]*t01[Y] ); // tgt segment size
-        if ( tSize < _precision ) return 0; // degenerated segment
+        if ( tSize < _precision )
+          return false; // degenerated segment
         t01[X] /= tSize, t01[Y] /= tSize; // normalize t01
 
-        double* s0 = Coords_S;
-        double* s1 = Coords_S + 2;
+        const double* s0 = coordsS;
+        const double* s1 = coordsS + 2;
         double t0s0[2] = { s0[X]-t0[X], s0[Y]-t0[Y] };
         double t0s1[2] = { s1[X]-t0[X], s1[Y]-t0[Y] };
         double nt01_x_t0s0 = t0s0[X] * t01[Y] - t0s0[Y] * t01[X]; // t0s0 dot norm of t01
@@ -338,7 +361,7 @@ namespace INTERP_KERNEL
         bool s0_out_of_tol = ( dist_ts0 > _tolerance );
         bool s1_out_of_tol = ( dist_ts1 > _tolerance );
         if ( nt01_x_t0s0 * nt01_x_t0s1 > 0 && ( s0_out_of_tol || s1_out_of_tol ))
-          return 0; // tgt segment is to far from src segment
+          return false; // tgt segment is to far from src segment
 
         double S0[2] = { s0[X], s0[Y] };
         double S1[2] = { s1[X], s1[Y] };
@@ -361,7 +384,8 @@ namespace INTERP_KERNEL
 
         double s01[2] = { S1[X]-S0[X], S1[Y]-S0[Y] }; // src segment direction
         double sSize = sqrt( s01[X]*s01[X] + s01[Y]*s01[Y] ); // src segment size
-        if ( sSize < _precision ) return 0; // degenerated segment
+        if ( sSize < _precision )
+          return false; // degenerated segment
         s01[X] /= sSize, s01[Y] /= sSize; // normalize s01
 
         // make t01 and s01 codirected
@@ -375,7 +399,7 @@ namespace INTERP_KERNEL
         };
         double medianSize = sqrt( medianDir[X]*medianDir[X] + medianDir[Y]*medianDir[Y] );
         if ( medianSize < std::numeric_limits<double>::min() )
-          return 0; // strange...
+          return false; // strange...
         medianDir[X] /= medianSize, medianDir[Y] /= medianSize;
 
         xt0 = t0[X] * medianDir[X] + t0[Y] * medianDir[Y];
@@ -384,7 +408,18 @@ namespace INTERP_KERNEL
         xs1 = S1[X] * medianDir[X] + S1[Y] * medianDir[Y];
 
       } // if ( SPACEDIM == 2 )
-
+    return true;
+  }
+  
+  /*!
+   * \brief Return length of intersection of two segments
+   */
+  template<class MyMeshType, class MyMatrix>
+  double CurveIntersector<MyMeshType,MyMatrix>::intersectSegmentsInternal(const double *coordsT, const double *coordsS, double& xs0, double& xs1, double& xt0, double& xt1) const
+  {
+    if(!projectionThis(coordsT,coordsS,xs0,xs1,xt0,xt1))
+      return 0.;
+    
     if ( xt0 > xt1 ) std::swap( xt0, xt1 );
     if ( xs0 > xs1 ) std::swap( xs0, xs1 );
 
@@ -393,6 +428,46 @@ namespace INTERP_KERNEL
     return ( x0 < x1 ) ? ( x1 - x0 ) : 0.;
   }
 
-}
+  template<class MyMeshType>
+  class DummyMyMeshType1D
+  {
+  public:
+    static const int MY_SPACEDIM=1;
+    static const int MY_MESHDIM=8;
+    typedef mcIdType MyConnType;
+    static const INTERP_KERNEL::NumberingPolicy My_numPol=MyMeshType::My_numPol;
+    // begin
+    // useless, but for windows compilation ...
+    const double *getCoordinatesPtr() const { return nullptr; }
+    const MyConnType *getConnectivityPtr() const { return nullptr; }
+    const MyConnType *getConnectivityIndexPtr() const { return nullptr; }
+    INTERP_KERNEL::NormalizedCellType getTypeOfElement(MyConnType) const { return (INTERP_KERNEL::NormalizedCellType)0; }
+    // end
+  };
 
-#endif
+  /*!
+   * This method determines if a target point ( \a coordsT ) is in source seg2 contained in \a coordsS. To do so _precision attribute is used.
+   * If target point is in, \a xs0, \a xs1 and \a xt are set to 1D referential for a further barycentric computation.
+   * This method deals with SPACEDIM == 2 (see projectionThis).
+   */
+  template<class MyMeshType, class MyMatrix>
+  bool CurveIntersector<MyMeshType,MyMatrix>::isPtIncludedInSeg(const double *coordsT, const double *coordsS, double& xs0, double& xs1, double& xt) const
+  {
+    if(!projectionThis(coordsT,coordsS,xs0,xs1,xt))
+      return false;
+    constexpr ConnType TAB[2]={0,1};
+    const double coordsS_1D[2]={xs0,xs1};
+    const double *coordsT_1D(&xt);
+    return PointLocatorAlgos<DummyMyMeshType1D<MyMeshType>>::isElementContainsPoint(coordsT_1D,NORM_SEG2,coordsS_1D,TAB,2,this->_precision);
+  }
+  
+  /*!
+   * \brief Return length of intersection of two segments
+   */
+  template<class MyMeshType, class MyMatrix>
+  double CurveIntersector<MyMeshType,MyMatrix>::intersectSegments(const double *coordsT, const double *coordsS) const
+  {
+    double xs0,xs1,xt0,xt1;
+    return intersectSegmentsInternal(coordsT,coordsS,xs0,xs1,xt0,xt1);
+  }
+}