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_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
58 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
59 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
60 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
63 HYDROData_Bathymetry::HYDROData_Bathymetry()
64 : HYDROData_IAltitudeObject()
68 HYDROData_Bathymetry::~HYDROData_Bathymetry()
72 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
73 MapOfTreatedObjects& theTreatedObjects ) const
75 QStringList aResList = dumpObjectCreation( theTreatedObjects );
76 QString aBathymetryName = GetObjPyName();
78 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
79 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
81 TCollection_AsciiString aFilePath = GetFilePath();
82 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
83 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
84 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
85 aResList << QString( "" );
86 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
87 aResList << QString( "" );
92 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
94 RemoveAltitudePoints();
96 if ( thePoints.IsEmpty() )
100 Handle(TDataStd_RealArray) aCoordsArray =
101 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
103 AltitudePoints::Iterator anIter( thePoints );
104 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
106 const AltitudePoint& aPoint = anIter.Value();
108 aCoordsArray->SetValue( i * 3, aPoint.X() );
109 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
110 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
116 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
118 AltitudePoints aPoints;
120 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
121 if ( aLabel.IsNull() )
124 Handle(TDataStd_RealArray) aCoordsArray;
125 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
128 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
129 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
134 AltitudePoint aPoint;
135 aPoint.SetX( aCoordsArray->Value( i++ ) );
136 aPoint.SetY( aCoordsArray->Value( i++ ) );
137 aPoint.SetZ( aCoordsArray->Value( i++ ) );
139 if( IsConvertToGlobal )
140 aDoc->Transform( aPoint, false );
141 aPoints.Append( aPoint );
147 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
149 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
152 int labkey = myLab.Tag();
153 //int altkey = aLabel.Tag();
154 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
155 // if (myQuadtree->isEmpty() )
156 if (myQuadtrees.find(labkey) == myQuadtrees.end())
158 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
159 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
160 myQuadtrees[labkey] = aQuadtree;
161 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
165 Handle(TDataStd_RealArray) aCoordsArray;
166 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
169 Nodes_3D* aListOfNodes = new Nodes_3D();
172 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
177 double x = aCoordsArray->Value(i++);
178 double y = aCoordsArray->Value(i++);
179 double z = aCoordsArray->Value(i++);
180 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
182 aListOfNodes->push_back(aPoint);
184 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
185 aQuadtree->setNodesAndCompute(aListOfNodes);
189 return myQuadtrees[labkey];
192 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
194 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
197 int labkey = myLab.Tag();
198 //int altkey = aLabel.Tag();
199 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
200 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
202 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
204 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
207 Handle(TDataStd_RealArray) aCoordsArray;
208 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
211 vtkPoints *points = vtkPoints::New();
212 points->Allocate(aCoordsArray->Upper() +1);
213 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
217 double x = aCoordsArray->Value(i++);
218 double y = aCoordsArray->Value(i++);
219 double z = aCoordsArray->Value(i++);
220 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
221 //DEBTRACE(" " << index);
223 vtkPolyData* profile = vtkPolyData::New();
224 profile->SetPoints(points);
225 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
227 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
228 delaunay2D->SetInputData(profile);
229 delaunay2D->Update();
230 vtkPolyData* data = delaunay2D->GetOutput();
232 myDelaunay2D[labkey] = data;
236 return myDelaunay2D[labkey];
240 void HYDROData_Bathymetry::RemoveAltitudePoints()
242 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
243 if ( !aLabel.IsNull() )
245 aLabel.ForgetAllAttributes();
250 void interpolateAltitudeForPoints( const gp_XY& thePoint,
251 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
252 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
253 HYDROData_Bathymetry::AltitudePoint& theResPoint,
254 const bool& theIsVertical )
256 double aCoordX = thePoint.X();
257 double aCoordY = thePoint.Y();
261 aCoordX = theFirstPoint.X();
263 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
265 // Recalculate X coordinate by equation of line from two points
266 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
267 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
272 aCoordY = theFirstPoint.Y();
274 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
276 // Recalculate y by equation of line from two points
277 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
278 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
282 theResPoint.SetX( aCoordX );
283 theResPoint.SetY( aCoordY );
285 // Calculate coefficient for interpolation
286 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
287 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
289 double aInterCoeff = 0;
291 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
294 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
295 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
297 // Calculate interpolated value
298 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
300 theResPoint.SetZ( aResVal );
303 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
306 int nbPts = triangle->GetNumberOfIds();
309 DEBTRACE("not a triangle ?");
313 double v[3][3]; // v[i][j] = j coordinate of node i
314 for (int i=0; i<3; i++)
316 s[i] = triangle->GetId(i);
317 delaunay2D->GetPoint(s[i],v[i]);
319 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
320 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
321 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
322 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
324 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
325 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
326 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]);
329 DEBTRACE("flat triangle ?");
333 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
334 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
337 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
338 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
341 double l2 = 1 -l0 -l1;
342 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
344 if ((l0>=0) and (l0<=1) and (l1>=0) and (l1<=1) and (l2>=0) and (l2<=1))
346 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
352 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
354 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
355 double anInvalidAltitude = GetInvalidAltitude();
356 double aResAltitude = anInvalidAltitude;
358 // --- find the nearest point in the bathymetry cloud, with quadtree
360 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
363 DEBTRACE(" no Quadtree");
367 std::map<double, const gpi_XYZ*> dist2nodes;
368 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
369 while (dist2nodes.size() == 0)
371 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
372 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
373 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
375 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
376 aResAltitude = it->second->Z();
377 int nodeIndex = it->second->getIndex();
378 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
380 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
381 // interpolation is required.
382 // - get a Delaunay2D mesh on the bathymetry cloud,
383 // - get the triangle containing the point in the Delaunay2D mesh,
384 // - interpolate altitude
386 bool isBathyInterpolRequired = false;
388 isBathyInterpolRequired =true;
389 if (isBathyInterpolRequired)
391 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
392 vtkIdList* cells= vtkIdList::New();
394 vtkIdList* points= vtkIdList::New();
395 points->Allocate(64);
396 aDelaunay2D->GetPointCells(nodeIndex, cells);
397 vtkIdType nbCells = cells->GetNumberOfIds();
398 DEBTRACE(" triangles on nearest point: " << nbCells);
399 bool isInside = false;
400 for (int i=0; i<nbCells; i++)
402 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
404 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
408 DEBTRACE(" interpolated z: " << z);
412 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
417 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
419 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
422 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
424 TCollection_AsciiString aRes;
426 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
427 if ( !aLabel.IsNull() )
429 Handle(TDataStd_AsciiString) anAsciiStr;
430 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
431 aRes = anAsciiStr->Get();
437 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
438 const bool theIsUpdate )
440 bool anIsAltitudesInverted = IsAltitudesInverted();
441 if ( anIsAltitudesInverted == theIsInverted )
444 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
451 // Update altitude points
452 AltitudePoints anAltitudePoints = GetAltitudePoints();
453 if ( anAltitudePoints.IsEmpty() )
456 AltitudePoints::Iterator anIter( anAltitudePoints );
457 for ( ; anIter.More(); anIter.Next() )
459 AltitudePoint& aPoint = anIter.ChangeValue();
460 aPoint.SetZ( aPoint.Z() * -1 );
463 SetAltitudePoints( anAltitudePoints );
466 bool HYDROData_Bathymetry::IsAltitudesInverted() const
470 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
471 if ( !aLabel.IsNull() )
473 Handle(TDataStd_Integer) anIntVal;
474 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
475 aRes = (bool)anIntVal->Get();
481 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
483 // Try to open the file
484 QFile aFile( theFileName.ToCString() );
485 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
490 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
492 AltitudePoints aPoints;
494 // Try to import the file
495 if ( aFileSuf == "xyz" )
496 aRes = importFromXYZFile( aFile, aPoints );
497 else if ( aFileSuf == "asc" )
498 aRes = importFromASCFile( aFile, aPoints );
504 // Convert from global to local CS
505 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
506 AltitudePoints::Iterator anIter( aPoints );
507 for ( ; anIter.More(); anIter.Next() )
509 AltitudePoint& aPoint = anIter.ChangeValue();
510 aDoc->Transform( aPoint, true );
515 // Update file path and altitude points of this Bathymetry
516 SetFilePath( theFileName );
517 SetAltitudePoints( aPoints );
520 return aRes && !aPoints.IsEmpty();
523 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
524 AltitudePoints& thePoints ) const
526 if ( !theFile.isOpen() )
529 // Strings in file is written as:
530 // 1. X(float) Y(float) Z(float)
531 // 2. X(float) Y(float) Z(float)
539 bool anIsAltitudesInverted = IsAltitudesInverted();
540 while ( !theFile.atEnd() )
542 QString aLine = theFile.readLine().simplified();
543 if ( aLine.isEmpty() )
546 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
547 if ( aValues.length() < 3 )
550 AltitudePoint aPoint;
552 QString anX = aValues.value( 0 );
553 QString anY = aValues.value( 1 );
554 QString aZ = aValues.value( 2 );
556 bool isXOk = false, isYOk = false, isZOk = false;
558 aPoint.SetX( anX.toDouble( &isXOk ) );
559 aPoint.SetY( anY.toDouble( &isYOk ) );
560 aPoint.SetZ( aZ.toDouble( &isZOk ) );
562 if ( !isXOk || !isYOk || !isZOk )
565 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
566 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
567 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
570 // Invert the z value if requested
571 if ( anIsAltitudesInverted )
572 aPoint.SetZ( -aPoint.Z() );
574 thePoints.Append( aPoint );
579 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
580 aTimer.Show( stream );
586 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
587 AltitudePoints& thePoints ) const
589 if ( !theFile.isOpen() )
593 QStringList aStrList;
602 aLine = theFile.readLine().simplified();
603 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
604 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
606 aNCols = aStrList[1].toInt();
608 aLine = theFile.readLine().simplified();
609 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
610 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
612 aNRows = aStrList[1].toInt();
614 aLine = theFile.readLine().simplified();
615 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
616 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
618 anXllCorner = aStrList[1].toDouble();
620 aLine = theFile.readLine().simplified();
621 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
622 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
624 anYllCorner = aStrList[1].toDouble();
626 aLine = theFile.readLine().simplified();
627 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
628 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
630 aCellSize = aStrList[1].toDouble();
632 aLine = theFile.readLine().simplified();
633 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
634 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
636 aNoDataValue = aStrList[1].toDouble();
638 bool anIsAltitudesInverted = IsAltitudesInverted();
642 while ( !theFile.atEnd() )
644 aLine = theFile.readLine().simplified();
645 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
647 aStrLength = aStrList.length();
648 if ( aStrLength == 0 )
651 if ( aStrLength != aNRows )
654 for (int j = 0; j < aNCols; j++)
656 if (aStrList[j].toDouble() != aNoDataValue)
658 AltitudePoint aPoint;
659 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
660 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
661 aPoint.SetZ(aStrList[j].toDouble());
663 if ( anIsAltitudesInverted )
664 aPoint.SetZ( -aPoint.Z() );
666 thePoints.Append(aPoint);
678 Handle(HYDROData_PolylineXY) HYDROData_Bathymetry::CreateBoundaryPolyline() const
680 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
681 Handle(HYDROData_PolylineXY) aResult =
682 Handle(HYDROData_PolylineXY)::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
684 if( aResult.IsNull() )
688 QString aPolylinePref = GetName() + "_Boundary";
689 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
690 aResult->SetName( aPolylineName );
692 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
694 AltitudePoints aPoints = GetAltitudePoints();
696 AltitudePoints::Iterator anIter( aPoints );
697 for ( ; anIter.More(); anIter.Next() )
699 const AltitudePoint& aPoint = anIter.Value();
701 double x = aPoint.X(), y = aPoint.Y();
702 if( isFirst || x<Xmin )
704 if( isFirst || x>Xmax )
706 if( isFirst || y<Ymin )
708 if( isFirst || y>Ymax )
713 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
714 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
715 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
716 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
717 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
719 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
726 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
728 gp_XYZ aDelta( theDx, theDy, 0 );
729 AltitudePoints aPoints = GetAltitudePoints();
730 AltitudePoints::Iterator anIter( aPoints );
731 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
733 AltitudePoint& aPoint = anIter.ChangeValue();
736 SetAltitudePoints( aPoints );