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>
30 #include <TDataStd_ExtStringArray.hxx>
37 #include <QStringList>
41 #include <vtkPoints.h>
42 #include <vtkDelaunay2D.h>
43 #include <vtkPolyData.h>
44 #include <vtkSmartPointer.h>
45 #include <vtkIdList.h>
54 #include <OSD_Timer.hxx>
58 #include "HYDRO_trace.hxx"
60 IMPLEMENT_STANDARD_HANDLE(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
61 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
63 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
65 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
68 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
72 HYDROData_Bathymetry::HYDROData_Bathymetry()
73 : HYDROData_IAltitudeObject()
77 HYDROData_Bathymetry::~HYDROData_Bathymetry()
81 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
82 MapOfTreatedObjects& theTreatedObjects ) const
84 QStringList aResList = dumpObjectCreation( theTreatedObjects );
85 QString aBathymetryName = GetObjPyName();
87 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
88 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
90 TCollection_AsciiString aFilePath = GetFilePath();
91 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
92 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
93 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
94 aResList << QString( "" );
95 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
96 aResList << QString( "" );
101 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
103 RemoveAltitudePoints();
105 if ( thePoints.IsEmpty() )
109 Handle(TDataStd_RealArray) aCoordsArray =
110 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
112 AltitudePoints::Iterator anIter( thePoints );
113 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
115 const AltitudePoint& aPoint = anIter.Value();
117 aCoordsArray->SetValue( i * 3, aPoint.X() );
118 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
119 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
125 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
127 AltitudePoints aPoints;
129 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
130 if ( aLabel.IsNull() )
133 Handle(TDataStd_RealArray) aCoordsArray;
134 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
137 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
138 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
143 AltitudePoint aPoint;
144 aPoint.SetX( aCoordsArray->Value( i++ ) );
145 aPoint.SetY( aCoordsArray->Value( i++ ) );
146 aPoint.SetZ( aCoordsArray->Value( i++ ) );
148 if( IsConvertToGlobal )
149 aDoc->Transform( aPoint, false );
150 aPoints.Append( aPoint );
156 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
158 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
161 int labkey = myLab.Tag();
162 //int altkey = aLabel.Tag();
163 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
164 // if (myQuadtree->isEmpty() )
165 if (myQuadtrees.find(labkey) == myQuadtrees.end())
167 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
168 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
169 myQuadtrees[labkey] = aQuadtree;
170 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
174 Handle(TDataStd_RealArray) aCoordsArray;
175 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
178 Nodes_3D* aListOfNodes = new Nodes_3D();
181 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
186 double x = aCoordsArray->Value(i++);
187 double y = aCoordsArray->Value(i++);
188 double z = aCoordsArray->Value(i++);
189 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
191 aListOfNodes->push_back(aPoint);
193 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
194 aQuadtree->setNodesAndCompute(aListOfNodes);
198 return myQuadtrees[labkey];
202 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
204 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
207 int labkey = myLab.Tag();
208 //int altkey = aLabel.Tag();
209 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
210 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
212 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
214 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
217 Handle(TDataStd_RealArray) aCoordsArray;
218 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
221 vtkPoints *points = vtkPoints::New();
222 points->Allocate(aCoordsArray->Upper() +1);
223 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
227 double x = aCoordsArray->Value(i++);
228 double y = aCoordsArray->Value(i++);
229 double z = aCoordsArray->Value(i++);
230 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
231 //DEBTRACE(" " << index);
233 vtkPolyData* profile = vtkPolyData::New();
234 profile->SetPoints(points);
235 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
237 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
238 delaunay2D->SetInputData(profile);
239 delaunay2D->Update();
240 vtkPolyData* data = delaunay2D->GetOutput();
242 myDelaunay2D[labkey] = data;
246 return myDelaunay2D[labkey];
250 void HYDROData_Bathymetry::RemoveAltitudePoints()
252 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
253 if ( !aLabel.IsNull() )
255 aLabel.ForgetAllAttributes();
260 void interpolateAltitudeForPoints( const gp_XY& thePoint,
261 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
262 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
263 HYDROData_Bathymetry::AltitudePoint& theResPoint,
264 const bool& theIsVertical )
266 double aCoordX = thePoint.X();
267 double aCoordY = thePoint.Y();
271 aCoordX = theFirstPoint.X();
273 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
275 // Recalculate X coordinate by equation of line from two points
276 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
277 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
282 aCoordY = theFirstPoint.Y();
284 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
286 // Recalculate y by equation of line from two points
287 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
288 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
292 theResPoint.SetX( aCoordX );
293 theResPoint.SetY( aCoordY );
295 // Calculate coefficient for interpolation
296 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
297 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
299 double aInterCoeff = 0;
301 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
304 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
305 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
307 // Calculate interpolated value
308 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
310 theResPoint.SetZ( aResVal );
313 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
316 int nbPts = triangle->GetNumberOfIds();
319 DEBTRACE("not a triangle ?");
323 double v[3][3]; // v[i][j] = j coordinate of node i
324 for (int i=0; i<3; i++)
326 s[i] = triangle->GetId(i);
327 delaunay2D->GetPoint(s[i],v[i]);
329 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
330 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
331 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
332 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
334 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
335 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
336 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]);
339 DEBTRACE("flat triangle ?");
343 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
344 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
347 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
348 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
351 double l2 = 1 -l0 -l1;
352 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
354 if ((l0>=0) && (l0<=1) && (l1>=0) && (l1<=1) && (l2>=0) && (l2<=1))
356 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
363 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
365 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
366 double anInvalidAltitude = GetInvalidAltitude();
367 double aResAltitude = anInvalidAltitude;
369 // --- find the nearest point in the bathymetry cloud, with quadtree
371 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
374 DEBTRACE(" no Quadtree");
378 std::map<double, const gpi_XYZ*> dist2nodes;
379 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
380 while (dist2nodes.size() == 0)
382 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
383 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
384 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
386 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
387 aResAltitude = it->second->Z();
388 int nodeIndex = it->second->getIndex();
389 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
391 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
392 // interpolation is required.
393 // - get a Delaunay2D mesh on the bathymetry cloud,
394 // - get the triangle containing the point in the Delaunay2D mesh,
395 // - interpolate altitude
397 bool isBathyInterpolRequired = false;
399 isBathyInterpolRequired =true;
402 if (isBathyInterpolRequired)
404 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
405 vtkIdList* cells= vtkIdList::New();
407 vtkIdList* points= vtkIdList::New();
408 points->Allocate(64);
409 aDelaunay2D->GetPointCells(nodeIndex, cells);
410 vtkIdType nbCells = cells->GetNumberOfIds();
411 DEBTRACE(" triangles on nearest point: " << nbCells);
412 bool isInside = false;
413 for (int i=0; i<nbCells; i++)
415 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
417 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
421 DEBTRACE(" interpolated z: " << z);
425 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
431 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
433 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
436 void HYDROData_Bathymetry::SetFilePaths( const QStringList& theFilePaths )
439 Handle_TDataStd_ExtStringArray TExtStrArr = TDataStd_ExtStringArray::Set( myLab.FindChild( DataTag_FilePaths ), 1, theFilePaths.size() );
440 foreach (QString filepath, theFilePaths)
442 std::string sstr = filepath.toStdString();
443 const char* Val = sstr.c_str();
444 TExtStrArr->SetValue(i, TCollection_ExtendedString(Val));
449 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
451 TCollection_AsciiString aRes;
453 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
454 if ( !aLabel.IsNull() )
456 Handle(TDataStd_AsciiString) anAsciiStr;
457 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
458 aRes = anAsciiStr->Get();
462 aLabel = myLab.FindChild( DataTag_FilePaths, false );
463 if ( !aLabel.IsNull() )
465 Handle(TDataStd_ExtStringArray) anExtStrArr;
466 if ( aLabel.FindAttribute( TDataStd_ExtStringArray::GetID(), anExtStrArr ) )
467 aRes = anExtStrArr->Value(1); //try take the first; convert extstring to asciistring
474 QStringList HYDROData_Bathymetry::GetFilePaths() const
478 TDF_Label aLabel = myLab.FindChild( DataTag_FilePaths, false );
479 if ( !aLabel.IsNull() )
481 Handle(TDataStd_ExtStringArray) anExtStrArr;
482 if ( aLabel.FindAttribute( TDataStd_ExtStringArray::GetID(), anExtStrArr ) )
484 for (int i = anExtStrArr->Lower(); i <= anExtStrArr->Upper(); i++ )
486 Standard_ExtString str = anExtStrArr->Value(i).ToExtString();
487 TCollection_AsciiString aText (str);
488 aResL << QString(aText.ToCString());
492 else //backward compatibility
494 TDF_Label anOldLabel = myLab.FindChild( DataTag_FilePath, false );
495 if ( !anOldLabel.IsNull() )
497 Handle(TDataStd_AsciiString) anAsciiStr;
498 if ( anOldLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
499 aResL << QString(anAsciiStr->Get().ToCString());
506 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
507 const bool theIsUpdate )
509 bool anIsAltitudesInverted = IsAltitudesInverted();
510 if ( anIsAltitudesInverted == theIsInverted )
513 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
520 // Update altitude points
521 AltitudePoints anAltitudePoints = GetAltitudePoints();
522 if ( anAltitudePoints.IsEmpty() )
525 AltitudePoints::Iterator anIter( anAltitudePoints );
526 for ( ; anIter.More(); anIter.Next() )
528 AltitudePoint& aPoint = anIter.ChangeValue();
529 aPoint.SetZ( aPoint.Z() * -1 );
532 SetAltitudePoints( anAltitudePoints );
535 bool HYDROData_Bathymetry::IsAltitudesInverted() const
539 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
540 if ( !aLabel.IsNull() )
542 Handle(TDataStd_Integer) anIntVal;
543 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
544 aRes = (bool)anIntVal->Get();
550 bool HYDROData_Bathymetry::ImportFromFiles( const QStringList& theFileNames )
552 AltitudePoints AllPoints;
555 foreach (QString theFileName, theFileNames)
557 // Try to open the file
558 QFile aFile( theFileName );
559 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
562 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
564 AltitudePoints aPoints;
566 // Try to import the file
567 if ( aFileSuf == "xyz" )
568 Stat = importFromXYZFile( aFile, aPoints );
569 else if ( aFileSuf == "asc" )
570 Stat = importFromASCFile( aFile, aPoints );
573 continue; //ignore this points
578 AllPoints.Append(aPoints);
581 // Convert from global to local CS
582 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
583 AltitudePoints::Iterator anIter( AllPoints );
584 for ( ; anIter.More(); anIter.Next() )
586 AltitudePoint& aPoint = anIter.ChangeValue();
587 aDoc->Transform( aPoint, true );
592 // Update file path and altitude points of this Bathymetry
593 SetFilePaths (theFileNames );
594 SetAltitudePoints( AllPoints );
597 return Stat && !AllPoints.IsEmpty();
600 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
601 AltitudePoints& thePoints ) const
603 if ( !theFile.isOpen() )
606 // Strings in file is written as:
607 // 1. X(float) Y(float) Z(float)
608 // 2. X(float) Y(float) Z(float)
616 bool anIsAltitudesInverted = IsAltitudesInverted();
617 while ( !theFile.atEnd() )
619 QString aLine = theFile.readLine().simplified();
620 if ( aLine.isEmpty() )
623 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
624 if ( aValues.length() < 3 )
627 AltitudePoint aPoint;
629 QString anX = aValues.value( 0 );
630 QString anY = aValues.value( 1 );
631 QString aZ = aValues.value( 2 );
633 bool isXOk = false, isYOk = false, isZOk = false;
635 aPoint.SetX( anX.toDouble( &isXOk ) );
636 aPoint.SetY( anY.toDouble( &isYOk ) );
637 aPoint.SetZ( aZ.toDouble( &isZOk ) );
639 if ( !isXOk || !isYOk || !isZOk )
642 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
643 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
644 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
647 // Invert the z value if requested
648 if ( anIsAltitudesInverted )
649 aPoint.SetZ( -aPoint.Z() );
651 thePoints.Append( aPoint );
656 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
657 aTimer.Show( stream );
663 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
664 AltitudePoints& thePoints ) const
666 if ( !theFile.isOpen() )
670 QStringList aStrList;
679 aLine = theFile.readLine().simplified();
680 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
681 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
683 aNCols = aStrList[1].toInt();
685 aLine = theFile.readLine().simplified();
686 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
687 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
689 aNRows = aStrList[1].toInt();
691 aLine = theFile.readLine().simplified();
692 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
693 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
695 anXllCorner = aStrList[1].toDouble();
697 aLine = theFile.readLine().simplified();
698 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
699 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
701 anYllCorner = aStrList[1].toDouble();
703 aLine = theFile.readLine().simplified();
704 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
705 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
707 aCellSize = aStrList[1].toDouble();
709 aLine = theFile.readLine().simplified();
710 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
711 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
713 aNoDataValue = aStrList[1].toDouble();
715 bool anIsAltitudesInverted = IsAltitudesInverted();
719 while ( !theFile.atEnd() )
721 aLine = theFile.readLine().simplified();
722 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
724 aStrLength = aStrList.length();
725 if ( aStrLength == 0 )
728 if ( aStrLength != aNRows )
731 for (int j = 0; j < aNCols; j++)
733 if (aStrList[j].toDouble() != aNoDataValue)
735 AltitudePoint aPoint;
736 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
737 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
738 aPoint.SetZ(aStrList[j].toDouble());
740 if ( anIsAltitudesInverted )
741 aPoint.SetZ( -aPoint.Z() );
743 thePoints.Append(aPoint);
755 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
757 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
758 Handle_HYDROData_PolylineXY aResult =
759 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
761 if( aResult.IsNull() )
765 QString aPolylinePref = GetName() + "_Boundary";
766 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
767 aResult->SetName( aPolylineName );
769 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
771 AltitudePoints aPoints = GetAltitudePoints();
773 AltitudePoints::Iterator anIter( aPoints );
774 for ( ; anIter.More(); anIter.Next() )
776 const AltitudePoint& aPoint = anIter.Value();
778 double x = aPoint.X(), y = aPoint.Y();
779 if( isFirst || x<Xmin )
781 if( isFirst || x>Xmax )
783 if( isFirst || y<Ymin )
785 if( isFirst || y>Ymax )
790 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
791 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
792 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
793 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
794 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
796 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
803 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
805 gp_XYZ aDelta( theDx, theDy, 0 );
806 AltitudePoints aPoints = GetAltitudePoints();
807 AltitudePoints::Iterator anIter( aPoints );
808 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
810 AltitudePoint& aPoint = anIter.ChangeValue();
813 SetAltitudePoints( aPoints );