1 // Copyright (C) 2014-2015 EDF-R&D
2 // This library is free software; you can redistribute it and/or
3 // modify it under the terms of the GNU Lesser General Public
4 // License as published by the Free Software Foundation; either
5 // version 2.1 of the License, or (at your option) any later version.
7 // This library is distributed in the hope that it will be useful,
8 // but WITHOUT ANY WARRANTY; without even the implied warranty of
9 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
10 // Lesser General Public License for more details.
12 // You should have received a copy of the GNU Lesser General Public
13 // License along with this library; if not, write to the Free Software
14 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
16 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 #include "HYDROData_Bathymetry.h"
20 #include "HYDROData_Document.h"
21 #include "HYDROData_Tool.h"
22 #include "HYDROData_PolylineXY.h"
27 #include <TDataStd_RealArray.hxx>
28 #include <TDataStd_AsciiString.hxx>
29 #include <TDataStd_Integer.hxx>
36 #include <QStringList>
39 #include <vtkPoints.h>
40 #include <vtkDelaunay2D.h>
41 #include <vtkPolyData.h>
42 #include <vtkSmartPointer.h>
43 #include <vtkIdList.h>
52 #include <OSD_Timer.hxx>
56 #include "HYDRO_trace.hxx"
58 IMPLEMENT_STANDARD_HANDLE(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
59 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
61 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
63 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
64 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
68 HYDROData_Bathymetry::HYDROData_Bathymetry()
69 : HYDROData_IAltitudeObject()
73 HYDROData_Bathymetry::~HYDROData_Bathymetry()
77 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
78 MapOfTreatedObjects& theTreatedObjects ) const
80 QStringList aResList = dumpObjectCreation( theTreatedObjects );
81 QString aBathymetryName = GetObjPyName();
83 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
84 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
86 TCollection_AsciiString aFilePath = GetFilePath();
87 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
88 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
89 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
90 aResList << QString( "" );
91 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
92 aResList << QString( "" );
97 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
99 RemoveAltitudePoints();
101 if ( thePoints.IsEmpty() )
105 Handle(TDataStd_RealArray) aCoordsArray =
106 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
108 AltitudePoints::Iterator anIter( thePoints );
109 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
111 const AltitudePoint& aPoint = anIter.Value();
113 aCoordsArray->SetValue( i * 3, aPoint.X() );
114 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
115 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
121 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
123 AltitudePoints aPoints;
125 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
126 if ( aLabel.IsNull() )
129 Handle(TDataStd_RealArray) aCoordsArray;
130 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
133 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
134 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
139 AltitudePoint aPoint;
140 aPoint.SetX( aCoordsArray->Value( i++ ) );
141 aPoint.SetY( aCoordsArray->Value( i++ ) );
142 aPoint.SetZ( aCoordsArray->Value( i++ ) );
144 if( IsConvertToGlobal )
145 aDoc->Transform( aPoint, false );
146 aPoints.Append( aPoint );
152 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
154 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
157 int labkey = myLab.Tag();
158 //int altkey = aLabel.Tag();
159 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
160 // if (myQuadtree->isEmpty() )
161 if (myQuadtrees.find(labkey) == myQuadtrees.end())
163 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
164 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
165 myQuadtrees[labkey] = aQuadtree;
166 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
170 Handle(TDataStd_RealArray) aCoordsArray;
171 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
174 Nodes_3D* aListOfNodes = new Nodes_3D();
177 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
182 double x = aCoordsArray->Value(i++);
183 double y = aCoordsArray->Value(i++);
184 double z = aCoordsArray->Value(i++);
185 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
187 aListOfNodes->push_back(aPoint);
189 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
190 aQuadtree->setNodesAndCompute(aListOfNodes);
194 return myQuadtrees[labkey];
197 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
199 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
202 int labkey = myLab.Tag();
203 //int altkey = aLabel.Tag();
204 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
205 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
207 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
209 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
212 Handle(TDataStd_RealArray) aCoordsArray;
213 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
216 vtkPoints *points = vtkPoints::New();
217 points->Allocate(aCoordsArray->Upper() +1);
218 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
222 double x = aCoordsArray->Value(i++);
223 double y = aCoordsArray->Value(i++);
224 double z = aCoordsArray->Value(i++);
225 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
226 //DEBTRACE(" " << index);
228 vtkPolyData* profile = vtkPolyData::New();
229 profile->SetPoints(points);
230 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
232 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
233 delaunay2D->SetInputData(profile);
234 delaunay2D->Update();
235 vtkPolyData* data = delaunay2D->GetOutput();
237 myDelaunay2D[labkey] = data;
241 return myDelaunay2D[labkey];
245 void HYDROData_Bathymetry::RemoveAltitudePoints()
247 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
248 if ( !aLabel.IsNull() )
250 aLabel.ForgetAllAttributes();
255 void interpolateAltitudeForPoints( const gp_XY& thePoint,
256 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
257 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
258 HYDROData_Bathymetry::AltitudePoint& theResPoint,
259 const bool& theIsVertical )
261 double aCoordX = thePoint.X();
262 double aCoordY = thePoint.Y();
266 aCoordX = theFirstPoint.X();
268 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
270 // Recalculate X coordinate by equation of line from two points
271 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
272 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
277 aCoordY = theFirstPoint.Y();
279 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
281 // Recalculate y by equation of line from two points
282 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
283 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
287 theResPoint.SetX( aCoordX );
288 theResPoint.SetY( aCoordY );
290 // Calculate coefficient for interpolation
291 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
292 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
294 double aInterCoeff = 0;
296 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
299 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
300 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
302 // Calculate interpolated value
303 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
305 theResPoint.SetZ( aResVal );
308 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
311 int nbPts = triangle->GetNumberOfIds();
314 DEBTRACE("not a triangle ?");
318 double v[3][3]; // v[i][j] = j coordinate of node i
319 for (int i=0; i<3; i++)
321 s[i] = triangle->GetId(i);
322 delaunay2D->GetPoint(s[i],v[i]);
324 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
325 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
326 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
327 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
329 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
330 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
331 double det = (v[1][1]-v[2][1])*(v[0][0]-v[2][0]) + (v[2][0]-v[1][0])*(v[0][1]-v[2][1]);
334 DEBTRACE("flat triangle ?");
338 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
339 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
342 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
343 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
346 double l2 = 1 -l0 -l1;
347 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
349 if ((l0>=0) && (l0<=1) && (l1>=0) && (l1<=1) && (l2>=0) && (l2<=1))
351 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
357 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
359 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
360 double anInvalidAltitude = GetInvalidAltitude();
361 double aResAltitude = anInvalidAltitude;
363 // --- find the nearest point in the bathymetry cloud, with quadtree
365 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
368 DEBTRACE(" no Quadtree");
372 std::map<double, const gpi_XYZ*> dist2nodes;
373 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
374 while (dist2nodes.size() == 0)
376 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
377 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
378 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
380 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
381 aResAltitude = it->second->Z();
382 int nodeIndex = it->second->getIndex();
383 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
385 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
386 // interpolation is required.
387 // - get a Delaunay2D mesh on the bathymetry cloud,
388 // - get the triangle containing the point in the Delaunay2D mesh,
389 // - interpolate altitude
391 bool isBathyInterpolRequired = false;
393 isBathyInterpolRequired =true;
394 if (isBathyInterpolRequired)
396 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
397 vtkIdList* cells= vtkIdList::New();
399 vtkIdList* points= vtkIdList::New();
400 points->Allocate(64);
401 aDelaunay2D->GetPointCells(nodeIndex, cells);
402 vtkIdType nbCells = cells->GetNumberOfIds();
403 DEBTRACE(" triangles on nearest point: " << nbCells);
404 bool isInside = false;
405 for (int i=0; i<nbCells; i++)
407 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
409 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
413 DEBTRACE(" interpolated z: " << z);
417 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
422 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
424 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
427 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
429 TCollection_AsciiString aRes;
431 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
432 if ( !aLabel.IsNull() )
434 Handle(TDataStd_AsciiString) anAsciiStr;
435 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
436 aRes = anAsciiStr->Get();
442 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
443 const bool theIsUpdate )
445 bool anIsAltitudesInverted = IsAltitudesInverted();
446 if ( anIsAltitudesInverted == theIsInverted )
449 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
456 // Update altitude points
457 AltitudePoints anAltitudePoints = GetAltitudePoints();
458 if ( anAltitudePoints.IsEmpty() )
461 AltitudePoints::Iterator anIter( anAltitudePoints );
462 for ( ; anIter.More(); anIter.Next() )
464 AltitudePoint& aPoint = anIter.ChangeValue();
465 aPoint.SetZ( aPoint.Z() * -1 );
468 SetAltitudePoints( anAltitudePoints );
471 bool HYDROData_Bathymetry::IsAltitudesInverted() const
475 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
476 if ( !aLabel.IsNull() )
478 Handle(TDataStd_Integer) anIntVal;
479 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
480 aRes = (bool)anIntVal->Get();
486 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
488 // Try to open the file
489 QFile aFile( theFileName.ToCString() );
490 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
495 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
497 AltitudePoints aPoints;
499 // Try to import the file
500 if ( aFileSuf == "xyz" )
501 aRes = importFromXYZFile( aFile, aPoints );
502 else if ( aFileSuf == "asc" )
503 aRes = importFromASCFile( aFile, aPoints );
509 // Convert from global to local CS
510 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
511 AltitudePoints::Iterator anIter( aPoints );
512 for ( ; anIter.More(); anIter.Next() )
514 AltitudePoint& aPoint = anIter.ChangeValue();
515 aDoc->Transform( aPoint, true );
520 // Update file path and altitude points of this Bathymetry
521 SetFilePath( theFileName );
522 SetAltitudePoints( aPoints );
525 return aRes && !aPoints.IsEmpty();
528 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
529 AltitudePoints& thePoints ) const
531 if ( !theFile.isOpen() )
534 // Strings in file is written as:
535 // 1. X(float) Y(float) Z(float)
536 // 2. X(float) Y(float) Z(float)
544 bool anIsAltitudesInverted = IsAltitudesInverted();
545 while ( !theFile.atEnd() )
547 QString aLine = theFile.readLine().simplified();
548 if ( aLine.isEmpty() )
551 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
552 if ( aValues.length() < 3 )
555 AltitudePoint aPoint;
557 QString anX = aValues.value( 0 );
558 QString anY = aValues.value( 1 );
559 QString aZ = aValues.value( 2 );
561 bool isXOk = false, isYOk = false, isZOk = false;
563 aPoint.SetX( anX.toDouble( &isXOk ) );
564 aPoint.SetY( anY.toDouble( &isYOk ) );
565 aPoint.SetZ( aZ.toDouble( &isZOk ) );
567 if ( !isXOk || !isYOk || !isZOk )
570 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
571 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
572 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
575 // Invert the z value if requested
576 if ( anIsAltitudesInverted )
577 aPoint.SetZ( -aPoint.Z() );
579 thePoints.Append( aPoint );
584 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
585 aTimer.Show( stream );
591 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
592 AltitudePoints& thePoints ) const
594 if ( !theFile.isOpen() )
598 QStringList aStrList;
607 aLine = theFile.readLine().simplified();
608 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
609 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
611 aNCols = aStrList[1].toInt();
613 aLine = theFile.readLine().simplified();
614 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
615 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
617 aNRows = aStrList[1].toInt();
619 aLine = theFile.readLine().simplified();
620 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
621 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
623 anXllCorner = aStrList[1].toDouble();
625 aLine = theFile.readLine().simplified();
626 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
627 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
629 anYllCorner = aStrList[1].toDouble();
631 aLine = theFile.readLine().simplified();
632 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
633 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
635 aCellSize = aStrList[1].toDouble();
637 aLine = theFile.readLine().simplified();
638 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
639 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
641 aNoDataValue = aStrList[1].toDouble();
643 bool anIsAltitudesInverted = IsAltitudesInverted();
647 while ( !theFile.atEnd() )
649 aLine = theFile.readLine().simplified();
650 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
652 aStrLength = aStrList.length();
653 if ( aStrLength == 0 )
656 if ( aStrLength != aNRows )
659 for (int j = 0; j < aNCols; j++)
661 if (aStrList[j].toDouble() != aNoDataValue)
663 AltitudePoint aPoint;
664 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
665 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
666 aPoint.SetZ(aStrList[j].toDouble());
668 if ( anIsAltitudesInverted )
669 aPoint.SetZ( -aPoint.Z() );
671 thePoints.Append(aPoint);
683 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
685 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
686 Handle_HYDROData_PolylineXY aResult =
687 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
689 if( aResult.IsNull() )
693 QString aPolylinePref = GetName() + "_Boundary";
694 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
695 aResult->SetName( aPolylineName );
697 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
699 AltitudePoints aPoints = GetAltitudePoints();
701 AltitudePoints::Iterator anIter( aPoints );
702 for ( ; anIter.More(); anIter.Next() )
704 const AltitudePoint& aPoint = anIter.Value();
706 double x = aPoint.X(), y = aPoint.Y();
707 if( isFirst || x<Xmin )
709 if( isFirst || x>Xmax )
711 if( isFirst || y<Ymin )
713 if( isFirst || y>Ymax )
718 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
719 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
720 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
721 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
722 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
724 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
731 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
733 gp_XYZ aDelta( theDx, theDy, 0 );
734 AltitudePoints aPoints = GetAltitudePoints();
735 AltitudePoints::Iterator anIter( aPoints );
736 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
738 AltitudePoint& aPoint = anIter.ChangeValue();
741 SetAltitudePoints( aPoints );