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 const int BLOCK_SIZE = 1000;
61 IMPLEMENT_STANDARD_HANDLE(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
62 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
64 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
66 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
69 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
73 HYDROData_Bathymetry::HYDROData_Bathymetry()
74 : HYDROData_IAltitudeObject()
78 HYDROData_Bathymetry::~HYDROData_Bathymetry()
82 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
83 MapOfTreatedObjects& theTreatedObjects ) const
85 QStringList aResList = dumpObjectCreation( theTreatedObjects );
86 QString aBathymetryName = GetObjPyName();
88 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
89 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
91 TCollection_AsciiString aFilePath = GetFilePath();
92 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
93 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
94 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
95 aResList << QString( "" );
96 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
97 aResList << QString( "" );
102 void HYDROData_Bathymetry::SetAltitudePoints( const HYDROData_Bathymetry::AltitudePoints& thePoints )
104 RemoveAltitudePoints();
106 if ( thePoints.empty() )
110 Handle(TDataStd_RealArray) aCoordsArray =
111 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.size() * 3 - 1 );
113 HYDROData_Bathymetry::AltitudePoints::const_iterator anIter = thePoints.begin(), aLast = thePoints.end();
114 for ( int i = 0 ; anIter!=aLast; ++i, ++anIter )
116 const HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
118 aCoordsArray->SetValue( i * 3, aPoint.X );
119 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y );
120 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z );
126 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
128 HYDROData_Bathymetry::AltitudePoints aPoints;
130 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
131 if ( aLabel.IsNull() )
134 Handle(TDataStd_RealArray) aCoordsArray;
135 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
138 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
139 int q = ( aCoordsArray->Upper() - aCoordsArray->Lower() + 1 ) / 3;
140 aPoints.reserve( q );
141 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
146 HYDROData_Bathymetry::AltitudePoint aPoint;
147 aPoint.X = aCoordsArray->Value( i++ );
148 aPoint.Y = aCoordsArray->Value( i++ );
149 aPoint.Z = aCoordsArray->Value( i++ );
151 if( IsConvertToGlobal )
152 aDoc->Transform( aPoint.X, aPoint.Y, aPoint.Z, false );
153 aPoints.push_back( aPoint );
159 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
161 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
164 int labkey = myLab.Tag();
165 //int altkey = aLabel.Tag();
166 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
167 // if (myQuadtree->isEmpty() )
168 if (myQuadtrees.find(labkey) == myQuadtrees.end())
170 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
171 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
172 myQuadtrees[labkey] = aQuadtree;
173 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
177 Handle(TDataStd_RealArray) aCoordsArray;
178 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
181 Nodes_3D* aListOfNodes = new Nodes_3D();
184 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
189 double x = aCoordsArray->Value(i++);
190 double y = aCoordsArray->Value(i++);
191 double z = aCoordsArray->Value(i++);
192 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
194 aListOfNodes->push_back(aPoint);
196 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
197 aQuadtree->setNodesAndCompute(aListOfNodes);
201 return myQuadtrees[labkey];
205 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
207 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
210 int labkey = myLab.Tag();
211 //int altkey = aLabel.Tag();
212 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
213 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
215 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
217 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
220 Handle(TDataStd_RealArray) aCoordsArray;
221 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
224 vtkPoints *points = vtkPoints::New();
225 points->Allocate(aCoordsArray->Upper() +1);
226 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
230 double x = aCoordsArray->Value(i++);
231 double y = aCoordsArray->Value(i++);
232 double z = aCoordsArray->Value(i++);
233 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
234 //DEBTRACE(" " << index);
236 vtkPolyData* profile = vtkPolyData::New();
237 profile->SetPoints(points);
238 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
240 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
241 delaunay2D->SetInputData(profile);
242 delaunay2D->Update();
243 vtkPolyData* data = delaunay2D->GetOutput();
245 myDelaunay2D[labkey] = data;
249 return myDelaunay2D[labkey];
253 void HYDROData_Bathymetry::RemoveAltitudePoints()
255 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
256 if ( !aLabel.IsNull() )
258 aLabel.ForgetAllAttributes();
263 void interpolateAltitudeForPoints( const gp_XY& thePoint,
264 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
265 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
266 HYDROData_Bathymetry::AltitudePoint& theResPoint,
267 const bool& theIsVertical )
269 double aCoordX = thePoint.X();
270 double aCoordY = thePoint.Y();
274 aCoordX = theFirstPoint.X;
276 if ( !ValuesEquals( theFirstPoint.X, theSecPoint.X ) )
278 // Recalculate X coordinate by equation of line from two points
279 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y ) * ( theSecPoint.X - theFirstPoint.X ) ) /
280 ( theSecPoint.Y - theFirstPoint.Y ) ) + theFirstPoint.X;
285 aCoordY = theFirstPoint.Y;
287 if ( !ValuesEquals( theFirstPoint.Y, theSecPoint.Y ) )
289 // Recalculate y by equation of line from two points
290 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X ) * ( theSecPoint.Y - theFirstPoint.Y ) ) /
291 ( theSecPoint.X - theFirstPoint.X ) ) + theFirstPoint.Y;
295 theResPoint.X = aCoordX;
296 theResPoint.Y = aCoordY;
298 // Calculate coefficient for interpolation
299 double aLength = Sqrt( Pow( theSecPoint.Y - theFirstPoint.Y, 2 ) +
300 Pow( theSecPoint.X - theFirstPoint.X, 2 ) );
302 double aInterCoeff = 0;
304 aInterCoeff = ( theSecPoint.Z - theFirstPoint.Z ) / aLength;
307 double aNewLength = Sqrt( Pow( theResPoint.Y - theFirstPoint.Y, 2 ) +
308 Pow( theResPoint.X - theFirstPoint.X, 2 ) );
310 // Calculate interpolated value
311 double aResVal = theFirstPoint.Z + aInterCoeff * aNewLength;
313 theResPoint.Z = aResVal;
316 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
319 int nbPts = triangle->GetNumberOfIds();
322 DEBTRACE("not a triangle ?");
326 double v[3][3]; // v[i][j] = j coordinate of node i
327 for (int i=0; i<3; i++)
329 s[i] = triangle->GetId(i);
330 delaunay2D->GetPoint(s[i],v[i]);
332 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
333 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
334 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
335 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
337 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
338 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
339 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]);
342 DEBTRACE("flat triangle ?");
346 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
347 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
350 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
351 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
354 double l2 = 1 -l0 -l1;
355 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
357 if ((l0>=0) && (l0<=1) && (l1>=0) && (l1<=1) && (l2>=0) && (l2<=1))
359 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
366 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
368 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
369 double anInvalidAltitude = GetInvalidAltitude();
370 double aResAltitude = anInvalidAltitude;
372 // --- find the nearest point in the bathymetry cloud, with quadtree
374 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
377 DEBTRACE(" no Quadtree");
381 std::map<double, const gpi_XYZ*> dist2nodes;
382 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
383 while (dist2nodes.size() == 0)
385 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
386 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
387 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
389 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
390 aResAltitude = it->second->Z();
391 int nodeIndex = it->second->getIndex();
392 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
394 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
395 // interpolation is required.
396 // - get a Delaunay2D mesh on the bathymetry cloud,
397 // - get the triangle containing the point in the Delaunay2D mesh,
398 // - interpolate altitude
400 bool isBathyInterpolRequired = false;
402 isBathyInterpolRequired =true;
405 if (isBathyInterpolRequired)
407 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
408 vtkIdList* cells= vtkIdList::New();
410 vtkIdList* points= vtkIdList::New();
411 points->Allocate(64);
412 aDelaunay2D->GetPointCells(nodeIndex, cells);
413 vtkIdType nbCells = cells->GetNumberOfIds();
414 DEBTRACE(" triangles on nearest point: " << nbCells);
415 bool isInside = false;
416 for (int i=0; i<nbCells; i++)
418 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
420 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
424 DEBTRACE(" interpolated z: " << z);
428 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
434 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
436 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
439 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
441 TCollection_AsciiString aRes;
443 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
444 if ( !aLabel.IsNull() )
446 Handle(TDataStd_AsciiString) anAsciiStr;
447 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
448 aRes = anAsciiStr->Get();
454 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
455 const bool theIsUpdate )
457 bool anIsAltitudesInverted = IsAltitudesInverted();
458 if ( anIsAltitudesInverted == theIsInverted )
461 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
468 // Update altitude points
469 HYDROData_Bathymetry::AltitudePoints anAltitudePoints = GetAltitudePoints();
470 if ( anAltitudePoints.empty() )
473 HYDROData_Bathymetry::AltitudePoints::iterator anIter = anAltitudePoints.begin(), aLast = anAltitudePoints.end();
474 for ( ; anIter!=aLast; ++anIter )
476 HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
480 SetAltitudePoints( anAltitudePoints );
483 bool HYDROData_Bathymetry::IsAltitudesInverted() const
487 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
488 if ( !aLabel.IsNull() )
490 Handle(TDataStd_Integer) anIntVal;
491 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
492 aRes = (bool)anIntVal->Get();
498 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
500 // Try to open the file
501 QFile aFile( theFileName.ToCString() );
502 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
507 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
509 HYDROData_Bathymetry::AltitudePoints aPoints;
511 // Try to import the file
512 if ( aFileSuf == "xyz" )
513 aRes = importFromXYZFile( aFile, aPoints );
514 else if ( aFileSuf == "asc" )
515 aRes = importFromASCFile( aFile, aPoints );
521 // Convert from global to local CS
522 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
523 HYDROData_Bathymetry::AltitudePoints::iterator anIter = aPoints.begin(), aLast = aPoints.end();
524 for ( ; anIter!=aLast; ++anIter )
526 HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
527 aDoc->Transform( aPoint.X, aPoint.Y, aPoint.Z, true );
532 // Update file path and altitude points of this Bathymetry
533 SetFilePath( theFileName );
534 SetAltitudePoints( aPoints );
537 return aRes && !aPoints.empty();
540 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
541 HYDROData_Bathymetry::AltitudePoints& thePoints ) const
543 if ( !theFile.isOpen() )
546 // Strings in file is written as:
547 // 1. X(float) Y(float) Z(float)
548 // 2. X(float) Y(float) Z(float)
556 bool anIsAltitudesInverted = IsAltitudesInverted();
557 while ( !theFile.atEnd() )
559 QString aLine = theFile.readLine().simplified();
560 if ( aLine.isEmpty() )
563 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
564 if ( aValues.length() < 3 )
567 HYDROData_Bathymetry::AltitudePoint aPoint;
569 QString anX = aValues.value( 0 );
570 QString anY = aValues.value( 1 );
571 QString aZ = aValues.value( 2 );
573 bool isXOk = false, isYOk = false, isZOk = false;
575 aPoint.X = anX.toDouble( &isXOk );
576 aPoint.Y = anY.toDouble( &isYOk );
577 aPoint.Z = aZ.toDouble( &isZOk );
579 if ( !isXOk || !isYOk || !isZOk )
582 if ( HYDROData_Tool::IsNan( aPoint.X ) || HYDROData_Tool::IsInf( aPoint.X ) ||
583 HYDROData_Tool::IsNan( aPoint.Y ) || HYDROData_Tool::IsInf( aPoint.Y ) ||
584 HYDROData_Tool::IsNan( aPoint.Z ) || HYDROData_Tool::IsInf( aPoint.Z ) )
587 // Invert the z value if requested
588 if ( anIsAltitudesInverted )
589 aPoint.Z = -aPoint.Z;
591 if( thePoints.size()>=thePoints.capacity() )
592 thePoints.reserve( thePoints.size()+BLOCK_SIZE );
594 thePoints.push_back( aPoint );
599 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
600 aTimer.Show( stream );
606 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
607 HYDROData_Bathymetry::AltitudePoints& thePoints ) const
609 if ( !theFile.isOpen() )
613 QStringList aStrList;
622 aLine = theFile.readLine().simplified();
623 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
624 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
626 aNCols = aStrList[1].toInt();
628 aLine = theFile.readLine().simplified();
629 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
630 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
632 aNRows = aStrList[1].toInt();
634 aLine = theFile.readLine().simplified();
635 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
636 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
638 anXllCorner = aStrList[1].toDouble();
640 aLine = theFile.readLine().simplified();
641 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
642 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
644 anYllCorner = aStrList[1].toDouble();
646 aLine = theFile.readLine().simplified();
647 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
648 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
650 aCellSize = aStrList[1].toDouble();
652 aLine = theFile.readLine().simplified();
653 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
654 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
656 aNoDataValue = aStrList[1].toDouble();
658 bool anIsAltitudesInverted = IsAltitudesInverted();
662 while ( !theFile.atEnd() )
664 aLine = theFile.readLine().simplified();
665 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
667 aStrLength = aStrList.length();
668 if ( aStrLength == 0 )
671 if ( aStrLength != aNRows )
674 for (int j = 0; j < aNCols; j++)
676 if (aStrList[j].toDouble() != aNoDataValue)
678 HYDROData_Bathymetry::AltitudePoint aPoint;
679 aPoint.X = anXllCorner + aCellSize*(j + 0.5);
680 aPoint.Y = anYllCorner + aCellSize*(aNRows - i + 0.5);
681 aPoint.Z = aStrList[j].toDouble();
683 if ( anIsAltitudesInverted )
684 aPoint.Z = -aPoint.Z;
686 if( thePoints.size()>=thePoints.capacity() )
687 thePoints.reserve( thePoints.size()+BLOCK_SIZE );
688 thePoints.push_back(aPoint);
698 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
700 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
701 Handle_HYDROData_PolylineXY aResult =
702 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
704 if( aResult.IsNull() )
708 QString aPolylinePref = GetName() + "_Boundary";
709 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
710 aResult->SetName( aPolylineName );
712 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
714 HYDROData_Bathymetry::AltitudePoints aPoints = GetAltitudePoints();
716 HYDROData_Bathymetry::AltitudePoints::const_iterator anIter = aPoints.begin(), aLast = aPoints.end();
717 for ( ; anIter!=aLast; ++anIter )
719 const HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
721 double x = aPoint.X, y = aPoint.Y;
722 if( isFirst || x<Xmin )
724 if( isFirst || x>Xmax )
726 if( isFirst || y<Ymin )
728 if( isFirst || y>Ymax )
733 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
734 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
735 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
736 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
737 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
739 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
746 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
748 gp_XYZ aDelta( theDx, theDy, 0 );
749 HYDROData_Bathymetry::AltitudePoints aPoints = GetAltitudePoints();
750 HYDROData_Bathymetry::AltitudePoints::iterator anIter = aPoints.begin(), aLast = aPoints.end();
751 for ( int i = 0; anIter!=aLast; ++i, ++anIter )
753 HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
754 aPoint.X += aDelta.X();
755 aPoint.Y += aDelta.Y();
756 aPoint.Z += aDelta.Z();
758 SetAltitudePoints( aPoints );