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>
31 #include <TDataStd_ExtStringArray.hxx>
38 #include <QStringList>
42 #include <vtkPoints.h>
43 #include <vtkDelaunay2D.h>
44 #include <vtkPolyData.h>
45 #include <vtkSmartPointer.h>
46 #include <vtkIdList.h>
55 #include <OSD_Timer.hxx>
59 #include "HYDRO_trace.hxx"
61 const int BLOCK_SIZE = 1000;
63 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
65 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
67 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
70 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
73 inline double sqr( double x )
78 HYDROData_Bathymetry::AltitudePoint::AltitudePoint( double x, double y, double z )
83 double HYDROData_Bathymetry::AltitudePoint::SquareDistance( const HYDROData_Bathymetry::AltitudePoint& p ) const
92 HYDROData_Bathymetry::HYDROData_Bathymetry()
93 : HYDROData_IAltitudeObject()
97 HYDROData_Bathymetry::~HYDROData_Bathymetry()
101 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
102 MapOfTreatedObjects& theTreatedObjects ) const
104 QStringList aResList = dumpObjectCreation( theTreatedObjects );
105 QString aBathymetryName = GetObjPyName();
107 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
108 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
110 TCollection_AsciiString aFilePath = GetFilePath();
111 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
112 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
113 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
114 aResList << QString( "" );
115 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
116 aResList << QString( "" );
121 void HYDROData_Bathymetry::SetAltitudePoints( const HYDROData_Bathymetry::AltitudePoints& thePoints )
123 RemoveAltitudePoints();
125 if ( thePoints.empty() )
129 Handle(TDataStd_RealArray) aCoordsArray =
130 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.size() * 3 - 1 );
132 HYDROData_Bathymetry::AltitudePoints::const_iterator anIter = thePoints.begin(), aLast = thePoints.end();
133 for ( int i = 0 ; anIter!=aLast; ++i, ++anIter )
135 const HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
137 aCoordsArray->SetValue( i * 3, aPoint.X );
138 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y );
139 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z );
145 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
147 HYDROData_Bathymetry::AltitudePoints aPoints;
149 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
150 if ( aLabel.IsNull() )
153 Handle(TDataStd_RealArray) aCoordsArray;
154 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
157 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
158 int q = ( aCoordsArray->Upper() - aCoordsArray->Lower() + 1 ) / 3;
159 aPoints.reserve( q );
160 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
165 HYDROData_Bathymetry::AltitudePoint aPoint;
166 aPoint.X = aCoordsArray->Value( i++ );
167 aPoint.Y = aCoordsArray->Value( i++ );
168 aPoint.Z = aCoordsArray->Value( i++ );
170 if( IsConvertToGlobal )
171 aDoc->Transform( aPoint.X, aPoint.Y, aPoint.Z, false );
172 aPoints.push_back( aPoint );
178 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
180 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
183 int labkey = myLab.Tag();
184 //int altkey = aLabel.Tag();
185 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
186 // if (myQuadtree->isEmpty() )
187 if (myQuadtrees.find(labkey) == myQuadtrees.end())
189 //DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
190 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
191 myQuadtrees[labkey] = aQuadtree;
192 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
196 Handle(TDataStd_RealArray) aCoordsArray;
197 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
200 Nodes_3D* aListOfNodes = new Nodes_3D();
203 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
208 double x = aCoordsArray->Value(i++);
209 double y = aCoordsArray->Value(i++);
210 double z = aCoordsArray->Value(i++);
211 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
213 aListOfNodes->push_back(aPoint);
215 //DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
216 aQuadtree->setNodesAndCompute(aListOfNodes);
220 return myQuadtrees[labkey];
224 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
226 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
229 int labkey = myLab.Tag();
230 //int altkey = aLabel.Tag();
231 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
232 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
234 //DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
236 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
239 Handle(TDataStd_RealArray) aCoordsArray;
240 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
243 vtkPoints *points = vtkPoints::New();
244 points->Allocate(aCoordsArray->Upper() +1);
245 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
249 double x = aCoordsArray->Value(i++);
250 double y = aCoordsArray->Value(i++);
251 double z = aCoordsArray->Value(i++);
252 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
253 //DEBTRACE(" " << index);
255 vtkPolyData* profile = vtkPolyData::New();
256 profile->SetPoints(points);
257 //DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
259 vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
260 delaunay2D->SetInputData(profile);
261 delaunay2D->Update();
262 vtkPolyData* data = delaunay2D->GetOutput();
264 myDelaunay2D[labkey] = data;
268 return myDelaunay2D[labkey];
272 void HYDROData_Bathymetry::RemoveAltitudePoints()
274 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
275 if ( !aLabel.IsNull() )
277 aLabel.ForgetAllAttributes();
282 void interpolateAltitudeForPoints( const gp_XY& thePoint,
283 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
284 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
285 HYDROData_Bathymetry::AltitudePoint& theResPoint,
286 const bool& theIsVertical )
288 double aCoordX = thePoint.X();
289 double aCoordY = thePoint.Y();
293 aCoordX = theFirstPoint.X;
295 if ( !ValuesEquals( theFirstPoint.X, theSecPoint.X ) )
297 // Recalculate X coordinate by equation of line from two points
298 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y ) * ( theSecPoint.X - theFirstPoint.X ) ) /
299 ( theSecPoint.Y - theFirstPoint.Y ) ) + theFirstPoint.X;
304 aCoordY = theFirstPoint.Y;
306 if ( !ValuesEquals( theFirstPoint.Y, theSecPoint.Y ) )
308 // Recalculate y by equation of line from two points
309 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X ) * ( theSecPoint.Y - theFirstPoint.Y ) ) /
310 ( theSecPoint.X - theFirstPoint.X ) ) + theFirstPoint.Y;
314 theResPoint.X = aCoordX;
315 theResPoint.Y = aCoordY;
317 // Calculate coefficient for interpolation
318 double aLength = Sqrt( Pow( theSecPoint.Y - theFirstPoint.Y, 2 ) +
319 Pow( theSecPoint.X - theFirstPoint.X, 2 ) );
321 double aInterCoeff = 0;
323 aInterCoeff = ( theSecPoint.Z - theFirstPoint.Z ) / aLength;
326 double aNewLength = Sqrt( Pow( theResPoint.Y - theFirstPoint.Y, 2 ) +
327 Pow( theResPoint.X - theFirstPoint.X, 2 ) );
329 // Calculate interpolated value
330 double aResVal = theFirstPoint.Z + aInterCoeff * aNewLength;
332 theResPoint.Z = aResVal;
335 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
338 int nbPts = triangle->GetNumberOfIds();
341 //DEBTRACE("not a triangle ?");
345 double v[3][3]; // v[i][j] = j coordinate of node i
346 for (int i=0; i<3; i++)
348 s[i] = triangle->GetId(i);
349 delaunay2D->GetPoint(s[i],v[i]);
351 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
352 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
353 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
354 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
356 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
357 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
358 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]);
361 //DEBTRACE("flat triangle ?");
365 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
366 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
369 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
370 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
373 double l2 = 1 -l0 -l1;
374 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
376 if ((l0>=0) && (l0<=1) && (l1>=0) && (l1<=1) && (l2>=0) && (l2<=1))
378 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
385 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
387 //DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
388 double anInvalidAltitude = GetInvalidAltitude();
389 double aResAltitude = anInvalidAltitude;
391 // --- find the nearest point in the bathymetry cloud, with quadtree
393 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
396 //DEBTRACE(" no Quadtree");
400 std::map<double, const gpi_XYZ*> dist2nodes;
401 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
402 while (dist2nodes.size() == 0)
404 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
405 //DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
406 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
408 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
409 aResAltitude = it->second->Z();
410 int nodeIndex = it->second->getIndex();
411 //DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
413 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
414 // interpolation is required.
415 // - get a Delaunay2D mesh on the bathymetry cloud,
416 // - get the triangle containing the point in the Delaunay2D mesh,
417 // - interpolate altitude
419 bool isBathyInterpolRequired = false;
421 isBathyInterpolRequired =true;
424 if (isBathyInterpolRequired)
426 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
427 vtkIdList* cells= vtkIdList::New();
429 vtkIdList* points= vtkIdList::New();
430 points->Allocate(64);
431 aDelaunay2D->GetPointCells(nodeIndex, cells);
432 vtkIdType nbCells = cells->GetNumberOfIds();
433 //DEBTRACE(" triangles on nearest point: " << nbCells);
434 bool isInside = false;
435 for (int i=0; i<nbCells; i++)
437 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
439 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
443 //DEBTRACE(" interpolated z: " << z);
449 // DEBTRACE(" point outside triangles, nearest z kept");
456 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
458 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
461 void HYDROData_Bathymetry::SetFilePaths( const QStringList& theFilePaths )
464 Handle_TDataStd_ExtStringArray TExtStrArr = TDataStd_ExtStringArray::Set( myLab.FindChild( DataTag_FilePaths ), 1, theFilePaths.size() );
465 foreach (QString filepath, theFilePaths)
467 std::string sstr = filepath.toStdString();
468 const char* Val = sstr.c_str();
469 TExtStrArr->SetValue(i, TCollection_ExtendedString(Val));
474 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
476 TCollection_AsciiString aRes;
478 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
479 if ( !aLabel.IsNull() )
481 Handle(TDataStd_AsciiString) anAsciiStr;
482 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
483 aRes = anAsciiStr->Get();
487 aLabel = myLab.FindChild( DataTag_FilePaths, false );
488 if ( !aLabel.IsNull() )
490 Handle(TDataStd_ExtStringArray) anExtStrArr;
491 if ( aLabel.FindAttribute( TDataStd_ExtStringArray::GetID(), anExtStrArr ) )
492 aRes = anExtStrArr->Value(1); //try take the first; convert extstring to asciistring
499 QStringList HYDROData_Bathymetry::GetFilePaths() const
503 TDF_Label aLabel = myLab.FindChild( DataTag_FilePaths, false );
504 if ( !aLabel.IsNull() )
506 Handle(TDataStd_ExtStringArray) anExtStrArr;
507 if ( aLabel.FindAttribute( TDataStd_ExtStringArray::GetID(), anExtStrArr ) )
509 for (int i = anExtStrArr->Lower(); i <= anExtStrArr->Upper(); i++ )
511 Standard_ExtString str = anExtStrArr->Value(i).ToExtString();
512 TCollection_AsciiString aText (str);
513 aResL << QString(aText.ToCString());
517 else //backward compatibility
519 TDF_Label anOldLabel = myLab.FindChild( DataTag_FilePath, false );
520 if ( !anOldLabel.IsNull() )
522 Handle(TDataStd_AsciiString) anAsciiStr;
523 if ( anOldLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
524 aResL << QString(anAsciiStr->Get().ToCString());
531 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
532 const bool theIsUpdate )
534 bool anIsAltitudesInverted = IsAltitudesInverted();
535 if ( anIsAltitudesInverted == theIsInverted )
538 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
545 // Update altitude points
546 HYDROData_Bathymetry::AltitudePoints anAltitudePoints = GetAltitudePoints();
547 if ( anAltitudePoints.empty() )
550 HYDROData_Bathymetry::AltitudePoints::iterator anIter = anAltitudePoints.begin(), aLast = anAltitudePoints.end();
551 for ( ; anIter!=aLast; ++anIter )
553 HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
557 SetAltitudePoints( anAltitudePoints );
560 bool HYDROData_Bathymetry::IsAltitudesInverted() const
564 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
565 if ( !aLabel.IsNull() )
567 Handle(TDataStd_Integer) anIntVal;
568 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
569 aRes = (bool)anIntVal->Get();
575 bool HYDROData_Bathymetry::ImportFromFile( const QString& theFileName )
577 return ImportFromFiles(QStringList(theFileName));
580 bool HYDROData_Bathymetry::ImportFromFiles( const QStringList& theFileNames )
582 AltitudePoints AllPoints;
585 foreach (QString theFileName, theFileNames)
587 // Try to open the file
588 QFile aFile( theFileName );
589 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
592 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
594 HYDROData_Bathymetry::AltitudePoints aPoints;
596 // Try to import the file
597 if ( aFileSuf == "xyz" )
598 Stat = importFromXYZFile( aFile, aPoints );
599 else if ( aFileSuf == "asc" )
600 Stat = importFromASCFile( aFile, aPoints );
603 continue; //ignore this points
608 AllPoints.insert(AllPoints.end(), aPoints.begin(), aPoints.end());
611 // Convert from global to local CS
612 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
613 HYDROData_Bathymetry::AltitudePoints::iterator anIter = AllPoints.begin(), aLast = AllPoints.end();
614 for ( ; anIter!=aLast; ++anIter )
616 HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
617 aDoc->Transform( aPoint.X, aPoint.Y, aPoint.Z, true );
622 // Update file path and altitude points of this Bathymetry
623 SetFilePaths (theFileNames );
624 SetAltitudePoints( AllPoints );
627 return Stat && !AllPoints.empty();
630 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
631 HYDROData_Bathymetry::AltitudePoints& thePoints ) const
633 if ( !theFile.isOpen() )
636 // Strings in file is written as:
637 // 1. X(float) Y(float) Z(float)
638 // 2. X(float) Y(float) Z(float)
646 bool anIsAltitudesInverted = IsAltitudesInverted();
647 while ( !theFile.atEnd() )
649 QString aLine = theFile.readLine().simplified();
650 if ( aLine.isEmpty() )
653 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
654 if ( aValues.length() < 3 )
657 HYDROData_Bathymetry::AltitudePoint aPoint;
659 QString anX = aValues.value( 0 );
660 QString anY = aValues.value( 1 );
661 QString aZ = aValues.value( 2 );
663 bool isXOk = false, isYOk = false, isZOk = false;
665 aPoint.X = anX.toDouble( &isXOk );
666 aPoint.Y = anY.toDouble( &isYOk );
667 aPoint.Z = aZ.toDouble( &isZOk );
669 if ( !isXOk || !isYOk || !isZOk )
672 if ( HYDROData_Tool::IsNan( aPoint.X ) || HYDROData_Tool::IsInf( aPoint.X ) ||
673 HYDROData_Tool::IsNan( aPoint.Y ) || HYDROData_Tool::IsInf( aPoint.Y ) ||
674 HYDROData_Tool::IsNan( aPoint.Z ) || HYDROData_Tool::IsInf( aPoint.Z ) )
677 // Invert the z value if requested
678 if ( anIsAltitudesInverted )
679 aPoint.Z = -aPoint.Z;
681 if( thePoints.size()>=thePoints.capacity() )
682 thePoints.reserve( thePoints.size()+BLOCK_SIZE );
684 thePoints.push_back( aPoint );
689 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
690 aTimer.Show( stream );
696 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
697 HYDROData_Bathymetry::AltitudePoints& thePoints ) const
699 if ( !theFile.isOpen() )
703 QStringList aStrList;
712 aLine = theFile.readLine().simplified();
713 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
714 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
716 aNCols = aStrList[1].toInt();
718 aLine = theFile.readLine().simplified();
719 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
720 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
722 aNRows = aStrList[1].toInt();
724 aLine = theFile.readLine().simplified();
725 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
726 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
728 anXllCorner = aStrList[1].toDouble();
730 aLine = theFile.readLine().simplified();
731 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
732 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
734 anYllCorner = aStrList[1].toDouble();
736 aLine = theFile.readLine().simplified();
737 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
738 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
740 aCellSize = aStrList[1].toDouble();
742 aLine = theFile.readLine().simplified();
743 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
744 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
746 aNoDataValue = aStrList[1].toDouble();
748 bool anIsAltitudesInverted = IsAltitudesInverted();
752 while ( !theFile.atEnd() )
754 aLine = theFile.readLine().simplified();
755 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
757 aStrLength = aStrList.length();
758 if ( aStrLength == 0 )
761 if ( aStrLength != aNRows )
764 for (int j = 0; j < aNCols; j++)
766 if (aStrList[j].toDouble() != aNoDataValue)
768 HYDROData_Bathymetry::AltitudePoint aPoint;
769 aPoint.X = anXllCorner + aCellSize*(j + 0.5);
770 aPoint.Y = anYllCorner + aCellSize*(aNRows - i + 0.5);
771 aPoint.Z = aStrList[j].toDouble();
773 if ( anIsAltitudesInverted )
774 aPoint.Z = -aPoint.Z;
776 if( thePoints.size()>=thePoints.capacity() )
777 thePoints.reserve( thePoints.size()+BLOCK_SIZE );
778 thePoints.push_back(aPoint);
788 Handle(HYDROData_PolylineXY) HYDROData_Bathymetry::CreateBoundaryPolyline() const
790 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
791 Handle(HYDROData_PolylineXY) aResult =
792 Handle(HYDROData_PolylineXY)::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
794 if( aResult.IsNull() )
798 QString aPolylinePref = GetName() + "_Boundary";
799 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
800 aResult->SetName( aPolylineName );
802 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
804 HYDROData_Bathymetry::AltitudePoints aPoints = GetAltitudePoints();
806 HYDROData_Bathymetry::AltitudePoints::const_iterator anIter = aPoints.begin(), aLast = aPoints.end();
807 for ( ; anIter!=aLast; ++anIter )
809 const HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
811 double x = aPoint.X, y = aPoint.Y;
812 if( isFirst || x<Xmin )
814 if( isFirst || x>Xmax )
816 if( isFirst || y<Ymin )
818 if( isFirst || y>Ymax )
823 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
824 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
825 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
826 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
827 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
829 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
836 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
838 gp_XYZ aDelta( theDx, theDy, 0 );
839 HYDROData_Bathymetry::AltitudePoints aPoints = GetAltitudePoints();
840 HYDROData_Bathymetry::AltitudePoints::iterator anIter = aPoints.begin(), aLast = aPoints.end();
841 for ( int i = 0; anIter!=aLast; ++i, ++anIter )
843 HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
844 aPoint.X += aDelta.X();
845 aPoint.Y += aDelta.Y();
846 aPoint.Z += aDelta.Z();
848 SetAltitudePoints( aPoints );