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"
23 #include "HYDROData_QuadtreeNode.hxx"
28 #include <TDataStd_RealArray.hxx>
29 #include <TDataStd_AsciiString.hxx>
30 #include <TDataStd_Integer.hxx>
37 #include <QStringList>
40 #include <vtkPoints.h>
41 #include <vtkDelaunay2D.h>
42 #include <vtkPolyData.h>
43 #include <vtkSmartPointer.h>
44 #include <vtkIdList.h>
53 #include <OSD_Timer.hxx>
57 #include "HYDRO_trace.hxx"
59 IMPLEMENT_STANDARD_HANDLE(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
60 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
62 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
64 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
67 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
71 HYDROData_Bathymetry::HYDROData_Bathymetry()
72 : HYDROData_IAltitudeObject()
76 HYDROData_Bathymetry::~HYDROData_Bathymetry()
80 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
81 MapOfTreatedObjects& theTreatedObjects ) const
83 QStringList aResList = dumpObjectCreation( theTreatedObjects );
84 QString aBathymetryName = GetObjPyName();
86 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
87 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
89 TCollection_AsciiString aFilePath = GetFilePath();
90 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
91 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
92 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
93 aResList << QString( "" );
94 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
95 aResList << QString( "" );
100 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
102 RemoveAltitudePoints();
104 if ( thePoints.IsEmpty() )
108 Handle(TDataStd_RealArray) aCoordsArray =
109 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
111 AltitudePoints::Iterator anIter( thePoints );
112 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
114 const AltitudePoint& aPoint = anIter.Value();
116 aCoordsArray->SetValue( i * 3, aPoint.X() );
117 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
118 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
124 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
126 AltitudePoints aPoints;
128 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
129 if ( aLabel.IsNull() )
132 Handle(TDataStd_RealArray) aCoordsArray;
133 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
136 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
137 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
142 AltitudePoint aPoint;
143 aPoint.SetX( aCoordsArray->Value( i++ ) );
144 aPoint.SetY( aCoordsArray->Value( i++ ) );
145 aPoint.SetZ( aCoordsArray->Value( i++ ) );
147 if( IsConvertToGlobal )
148 aDoc->Transform( aPoint, false );
149 aPoints.Append( aPoint );
155 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
157 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
160 int labkey = myLab.Tag();
161 //int altkey = aLabel.Tag();
162 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
163 // if (myQuadtree->isEmpty() )
164 if (myQuadtrees.find(labkey) == myQuadtrees.end())
166 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
167 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
168 myQuadtrees[labkey] = aQuadtree;
169 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
173 Handle(TDataStd_RealArray) aCoordsArray;
174 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
177 Nodes_3D* aListOfNodes = new Nodes_3D();
180 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
185 double x = aCoordsArray->Value(i++);
186 double y = aCoordsArray->Value(i++);
187 double z = aCoordsArray->Value(i++);
188 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
190 aListOfNodes->push_back(aPoint);
192 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
193 aQuadtree->setNodesAndCompute(aListOfNodes);
197 return myQuadtrees[labkey];
201 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
203 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
206 int labkey = myLab.Tag();
207 //int altkey = aLabel.Tag();
208 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
209 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
211 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
213 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
216 Handle(TDataStd_RealArray) aCoordsArray;
217 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
220 vtkPoints *points = vtkPoints::New();
221 points->Allocate(aCoordsArray->Upper() +1);
222 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
226 double x = aCoordsArray->Value(i++);
227 double y = aCoordsArray->Value(i++);
228 double z = aCoordsArray->Value(i++);
229 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
230 //DEBTRACE(" " << index);
232 vtkPolyData* profile = vtkPolyData::New();
233 profile->SetPoints(points);
234 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
236 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
237 delaunay2D->SetInputData(profile);
238 delaunay2D->Update();
239 vtkPolyData* data = delaunay2D->GetOutput();
241 myDelaunay2D[labkey] = data;
245 return myDelaunay2D[labkey];
249 void HYDROData_Bathymetry::RemoveAltitudePoints()
251 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
252 if ( !aLabel.IsNull() )
254 aLabel.ForgetAllAttributes();
259 void interpolateAltitudeForPoints( const gp_XY& thePoint,
260 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
261 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
262 HYDROData_Bathymetry::AltitudePoint& theResPoint,
263 const bool& theIsVertical )
265 double aCoordX = thePoint.X();
266 double aCoordY = thePoint.Y();
270 aCoordX = theFirstPoint.X();
272 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
274 // Recalculate X coordinate by equation of line from two points
275 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
276 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
281 aCoordY = theFirstPoint.Y();
283 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
285 // Recalculate y by equation of line from two points
286 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
287 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
291 theResPoint.SetX( aCoordX );
292 theResPoint.SetY( aCoordY );
294 // Calculate coefficient for interpolation
295 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
296 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
298 double aInterCoeff = 0;
300 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
303 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
304 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
306 // Calculate interpolated value
307 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
309 theResPoint.SetZ( aResVal );
312 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
315 int nbPts = triangle->GetNumberOfIds();
318 DEBTRACE("not a triangle ?");
322 double v[3][3]; // v[i][j] = j coordinate of node i
323 for (int i=0; i<3; i++)
325 s[i] = triangle->GetId(i);
326 delaunay2D->GetPoint(s[i],v[i]);
328 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
329 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
330 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
331 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
333 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
334 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
335 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]);
338 DEBTRACE("flat triangle ?");
342 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
343 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
346 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
347 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
350 double l2 = 1 -l0 -l1;
351 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
353 if ((l0>=0) && (l0<=1) && (l1>=0) && (l1<=1) && (l2>=0) && (l2<=1))
355 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
362 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
364 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
365 double anInvalidAltitude = GetInvalidAltitude();
366 double aResAltitude = anInvalidAltitude;
368 // --- find the nearest point in the bathymetry cloud, with quadtree
370 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
373 DEBTRACE(" no Quadtree");
377 std::map<double, const gpi_XYZ*> dist2nodes;
378 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
379 while (dist2nodes.size() == 0)
381 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
382 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
383 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
385 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
386 aResAltitude = it->second->Z();
387 int nodeIndex = it->second->getIndex();
388 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
390 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
391 // interpolation is required.
392 // - get a Delaunay2D mesh on the bathymetry cloud,
393 // - get the triangle containing the point in the Delaunay2D mesh,
394 // - interpolate altitude
396 bool isBathyInterpolRequired = false;
398 isBathyInterpolRequired =true;
401 if (isBathyInterpolRequired)
403 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
404 vtkIdList* cells= vtkIdList::New();
406 vtkIdList* points= vtkIdList::New();
407 points->Allocate(64);
408 aDelaunay2D->GetPointCells(nodeIndex, cells);
409 vtkIdType nbCells = cells->GetNumberOfIds();
410 DEBTRACE(" triangles on nearest point: " << nbCells);
411 bool isInside = false;
412 for (int i=0; i<nbCells; i++)
414 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
416 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
420 DEBTRACE(" interpolated z: " << z);
424 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
430 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
432 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
435 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
437 TCollection_AsciiString aRes;
439 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
440 if ( !aLabel.IsNull() )
442 Handle(TDataStd_AsciiString) anAsciiStr;
443 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
444 aRes = anAsciiStr->Get();
450 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
451 const bool theIsUpdate )
453 bool anIsAltitudesInverted = IsAltitudesInverted();
454 if ( anIsAltitudesInverted == theIsInverted )
457 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
464 // Update altitude points
465 AltitudePoints anAltitudePoints = GetAltitudePoints();
466 if ( anAltitudePoints.IsEmpty() )
469 AltitudePoints::Iterator anIter( anAltitudePoints );
470 for ( ; anIter.More(); anIter.Next() )
472 AltitudePoint& aPoint = anIter.ChangeValue();
473 aPoint.SetZ( aPoint.Z() * -1 );
476 SetAltitudePoints( anAltitudePoints );
479 bool HYDROData_Bathymetry::IsAltitudesInverted() const
483 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
484 if ( !aLabel.IsNull() )
486 Handle(TDataStd_Integer) anIntVal;
487 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
488 aRes = (bool)anIntVal->Get();
494 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
496 // Try to open the file
497 QFile aFile( theFileName.ToCString() );
498 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
503 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
505 AltitudePoints aPoints;
507 // Try to import the file
508 if ( aFileSuf == "xyz" )
509 aRes = importFromXYZFile( aFile, aPoints );
510 else if ( aFileSuf == "asc" )
511 aRes = importFromASCFile( aFile, aPoints );
517 // Convert from global to local CS
518 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
519 AltitudePoints::Iterator anIter( aPoints );
520 for ( ; anIter.More(); anIter.Next() )
522 AltitudePoint& aPoint = anIter.ChangeValue();
523 aDoc->Transform( aPoint, true );
528 // Update file path and altitude points of this Bathymetry
529 SetFilePath( theFileName );
530 SetAltitudePoints( aPoints );
533 return aRes && !aPoints.IsEmpty();
536 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
537 AltitudePoints& thePoints ) const
539 if ( !theFile.isOpen() )
542 // Strings in file is written as:
543 // 1. X(float) Y(float) Z(float)
544 // 2. X(float) Y(float) Z(float)
552 bool anIsAltitudesInverted = IsAltitudesInverted();
553 while ( !theFile.atEnd() )
555 QString aLine = theFile.readLine().simplified();
556 if ( aLine.isEmpty() )
559 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
560 if ( aValues.length() < 3 )
563 AltitudePoint aPoint;
565 QString anX = aValues.value( 0 );
566 QString anY = aValues.value( 1 );
567 QString aZ = aValues.value( 2 );
569 bool isXOk = false, isYOk = false, isZOk = false;
571 aPoint.SetX( anX.toDouble( &isXOk ) );
572 aPoint.SetY( anY.toDouble( &isYOk ) );
573 aPoint.SetZ( aZ.toDouble( &isZOk ) );
575 if ( !isXOk || !isYOk || !isZOk )
578 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
579 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
580 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
583 // Invert the z value if requested
584 if ( anIsAltitudesInverted )
585 aPoint.SetZ( -aPoint.Z() );
587 thePoints.Append( aPoint );
592 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
593 aTimer.Show( stream );
599 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
600 AltitudePoints& thePoints ) const
602 if ( !theFile.isOpen() )
606 QStringList aStrList;
615 aLine = theFile.readLine().simplified();
616 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
617 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
619 aNCols = aStrList[1].toInt();
621 aLine = theFile.readLine().simplified();
622 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
623 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
625 aNRows = aStrList[1].toInt();
627 aLine = theFile.readLine().simplified();
628 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
629 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
631 anXllCorner = aStrList[1].toDouble();
633 aLine = theFile.readLine().simplified();
634 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
635 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
637 anYllCorner = aStrList[1].toDouble();
639 aLine = theFile.readLine().simplified();
640 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
641 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
643 aCellSize = aStrList[1].toDouble();
645 aLine = theFile.readLine().simplified();
646 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
647 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
649 aNoDataValue = aStrList[1].toDouble();
651 bool anIsAltitudesInverted = IsAltitudesInverted();
655 while ( !theFile.atEnd() )
657 aLine = theFile.readLine().simplified();
658 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
660 aStrLength = aStrList.length();
661 if ( aStrLength == 0 )
664 if ( aStrLength != aNRows )
667 for (int j = 0; j < aNCols; j++)
669 if (aStrList[j].toDouble() != aNoDataValue)
671 AltitudePoint aPoint;
672 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
673 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
674 aPoint.SetZ(aStrList[j].toDouble());
676 if ( anIsAltitudesInverted )
677 aPoint.SetZ( -aPoint.Z() );
679 thePoints.Append(aPoint);
691 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
693 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
694 Handle_HYDROData_PolylineXY aResult =
695 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
697 if( aResult.IsNull() )
701 QString aPolylinePref = GetName() + "_Boundary";
702 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
703 aResult->SetName( aPolylineName );
705 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
707 AltitudePoints aPoints = GetAltitudePoints();
709 AltitudePoints::Iterator anIter( aPoints );
710 for ( ; anIter.More(); anIter.Next() )
712 const AltitudePoint& aPoint = anIter.Value();
714 double x = aPoint.X(), y = aPoint.Y();
715 if( isFirst || x<Xmin )
717 if( isFirst || x>Xmax )
719 if( isFirst || y<Ymin )
721 if( isFirst || y>Ymax )
726 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
727 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
728 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
729 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
730 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
732 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
739 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
741 gp_XYZ aDelta( theDx, theDy, 0 );
742 AltitudePoints aPoints = GetAltitudePoints();
743 AltitudePoints::Iterator anIter( aPoints );
744 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
746 AltitudePoint& aPoint = anIter.ChangeValue();
749 SetAltitudePoints( aPoints );