]> SALOME platform Git repositories - modules/hydro.git/commitdiff
Salome HOME
implementation of algorithms for DTM
authorasl <asl@opencascade.com>
Wed, 7 Sep 2016 07:35:30 +0000 (10:35 +0300)
committerasl <asl@opencascade.com>
Wed, 7 Sep 2016 07:35:30 +0000 (10:35 +0300)
src/HYDROData/HYDROData_DTM.cxx
src/HYDROData/HYDROData_DTM.h
src/HYDRO_tests/test_HYDROData_DTM.cxx

index a80edf2202d3ddae55e07a21006e11a61d27cc06..384f672a1348247ef936ac23b5f9523756fb6397 100644 (file)
 #include <TColStd_Array1OfReal.hxx>
 #include <TColStd_Array1OfInteger.hxx>
 #include <TColgp_Array1OfPnt.hxx>
+#include <Geom2dAPI_InterCurveCurve.hxx>
+#include <limits>
 
 IMPLEMENT_STANDARD_HANDLE( HYDROData_DTM, HYDROData_Bathymetry )
 IMPLEMENT_STANDARD_RTTIEXT( HYDROData_DTM, HYDROData_Bathymetry )
 
+
+HYDROData_DTM::CurveUZ::CurveUZ( double theXCurv )
+  : myXcurv( theXCurv )
+{
+}
+
+HYDROData_DTM::CurveUZ::~CurveUZ()
+{
+}
+
+double HYDROData_DTM::CurveUZ::Xcurv() const
+{
+  return myXcurv;
+}
+
+HYDROData_DTM::CurveUZ HYDROData_DTM::CurveUZ::operator + ( const CurveUZ& c ) const
+{
+  HYDROData_DTM::CurveUZ res( Xcurv() + c.Xcurv() );
+  size_t n = size();
+  res.reserve( n );
+  for( int i=0; i<n; i++ )
+  {
+    PointUZ p;
+    p.U = operator[]( i ).U + c[i].U;
+    p.Z = operator[]( i ).Z;
+    res.push_back( p );
+  }
+  return res;
+}
+
+HYDROData_DTM::CurveUZ HYDROData_DTM::CurveUZ::operator * ( double d ) const
+{
+  HYDROData_DTM::CurveUZ res( Xcurv()*d );
+  size_t n = size();
+  res.reserve( n );
+  for( int i=0; i<n; i++ )
+  {
+    PointUZ p;
+    p.U = operator[]( i ).U * d;
+    p.Z = operator[]( i ).Z;
+    res.push_back( p );
+  }
+  return res;
+}
+
+
+
 HYDROData_DTM::HYDROData_DTM()
 {
 }
