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>
38 #include <vtkPoints.h>
39 #include <vtkDelaunay2D.h>
40 #include <vtkPolyData.h>
41 #include <vtkSmartPointer.h>
42 #include <vtkIdList.h>
50 #include <OSD_Timer.hxx>
54 #include "HYDRO_trace.hxx"
56 IMPLEMENT_STANDARD_HANDLE(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
57 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
59 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
60 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
61 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
64 HYDROData_Bathymetry::HYDROData_Bathymetry()
65 : HYDROData_IAltitudeObject()
69 HYDROData_Bathymetry::~HYDROData_Bathymetry()
73 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
74 MapOfTreatedObjects& theTreatedObjects ) const
76 QStringList aResList = dumpObjectCreation( theTreatedObjects );
77 QString aBathymetryName = GetObjPyName();
79 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
80 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
82 TCollection_AsciiString aFilePath = GetFilePath();
83 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
84 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
85 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
86 aResList << QString( "" );
87 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
88 aResList << QString( "" );
93 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
95 RemoveAltitudePoints();
97 if ( thePoints.IsEmpty() )
101 Handle(TDataStd_RealArray) aCoordsArray =
102 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
104 AltitudePoints::Iterator anIter( thePoints );
105 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
107 const AltitudePoint& aPoint = anIter.Value();
109 aCoordsArray->SetValue( i * 3, aPoint.X() );
110 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
111 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
117 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
119 AltitudePoints aPoints;
121 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
122 if ( aLabel.IsNull() )
125 Handle(TDataStd_RealArray) aCoordsArray;
126 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
129 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
130 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
135 AltitudePoint aPoint;
136 aPoint.SetX( aCoordsArray->Value( i++ ) );
137 aPoint.SetY( aCoordsArray->Value( i++ ) );
138 aPoint.SetZ( aCoordsArray->Value( i++ ) );
140 if( IsConvertToGlobal )
141 aDoc->Transform( aPoint, false );
142 aPoints.Append( aPoint );
148 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
150 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
153 int labkey = myLab.Tag();
154 //int altkey = aLabel.Tag();
155 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
156 // if (myQuadtree->isEmpty() )
157 if (myQuadtrees.find(labkey) == myQuadtrees.end())
159 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
160 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
161 myQuadtrees[labkey] = aQuadtree;
162 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
166 Handle(TDataStd_RealArray) aCoordsArray;
167 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
170 Nodes_3D* aListOfNodes = new Nodes_3D();
173 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
178 double x = aCoordsArray->Value(i++);
179 double y = aCoordsArray->Value(i++);
180 double z = aCoordsArray->Value(i++);
181 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
183 aListOfNodes->push_back(aPoint);
185 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
186 aQuadtree->setNodesAndCompute(aListOfNodes);
190 return myQuadtrees[labkey];
193 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
195 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
198 int labkey = myLab.Tag();
199 //int altkey = aLabel.Tag();
200 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
201 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
203 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
205 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
208 Handle(TDataStd_RealArray) aCoordsArray;
209 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
212 vtkPoints *points = vtkPoints::New();
213 points->Allocate(aCoordsArray->Upper() +1);
214 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
218 double x = aCoordsArray->Value(i++);
219 double y = aCoordsArray->Value(i++);
220 double z = aCoordsArray->Value(i++);
221 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
222 //DEBTRACE(" " << index);
224 vtkPolyData* profile = vtkPolyData::New();
225 profile->SetPoints(points);
226 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
228 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
229 delaunay2D->SetInputData(profile);
230 delaunay2D->Update();
231 vtkPolyData* data = delaunay2D->GetOutput();
233 myDelaunay2D[labkey] = data;
237 return myDelaunay2D[labkey];
241 void HYDROData_Bathymetry::RemoveAltitudePoints()
243 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
244 if ( !aLabel.IsNull() )
246 aLabel.ForgetAllAttributes();
251 void interpolateAltitudeForPoints( const gp_XY& thePoint,
252 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
253 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
254 HYDROData_Bathymetry::AltitudePoint& theResPoint,
255 const bool& theIsVertical )
257 double aCoordX = thePoint.X();
258 double aCoordY = thePoint.Y();
262 aCoordX = theFirstPoint.X();
264 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
266 // Recalculate X coordinate by equation of line from two points
267 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
268 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
273 aCoordY = theFirstPoint.Y();
275 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
277 // Recalculate y by equation of line from two points
278 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
279 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
283 theResPoint.SetX( aCoordX );
284 theResPoint.SetY( aCoordY );
286 // Calculate coefficient for interpolation
287 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
288 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
290 double aInterCoeff = 0;
292 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
295 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
296 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
298 // Calculate interpolated value
299 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
301 theResPoint.SetZ( aResVal );
304 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
307 int nbPts = triangle->GetNumberOfIds();
310 DEBTRACE("not a triangle ?");
314 double v[3][3]; // v[i][j] = j coordinate of node i
315 for (int i=0; i<3; i++)
317 s[i] = triangle->GetId(i);
318 delaunay2D->GetPoint(s[i],v[i]);
320 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
321 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
322 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
323 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
325 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
326 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
327 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]);
330 DEBTRACE("flat triangle ?");
334 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
335 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
338 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
339 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
342 double l2 = 1 -l0 -l1;
343 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
345 if ((l0>=0) and (l0<=1) and (l1>=0) and (l1<=1) and (l2>=0) and (l2<=1))
347 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
353 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint) const
355 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << ")");
356 double anInvalidAltitude = GetInvalidAltitude();
357 double aResAltitude = anInvalidAltitude;
359 // --- find the nearest point in the bathymetry cloud, with quadtree
361 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
364 DEBTRACE(" no Quadtree");
368 std::map<double, const gpi_XYZ*> dist2nodes;
369 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
370 while (dist2nodes.size() == 0)
372 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
373 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
374 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
376 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
377 aResAltitude = it->second->Z();
378 int nodeIndex = it->second->getIndex();
379 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
381 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
382 // interpolation is required.
383 // - get a Delaunay2D mesh on the bathymetry cloud,
384 // - get the triangle containing the point in the Delaunay2D mesh,
385 // - interpolate altitude
387 bool isBathyInterpolRequired = true;
388 if (isBathyInterpolRequired)
390 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
391 vtkIdList* cells= vtkIdList::New();
393 vtkIdList* points= vtkIdList::New();
394 points->Allocate(64);
395 aDelaunay2D->GetPointCells(nodeIndex, cells);
396 vtkIdType nbCells = cells->GetNumberOfIds();
397 DEBTRACE(" triangles on nearest point: " << nbCells);
398 bool isInside = false;
399 for (int i=0; i<nbCells; i++)
401 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
403 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
407 DEBTRACE(" interpolated z: " << z);
411 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
416 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
418 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
421 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
423 TCollection_AsciiString aRes;
425 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
426 if ( !aLabel.IsNull() )
428 Handle(TDataStd_AsciiString) anAsciiStr;
429 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
430 aRes = anAsciiStr->Get();
436 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
437 const bool theIsUpdate )
439 bool anIsAltitudesInverted = IsAltitudesInverted();
440 if ( anIsAltitudesInverted == theIsInverted )
443 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
450 // Update altitude points
451 AltitudePoints anAltitudePoints = GetAltitudePoints();
452 if ( anAltitudePoints.IsEmpty() )
455 AltitudePoints::Iterator anIter( anAltitudePoints );
456 for ( ; anIter.More(); anIter.Next() )
458 AltitudePoint& aPoint = anIter.ChangeValue();
459 aPoint.SetZ( aPoint.Z() * -1 );
462 SetAltitudePoints( anAltitudePoints );
465 bool HYDROData_Bathymetry::IsAltitudesInverted() const
469 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
470 if ( !aLabel.IsNull() )
472 Handle(TDataStd_Integer) anIntVal;
473 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
474 aRes = (bool)anIntVal->Get();
480 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
482 // Try to open the file
483 QFile aFile( theFileName.ToCString() );
484 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
489 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
491 AltitudePoints aPoints;
493 // Try to import the file
494 if ( aFileSuf == "xyz" )
495 aRes = importFromXYZFile( aFile, aPoints );
496 else if ( aFileSuf == "asc" )
497 aRes = importFromASCFile( aFile, aPoints );
503 // Convert from global to local CS
504 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
505 AltitudePoints::Iterator anIter( aPoints );
506 for ( ; anIter.More(); anIter.Next() )
508 AltitudePoint& aPoint = anIter.ChangeValue();
509 aDoc->Transform( aPoint, true );
514 // Update file path and altitude points of this Bathymetry
515 SetFilePath( theFileName );
516 SetAltitudePoints( aPoints );
519 return aRes && !aPoints.IsEmpty();
522 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
523 AltitudePoints& thePoints ) const
525 if ( !theFile.isOpen() )
528 // Strings in file is written as:
529 // 1. X(float) Y(float) Z(float)
530 // 2. X(float) Y(float) Z(float)
538 bool anIsAltitudesInverted = IsAltitudesInverted();
539 while ( !theFile.atEnd() )
541 QString aLine = theFile.readLine().simplified();
542 if ( aLine.isEmpty() )
545 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
546 if ( aValues.length() < 3 )
549 AltitudePoint aPoint;
551 QString anX = aValues.value( 0 );
552 QString anY = aValues.value( 1 );
553 QString aZ = aValues.value( 2 );
555 bool isXOk = false, isYOk = false, isZOk = false;
557 aPoint.SetX( anX.toDouble( &isXOk ) );
558 aPoint.SetY( anY.toDouble( &isYOk ) );
559 aPoint.SetZ( aZ.toDouble( &isZOk ) );
561 if ( !isXOk || !isYOk || !isZOk )
564 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
565 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
566 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
569 // Invert the z value if requested
570 if ( anIsAltitudesInverted )
571 aPoint.SetZ( -aPoint.Z() );
573 thePoints.Append( aPoint );
578 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
579 aTimer.Show( stream );
585 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
586 AltitudePoints& thePoints ) const
588 if ( !theFile.isOpen() )
592 QStringList aStrList;
601 aLine = theFile.readLine().simplified();
602 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
603 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
605 aNCols = aStrList[1].toInt();
607 aLine = theFile.readLine().simplified();
608 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
609 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
611 aNRows = aStrList[1].toInt();
613 aLine = theFile.readLine().simplified();
614 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
615 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
617 anXllCorner = aStrList[1].toDouble();
619 aLine = theFile.readLine().simplified();
620 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
621 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
623 anYllCorner = aStrList[1].toDouble();
625 aLine = theFile.readLine().simplified();
626 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
627 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
629 aCellSize = aStrList[1].toDouble();
631 aLine = theFile.readLine().simplified();
632 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
633 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
635 aNoDataValue = aStrList[1].toDouble();
637 bool anIsAltitudesInverted = IsAltitudesInverted();
641 while ( !theFile.atEnd() )
643 aLine = theFile.readLine().simplified();
644 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
646 aStrLength = aStrList.length();
647 if ( aStrLength == 0 )
650 if ( aStrLength != aNRows )
653 for (int j = 0; j < aNCols; j++)
655 if (aStrList[j].toDouble() != aNoDataValue)
657 AltitudePoint aPoint;
658 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
659 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
660 aPoint.SetZ(aStrList[j].toDouble());
662 if ( anIsAltitudesInverted )
663 aPoint.SetZ( -aPoint.Z() );
665 thePoints.Append(aPoint);
677 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
679 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
680 Handle_HYDROData_PolylineXY aResult =
681 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
683 if( aResult.IsNull() )
687 QString aPolylinePref = GetName() + "_Boundary";
688 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
689 aResult->SetName( aPolylineName );
691 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
693 AltitudePoints aPoints = GetAltitudePoints();
695 AltitudePoints::Iterator anIter( aPoints );
696 for ( ; anIter.More(); anIter.Next() )
698 const AltitudePoint& aPoint = anIter.Value();
700 double x = aPoint.X(), y = aPoint.Y();
701 if( isFirst || x<Xmin )
703 if( isFirst || x>Xmax )
705 if( isFirst || y<Ymin )
707 if( isFirst || y>Ymax )
712 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
713 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
714 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
715 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
716 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
718 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
725 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
727 gp_XYZ aDelta( theDx, theDy, 0 );
728 AltitudePoints aPoints = GetAltitudePoints();
729 AltitudePoints::Iterator anIter( aPoints );
730 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
732 AltitudePoint& aPoint = anIter.ChangeValue();
735 SetAltitudePoints( aPoints );