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;
66 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
70 HYDROData_Bathymetry::HYDROData_Bathymetry()
71 : HYDROData_IAltitudeObject()
75 HYDROData_Bathymetry::~HYDROData_Bathymetry()
79 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
80 MapOfTreatedObjects& theTreatedObjects ) const
82 QStringList aResList = dumpObjectCreation( theTreatedObjects );
83 QString aBathymetryName = GetObjPyName();
85 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
86 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
88 TCollection_AsciiString aFilePath = GetFilePath();
89 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
90 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
91 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
92 aResList << QString( "" );
93 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
94 aResList << QString( "" );
99 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
101 RemoveAltitudePoints();
103 if ( thePoints.IsEmpty() )
107 Handle(TDataStd_RealArray) aCoordsArray =
108 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
110 AltitudePoints::Iterator anIter( thePoints );
111 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
113 const AltitudePoint& aPoint = anIter.Value();
115 aCoordsArray->SetValue( i * 3, aPoint.X() );
116 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
117 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
123 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
125 AltitudePoints aPoints;
127 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
128 if ( aLabel.IsNull() )
131 Handle(TDataStd_RealArray) aCoordsArray;
132 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
135 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
136 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
141 AltitudePoint aPoint;
142 aPoint.SetX( aCoordsArray->Value( i++ ) );
143 aPoint.SetY( aCoordsArray->Value( i++ ) );
144 aPoint.SetZ( aCoordsArray->Value( i++ ) );
146 if( IsConvertToGlobal )
147 aDoc->Transform( aPoint, false );
148 aPoints.Append( aPoint );
154 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
156 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
159 int labkey = myLab.Tag();
160 //int altkey = aLabel.Tag();
161 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
162 // if (myQuadtree->isEmpty() )
163 if (myQuadtrees.find(labkey) == myQuadtrees.end())
165 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
166 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
167 myQuadtrees[labkey] = aQuadtree;
168 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
172 Handle(TDataStd_RealArray) aCoordsArray;
173 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
176 Nodes_3D* aListOfNodes = new Nodes_3D();
179 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
184 double x = aCoordsArray->Value(i++);
185 double y = aCoordsArray->Value(i++);
186 double z = aCoordsArray->Value(i++);
187 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
189 aListOfNodes->push_back(aPoint);
191 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
192 aQuadtree->setNodesAndCompute(aListOfNodes);
196 return myQuadtrees[labkey];
200 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
202 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
205 int labkey = myLab.Tag();
206 //int altkey = aLabel.Tag();
207 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
208 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
210 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
212 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
215 Handle(TDataStd_RealArray) aCoordsArray;
216 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
219 vtkPoints *points = vtkPoints::New();
220 points->Allocate(aCoordsArray->Upper() +1);
221 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
225 double x = aCoordsArray->Value(i++);
226 double y = aCoordsArray->Value(i++);
227 double z = aCoordsArray->Value(i++);
228 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
229 //DEBTRACE(" " << index);
231 vtkPolyData* profile = vtkPolyData::New();
232 profile->SetPoints(points);
233 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
235 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
236 delaunay2D->SetInputData(profile);
237 delaunay2D->Update();
238 vtkPolyData* data = delaunay2D->GetOutput();
240 myDelaunay2D[labkey] = data;
244 return myDelaunay2D[labkey];
248 void HYDROData_Bathymetry::RemoveAltitudePoints()
250 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
251 if ( !aLabel.IsNull() )
253 aLabel.ForgetAllAttributes();
258 void interpolateAltitudeForPoints( const gp_XY& thePoint,
259 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
260 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
261 HYDROData_Bathymetry::AltitudePoint& theResPoint,
262 const bool& theIsVertical )
264 double aCoordX = thePoint.X();
265 double aCoordY = thePoint.Y();
269 aCoordX = theFirstPoint.X();
271 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
273 // Recalculate X coordinate by equation of line from two points
274 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
275 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
280 aCoordY = theFirstPoint.Y();
282 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
284 // Recalculate y by equation of line from two points
285 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
286 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
290 theResPoint.SetX( aCoordX );
291 theResPoint.SetY( aCoordY );
293 // Calculate coefficient for interpolation
294 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
295 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
297 double aInterCoeff = 0;
299 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
302 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
303 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
305 // Calculate interpolated value
306 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
308 theResPoint.SetZ( aResVal );
311 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
314 int nbPts = triangle->GetNumberOfIds();
317 DEBTRACE("not a triangle ?");
321 double v[3][3]; // v[i][j] = j coordinate of node i
322 for (int i=0; i<3; i++)
324 s[i] = triangle->GetId(i);
325 delaunay2D->GetPoint(s[i],v[i]);
327 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
328 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
329 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
330 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
332 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
333 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
334 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]);
337 DEBTRACE("flat triangle ?");
341 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
342 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
345 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
346 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
349 double l2 = 1 -l0 -l1;
350 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
352 if ((l0>=0) && (l0<=1) && (l1>=0) && (l1<=1) && (l2>=0) && (l2<=1))
354 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
361 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
363 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
364 double anInvalidAltitude = GetInvalidAltitude();
365 double aResAltitude = anInvalidAltitude;
367 // --- find the nearest point in the bathymetry cloud, with quadtree
369 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
372 DEBTRACE(" no Quadtree");
376 std::map<double, const gpi_XYZ*> dist2nodes;
377 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
378 while (dist2nodes.size() == 0)
380 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
381 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
382 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
384 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
385 aResAltitude = it->second->Z();
386 int nodeIndex = it->second->getIndex();
387 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
389 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
390 // interpolation is required.
391 // - get a Delaunay2D mesh on the bathymetry cloud,
392 // - get the triangle containing the point in the Delaunay2D mesh,
393 // - interpolate altitude
395 bool isBathyInterpolRequired = false;
397 isBathyInterpolRequired =true;
400 if (isBathyInterpolRequired)
402 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
403 vtkIdList* cells= vtkIdList::New();
405 vtkIdList* points= vtkIdList::New();
406 points->Allocate(64);
407 aDelaunay2D->GetPointCells(nodeIndex, cells);
408 vtkIdType nbCells = cells->GetNumberOfIds();
409 DEBTRACE(" triangles on nearest point: " << nbCells);
410 bool isInside = false;
411 for (int i=0; i<nbCells; i++)
413 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
415 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
419 DEBTRACE(" interpolated z: " << z);
423 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
429 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
431 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
434 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
436 TCollection_AsciiString aRes;
438 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
439 if ( !aLabel.IsNull() )
441 Handle(TDataStd_AsciiString) anAsciiStr;
442 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
443 aRes = anAsciiStr->Get();
449 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
450 const bool theIsUpdate )
452 bool anIsAltitudesInverted = IsAltitudesInverted();
453 if ( anIsAltitudesInverted == theIsInverted )
456 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
463 // Update altitude points
464 AltitudePoints anAltitudePoints = GetAltitudePoints();
465 if ( anAltitudePoints.IsEmpty() )
468 AltitudePoints::Iterator anIter( anAltitudePoints );
469 for ( ; anIter.More(); anIter.Next() )
471 AltitudePoint& aPoint = anIter.ChangeValue();
472 aPoint.SetZ( aPoint.Z() * -1 );
475 SetAltitudePoints( anAltitudePoints );
478 bool HYDROData_Bathymetry::IsAltitudesInverted() const
482 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
483 if ( !aLabel.IsNull() )
485 Handle(TDataStd_Integer) anIntVal;
486 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
487 aRes = (bool)anIntVal->Get();
493 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
495 // Try to open the file
496 QFile aFile( theFileName.ToCString() );
497 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
502 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
504 AltitudePoints aPoints;
506 // Try to import the file
507 if ( aFileSuf == "xyz" )
508 aRes = importFromXYZFile( aFile, aPoints );
509 else if ( aFileSuf == "asc" )
510 aRes = importFromASCFile( aFile, aPoints );
516 // Convert from global to local CS
517 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
518 AltitudePoints::Iterator anIter( aPoints );
519 for ( ; anIter.More(); anIter.Next() )
521 AltitudePoint& aPoint = anIter.ChangeValue();
522 aDoc->Transform( aPoint, true );
527 // Update file path and altitude points of this Bathymetry
528 SetFilePath( theFileName );
529 SetAltitudePoints( aPoints );
532 return aRes && !aPoints.IsEmpty();
535 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
536 AltitudePoints& thePoints ) const
538 if ( !theFile.isOpen() )
541 // Strings in file is written as:
542 // 1. X(float) Y(float) Z(float)
543 // 2. X(float) Y(float) Z(float)
551 bool anIsAltitudesInverted = IsAltitudesInverted();
552 while ( !theFile.atEnd() )
554 QString aLine = theFile.readLine().simplified();
555 if ( aLine.isEmpty() )
558 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
559 if ( aValues.length() < 3 )
562 AltitudePoint aPoint;
564 QString anX = aValues.value( 0 );
565 QString anY = aValues.value( 1 );
566 QString aZ = aValues.value( 2 );
568 bool isXOk = false, isYOk = false, isZOk = false;
570 aPoint.SetX( anX.toDouble( &isXOk ) );
571 aPoint.SetY( anY.toDouble( &isYOk ) );
572 aPoint.SetZ( aZ.toDouble( &isZOk ) );
574 if ( !isXOk || !isYOk || !isZOk )
577 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
578 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
579 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
582 // Invert the z value if requested
583 if ( anIsAltitudesInverted )
584 aPoint.SetZ( -aPoint.Z() );
586 thePoints.Append( aPoint );
591 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
592 aTimer.Show( stream );
598 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
599 AltitudePoints& thePoints ) const
601 if ( !theFile.isOpen() )
605 QStringList aStrList;
614 aLine = theFile.readLine().simplified();
615 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
616 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
618 aNCols = aStrList[1].toInt();
620 aLine = theFile.readLine().simplified();
621 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
622 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
624 aNRows = aStrList[1].toInt();
626 aLine = theFile.readLine().simplified();
627 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
628 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
630 anXllCorner = aStrList[1].toDouble();
632 aLine = theFile.readLine().simplified();
633 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
634 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
636 anYllCorner = aStrList[1].toDouble();
638 aLine = theFile.readLine().simplified();
639 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
640 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
642 aCellSize = aStrList[1].toDouble();
644 aLine = theFile.readLine().simplified();
645 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
646 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
648 aNoDataValue = aStrList[1].toDouble();
650 bool anIsAltitudesInverted = IsAltitudesInverted();
654 while ( !theFile.atEnd() )
656 aLine = theFile.readLine().simplified();
657 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
659 aStrLength = aStrList.length();
660 if ( aStrLength == 0 )
663 if ( aStrLength != aNRows )
666 for (int j = 0; j < aNCols; j++)
668 if (aStrList[j].toDouble() != aNoDataValue)
670 AltitudePoint aPoint;
671 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
672 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
673 aPoint.SetZ(aStrList[j].toDouble());
675 if ( anIsAltitudesInverted )
676 aPoint.SetZ( -aPoint.Z() );
678 thePoints.Append(aPoint);
690 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
692 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
693 Handle_HYDROData_PolylineXY aResult =
694 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
696 if( aResult.IsNull() )
700 QString aPolylinePref = GetName() + "_Boundary";
701 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
702 aResult->SetName( aPolylineName );
704 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
706 AltitudePoints aPoints = GetAltitudePoints();
708 AltitudePoints::Iterator anIter( aPoints );
709 for ( ; anIter.More(); anIter.Next() )
711 const AltitudePoint& aPoint = anIter.Value();
713 double x = aPoint.X(), y = aPoint.Y();
714 if( isFirst || x<Xmin )
716 if( isFirst || x>Xmax )
718 if( isFirst || y<Ymin )
720 if( isFirst || y>Ymax )
725 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
726 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
727 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
728 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
729 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
731 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
738 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
740 gp_XYZ aDelta( theDx, theDy, 0 );
741 AltitudePoints aPoints = GetAltitudePoints();
742 AltitudePoints::Iterator anIter( aPoints );
743 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
745 AltitudePoint& aPoint = anIter.ChangeValue();
748 SetAltitudePoints( aPoints );