@@ -35,9 +84,9 @@ HYDROData_DTM::~HYDROData_DTM()
 
 
 void GetProperties( const Handle_HYDROData_Profile& theProfile,
-                    gp_Pnt& theLowestPoint,
-                    gp_Vec2d& theDir,
-                    bool isNormalDir )
+                    gp_Pnt& theLowestPoint, gp_Vec2d& theDir,
+                    bool isNormalDir,
+                    double& theZMin, double& theZMax )
 {
   theLowestPoint = theProfile->GetBottomPoint();
   
@@ -52,24 +101,34 @@ void GetProperties( const Handle_HYDROData_Profile& theProfile,
     theDir = gp_Vec2d( x, y );
 }
 
-inline gp_Pnt2d To2D( const gp_Pnt& thePnt, const gp_Trsf& theTr )
+inline gp_Pnt2d To2D( const gp_Pnt& thePnt, const gp_Trsf& theTr,
+                      double& theUMin, double& theUMax )
 {
   gp_Pnt p = thePnt.Transformed( theTr );
-  return gp_Pnt2d( p.X(), p.Z() );
+  double u = p.X();
+  double z = p.Z();
+  if( u<theUMin )
+    theUMin = u;
+  if( u>theUMax )
+    theUMax = u;
+  return gp_Pnt2d( u, z );
 }
 
-Handle(TColgp_HArray1OfPnt2d) To2D( const TColgp_Array1OfPnt& thePoints, const gp_Trsf& theTr )
+Handle(TColgp_HArray1OfPnt2d) To2D( const TColgp_Array1OfPnt& thePoints,
+                                    const gp_Trsf& theTr,
+                                    double& theUMin, double& theUMax )
 {
   int low = thePoints.Lower(), up = thePoints.Upper();
   Handle(TColgp_HArray1OfPnt2d) points = new TColgp_HArray1OfPnt2d( low, up );
   for( int i=low; i<=up; i++ )
-    points->SetValue( i, To2D( thePoints.Value( i ), theTr ) );
+    points->SetValue( i, To2D( thePoints.Value( i ), theTr, theUMin, theUMax ) );
   return points;
 }
 
 Handle(Geom2d_Curve) CurveTo2D( const Handle(Geom_Curve)& theCurve, 
                                 Standard_Real theFirst, Standard_Real theLast,
-                                const gp_Trsf& theTr )
+                                const gp_Trsf& theTr,
+                                double& theUMin, double& theUMax )
 {
   if( theCurve->IsKind( STANDARD_TYPE( Geom_Line ) ) )
   {
@@ -78,8 +137,8 @@ Handle(Geom2d_Curve) CurveTo2D( const Handle(Geom_Curve)& theCurve,
     theCurve->D0( theLast, aLastPnt );
 
     gp_Pnt2d
-      aFirst2d = To2D( aFirstPnt, theTr ),
-      aLast2d = To2D( aLastPnt, theTr );
+      aFirst2d = To2D( aFirstPnt, theTr, theUMin, theUMax ),
+      aLast2d = To2D( aLastPnt, theTr, theUMin, theUMax );
 
     gp_Vec2d dir( aFirst2d, aLast2d );
     Handle_Geom2d_Line aLine2d = new Geom2d_Line( aFirst2d, gp_Dir2d( dir.X(), dir.Y() ) );
@@ -90,7 +149,7 @@ Handle(Geom2d_Curve) CurveTo2D( const Handle(Geom_Curve)& theCurve,
   {
     Handle(Geom_BSplineCurve) aSpline = Handle(Geom_BSplineCurve)::DownCast( theCurve );
 
-    Handle(TColgp_HArray1OfPnt2d) poles = To2D( aSpline->Poles(), theTr );
+    Handle(TColgp_HArray1OfPnt2d) poles = To2D( aSpline->Poles(), theTr, theUMin, theUMax );
     const TColStd_Array1OfReal& knots = aSpline->Knots();
     const TColStd_Array1OfInteger& multiplicities = aSpline->Multiplicities();
     int aDegree = aSpline->Degree();
@@ -102,7 +161,8 @@ Handle(Geom2d_Curve) CurveTo2D( const Handle(Geom_Curve)& theCurve,
 }
 
 Handle_Geom2d_BSplineCurve HYDROData_DTM::CreateHydraulicAxis( 
-  const std::vector<Handle_HYDROData_Profile>& theProfiles )
+  const std::vector<Handle_HYDROData_Profile>& theProfiles,
+  std::vector<double>& theDistances )
 {
   size_t n = theProfiles.size();
   Handle_Geom2d_BSplineCurve aResult;
@@ -114,9 +174,12 @@ Handle_Geom2d_BSplineCurve HYDROData_DTM::CreateHydraulicAxis(
   for( size_t i = 1; i <= n; i++ )
   {
     Handle_HYDROData_Profile aProfile = theProfiles[i-1];
+    aProfile->Update();
+
     gp_Pnt aLowest;
     gp_Vec2d aTangent;
-    GetProperties( aProfile, aLowest, aTangent, true );
+    double zmin, zmax;
+    GetProperties( aProfile, aLowest, aTangent, true, zmin, zmax );
     aTangent.Normalize();
 
     points->SetValue( (int)i, gp_Pnt2d( aLowest.X(), aLowest.Y() ) );
@@ -128,19 +191,24 @@ Handle_Geom2d_BSplineCurve HYDROData_DTM::CreateHydraulicAxis(
   anInterpolator.Load( tangents, flags );
   anInterpolator.Perform();
   if( anInterpolator.IsDone() )
+  {
     aResult = anInterpolator.Curve();
+    //TODO: fill the distances vector
+  }
   return aResult;
 }
 
-std::vector<Handle_Geom2d_Curve> HYDROData_DTM::ProfileToParametric( const Handle_HYDROData_Profile& theProfile )
+std::vector<Handle_Geom2d_Curve> HYDROData_DTM::ProfileToParametric( 
+  const Handle_HYDROData_Profile& theProfile,
+  double& theUMin, double& theUMax )
 {
   std::vector<Handle_Geom2d_Curve> curves;
-  theProfile->Update();
   
   // Transformation of the coordinate systems
   gp_Pnt aLowest;
   gp_Vec2d aDir;
-  GetProperties( theProfile, aLowest, aDir, false );
+  double zmin, zmax;
+  GetProperties( theProfile, aLowest, aDir, false, zmin, zmax );
 
   gp_Ax3 aStd3d( gp_Pnt( 0, 0, 0 ), gp_Dir( 0, 0, 1 ), gp_Dir( 1, 0, 0 ) );
   gp_Ax3 aLocal( aLowest, gp_Dir( 0, 0, 1 ), gp_Dir( aDir.X(), aDir.Y(), 0 ) );
@@ -162,30 +230,195 @@ std::vector<Handle_Geom2d_Curve> HYDROData_DTM::ProfileToParametric( const Handl
     Handle(Geom_Curve) aCurve = BRep_Tool::Curve( anEdge, aLoc, aFirst, aLast );
 
     // Convert the curve to 2d CS
-    Handle(Geom2d_Curve) aCurve2d = CurveTo2D( aCurve, aFirst, aLast, aTransf );
+    Handle(Geom2d_Curve) aCurve2d = CurveTo2D( aCurve, aFirst, aLast, aTransf, theUMin, theUMax );
     if( !aCurve2d.IsNull() )
       curves.push_back( aCurve2d );
   }
   return curves;
 }
 
+
+double CalcGC( const std::vector<gp_Pnt2d>& intersections )
+{
+  double u = 0;
+  size_t n = intersections.size(); 
+  for( size_t i = 0; i < n; i++ )
+  {
+    u += intersections[i].X();
+  }
+  u /= n;
+  return u;
+}
+
+double CalcWidth( const std::vector<gp_Pnt2d>& intersections )
+{
+  double umin = std::numeric_limits<double>::max(),
+         umax = -umin;
+
+  size_t n = intersections.size();
+  if( n <= 1 )
+    return 0;
+
+  for( size_t i = 0; i < n; i++ )
+  {
+    double u = intersections[i].X();
+    if( u<umin )
+      umin = u;
+    if( u>umax )
+      umax = u;
+  }
+  return umax-umin;
+}
+
 void HYDROData_DTM::ProfileDiscretization( const Handle_HYDROData_Profile& theProfile, 
-                                           double theMinZ, double theMaxZ, double theDDZ,
+                                           double theXCurv, double theMinZ, double theMaxZ, double theDDZ,
                                            CurveUZ& theMidPointCurve,
-                                           CurveUZ& theWidthCurve )
+                                           CurveUZ& theWidthCurve,
+                                           double theTolerance )
 {
-  /*for( double z = theMinZ; z<=theMaxZ; z += theDDZ )
+  double aDblMax = std::numeric_limits<double>::max(),
+         aUMin = aDblMax,
+         aUMax = -aUMin;
+  
+  std::vector<Handle_Geom2d_Curve> curves = ProfileToParametric( theProfile, aUMin, aUMax );
+  size_t n = curves.size();
+
+  if( n==0 )
+    return;
+
+  // we add the "virtual" vertical lines to simulate the intersection with profile 
+  gp_Pnt2d aFirst, aLast;
+  curves[0]->D0( curves[0]->FirstParameter(), aFirst );
+  curves[n-1]->D0( curves[n-1]->LastParameter(), aLast );
+  Handle(Geom2d_Line) aV1 = new Geom2d_Line( aFirst, gp_Dir2d( 0, 1 ) );
+  Handle(Geom2d_TrimmedCurve) aT1 = new Geom2d_TrimmedCurve( aV1, 0.0, aDblMax );
+  
+  Handle(Geom2d_Line) aV2 = new Geom2d_Line( aLast, gp_Dir2d( 0, 1 ) );
+  Handle(Geom2d_TrimmedCurve) aT2 = new Geom2d_TrimmedCurve( aV2, 0.0, aDblMax );
+
+  curves.push_back( aT1 );
+  curves.push_back( aT2 );
+  
+  int psize = ( int )( ( theMaxZ-theMinZ ) / theDDZ + 1 );
+  theMidPointCurve = CurveUZ( theXCurv );
+  theMidPointCurve.reserve( psize );
+  theWidthCurve = CurveUZ( theXCurv );
+  theWidthCurve.reserve( psize );
+
+  n = curves.size();
+  // for each discrete value of z we search intersection with profile
+  for( double z = theMinZ; z <= theMaxZ; z += theDDZ )
   {
+    Handle(Geom2d_Line) aLine = new Geom2d_Line( gp_Pnt2d( 0, z ), gp_Dir2d( 1, 0 ) );
+    std::vector<gp_Pnt2d> intersections;
+    for( size_t i = 0; i < n; i++ )
+    {
+      Handle_Geom2d_Curve aCurve = curves[i];
+      Geom2dAPI_InterCurveCurve anIntersect( aCurve, aLine, theTolerance );
+      for( int k=0, m=anIntersect.NbPoints(); k<m; k++ )
+        intersections.push_back( anIntersect.Point( k ) );
+    }
+
+    if( intersections.size() >= 2 )
+    {
+      double u_mid = CalcGC( intersections );
+      PointUZ p_mid;
+      p_mid.U = u_mid;
+      p_mid.Z = z;
+      theMidPointCurve.push_back( p_mid );
 
-  }*/
+      double u_wid = CalcWidth( intersections );
+      PointUZ p_wid;
+      p_wid.U = u_wid;
+      p_wid.Z = z;
+      theWidthCurve.push_back( p_wid );
+    }
+  }
 }
 
 void HYDROData_DTM::Interpolate( const CurveUZ& theCurveA, const CurveUZ& theCurveB, 
-                                 int theNbSteps, std::vector<CurveUZ>& theInterpolation )
+                                 int theNbSteps, std::vector<CurveUZ>& theInterpolation,
+                                 bool isAddSecond )
 {
+  theInterpolation.clear();
+  int d = isAddSecond ? 2 : 1;
+  theInterpolation.reserve( theNbSteps+d );
+  double dt = 1.0 / double( theNbSteps + 1 );
+  double t = dt;
+  theInterpolation.push_back( theCurveA );
+  for( int i=0; i<theNbSteps; i++, t+=dt )
+  {
+    CurveUZ anInterp = theCurveA*(1-t) + theCurveB*t;
+    theInterpolation.push_back( anInterp );
+  }
+  if( isAddSecond )
+    theInterpolation.push_back( theCurveB );
 }
 
-void HYDROData_DTM::CurveTo3d( const CurveUZ& theCurve, const CurveUZ& theCurveB, 
-                               int theNbSteps, std::vector<CurveUZ>& theInterpolation )
+void HYDROData_DTM::CurveTo3d( const CurveUZ& theMidCurve, const CurveUZ& theWidthCurve,
+                               AltitudePoints& thePoints )
 {
+  //TODO
+}
+
+inline double max( double a, double b )
+{
+  if( a>b )
+    return a;
+  else
+    return b;
+}
+
+void HYDROData_DTM::Interpolate( const Handle_HYDROData_Profile& theProfileA,
+                                 double theXCurvA,
+                                 const Handle_HYDROData_Profile& theProfileB,
+                                 double theXCurvB,
+                                 double theDDZ, int theNbSteps, bool isAddSecond,
+                                 AltitudePoints& thePoints )
+{
+  double zminA, zmaxA, zminB, zmaxB;
+  gp_Pnt lowestA, lowestB;
+  gp_Vec2d dirA, dirB;
+
+  GetProperties( theProfileA, lowestA, dirA, false, zminA, zmaxA ); 
+  GetProperties( theProfileB, lowestB, dirB, false, zminB, zmaxB ); 
+
+  double zmin = max( zminA, zminB );
+  double zmax = max( zmaxA, zmaxB );
+
+  CurveUZ midA(0), midB(0);
+  CurveUZ widA(0), widB(0);
+
+  ProfileDiscretization( theProfileA, theXCurvA, zmin, zmax, theDDZ, midA, widA ); 
+  ProfileDiscretization( theProfileB, theXCurvB, zmin, zmax, theDDZ, midB, widB );
+
+  std::vector<CurveUZ> mid, wid;
+  Interpolate( midA, midB, theNbSteps, mid, isAddSecond );
+  Interpolate( widA, widB, theNbSteps, wid, isAddSecond );
+
+  size_t q = mid.size();
+  for( size_t i=0; i<q; i++ )
+    CurveTo3d( mid[i], wid[i], thePoints );
+}
+
+void HYDROData_DTM::Interpolate( const std::vector<Handle_HYDROData_Profile>& theProfiles,
+                                 double theDDZ, double theSpatialStep, AltitudePoints& thePoints )
+{
+  size_t n = theProfiles.size();
+  if( n<=1 )
+    return;
+
+  std::vector<double> distances;
+  Handle_Geom2d_BSplineCurve aHydraulicAxis = CreateHydraulicAxis( theProfiles, distances );
+  if( aHydraulicAxis.IsNull() )
+    return;
+
+  for( size_t i=0, n1=n-1; i<n1; i++ )
+  {
+    double aDistance = distances[i+1]-distances[i];
+    int aNbSteps = int(aDistance/theSpatialStep) + 1;
+    bool isAddSecond = i==n1-1;
+    Interpolate( theProfiles[i], distances[i], 
+      theProfiles[i+1], distances[i+1], theDDZ, aNbSteps, isAddSecond, thePoints );
+  }
 }
index a91357edce2107f30962c335f899839ac64120fc..49303b2a1a850529e9e7cb78f666865ee7be35b4 100644 (file)
@@ -60,27 +60,46 @@ protected:
   class CurveUZ : public std::vector<PointUZ>
   {
   public:
-    CurveUZ( double theXCurv );
+    CurveUZ( double theXcurv );
     ~CurveUZ();
 
+    double Xcurv() const;
+
+    CurveUZ operator + ( const CurveUZ& ) const;
+    CurveUZ operator * ( double ) const;
+
   private:
-    double Xcurv;
+    double myXcurv;
   };
 
-  static Handle_Geom2d_BSplineCurve CreateHydraulicAxis( const std::vector<Handle_HYDROData_Profile>& theProfiles );
+  static Handle_Geom2d_BSplineCurve CreateHydraulicAxis( 
+    const std::vector<Handle_HYDROData_Profile>& theProfiles,
+    std::vector<double>& theDistances );
 
-  static std::vector<Handle_Geom2d_Curve> ProfileToParametric( const Handle_HYDROData_Profile& theProfile );
+  static std::vector<Handle_Geom2d_Curve> ProfileToParametric( const Handle_HYDROData_Profile& theProfile,
+                                                               double& theUMin, double& theUMax );
 
-  static void ProfileDiscretization( const Handle_HYDROData_Profile& theProfile, 
-                                     double theMinZ, double theMaxZ, double theDDZ,
+  static void ProfileDiscretization( const Handle_HYDROData_Profile& theProfile,
+                                     double theXCurv, double theMinZ, double theMaxZ, double theDDZ,
                                      CurveUZ& theMidPointCurve,
-                                     CurveUZ& theWidthCurve );
+                                     CurveUZ& theWidthCurve,
+                                     double theTolerance = 1E-6 );
 
+  static void CurveTo3d( const CurveUZ& theMidCurve, const CurveUZ& theWidthCurve, AltitudePoints& thePoints );
+  
   static void Interpolate( const CurveUZ& theCurveA, const CurveUZ& theCurveB, 
-                           int theNbSteps, std::vector<CurveUZ>& theInterpolation );
-
-  static void CurveTo3d( const CurveUZ& theCurve, const CurveUZ& theCurveB, 
-                         int theNbSteps, std::vector<CurveUZ>& theInterpolation );
+                           int theNbSteps, std::vector<CurveUZ>& theInterpolation,
+                           bool isAddSecond );
+
+  static void Interpolate( const Handle_HYDROData_Profile& theProfileA,
+                           double theXCurvA,
+                           const Handle_HYDROData_Profile& theProfileB,
+                           double theXCurvB,
+                           double theDDZ, int theNbSteps, bool isAddSecond,
+                           AltitudePoints& thePoints );
+
+  static void Interpolate( const std::vector<Handle_HYDROData_Profile>& theProfiles,
+                           double theDDZ, double theSpatialStep, AltitudePoints& thePoints );
 };
 
 #endif
index d1189c74dbe9c894cdabf18ac23df92d84c3be3b..5786bef3d59e652882d332916d49297d412e2918 100644 (file)
@@ -66,11 +66,17 @@ void test_HYDROData_DTM::test_profile_conversion_to_2d()
   aProfile2->SetLeftPoint( gp_XY( 10, 10 ) );
   aProfile2->SetRightPoint( gp_XY( 20, 20 ) );
 
-  std::vector<Handle_Geom2d_Curve> curves1 = HYDROData_DTM::ProfileToParametric( aProfile1 );
-  std::vector<Handle_Geom2d_Curve> curves2 = HYDROData_DTM::ProfileToParametric( aProfile2 );
+  double aUMin1 = std::numeric_limits<double>::max(),
+         aUMax1 = -aUMin1,
+         aUMin2 = aUMin1,
+         aUMax2 = aUMax1;
+  std::vector<Handle_Geom2d_Curve> curves1 = HYDROData_DTM::ProfileToParametric( aProfile1, aUMin1, aUMax1 );
+  std::vector<Handle_Geom2d_Curve> curves2 = HYDROData_DTM::ProfileToParametric( aProfile2, aUMin2, aUMax2 );
 
   gp_Pnt2d aFirst, aLast;
   CPPUNIT_ASSERT_EQUAL( 3, (int)curves1.size() );
+  CPPUNIT_ASSERT_DOUBLES_EQUAL( -5.303, aUMin1, EPS );
+  CPPUNIT_ASSERT_DOUBLES_EQUAL(  8.839, aUMax1, EPS );
   curves1[0]->D0( curves1[0]->FirstParameter(), aFirst );
   curves1[0]->D0( curves1[0]->LastParameter(), aLast );
   CPPUNIT_ASSERT_DOUBLES_EQUAL( -5.303, aFirst.X(), EPS );
@@ -91,6 +97,8 @@ void test_HYDROData_DTM::test_profile_conversion_to_2d()
   CPPUNIT_ASSERT_DOUBLES_EQUAL(  4.0,   aLast.Y(),  EPS );
 
   CPPUNIT_ASSERT_EQUAL( 1, (int)curves2.size() );
+  CPPUNIT_ASSERT_DOUBLES_EQUAL( -5.303, aUMin2, EPS );
+  CPPUNIT_ASSERT_DOUBLES_EQUAL(  8.839, aUMax2, EPS );
   Handle(Geom2d_BSplineCurve) aBSpline = Handle(Geom2d_BSplineCurve)::DownCast( curves2[0] );
   CPPUNIT_ASSERT_EQUAL( false, (bool)aBSpline.IsNull() );
   const TColgp_Array1OfPnt2d& poles = aBSpline->Poles();