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>
43 #include <vtkInformation.h>
44 #include <vtkExecutive.h>
52 #include <OSD_Timer.hxx>
56 #include "HYDRO_trace.hxx"
58 IMPLEMENT_STANDARD_HANDLE(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
59 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
61 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
62 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
63 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
66 hydroDelaunay2D* hydroDelaunay2D::New()
68 DEBTRACE("hydroDelaunay2D::New()");
69 return new hydroDelaunay2D();
71 hydroDelaunay2D::hydroDelaunay2D() :
74 DEBTRACE("hydroDelaunay2D");
77 hydroDelaunay2D::~hydroDelaunay2D()
79 DEBTRACE("~hydroDelaunay2D");
82 int hydroDelaunay2D::RequestData(vtkInformation *request , vtkInformationVector **inv, vtkInformationVector *outv)
84 DEBTRACE("hydroDelaunay2D::RequestData");
85 vtkDelaunay2D::RequestData(request, inv, outv);
88 int hydroDelaunay2D::ProcessRequest(vtkInformation *request , vtkInformationVector **inv, vtkInformationVector *outv)
90 DEBTRACE("hydroDelaunay2D::ProcessRequest");
91 vtkDelaunay2D::ProcessRequest(request, inv, outv);
95 HYDROData_Bathymetry::HYDROData_Bathymetry()
96 : HYDROData_IAltitudeObject()
98 //DEBTRACE("HYDROData_Bathymetry constructor start " << this);
100 // myQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
101 //DEBTRACE("HYDROData_Bathymetry constructor end " << this);
104 HYDROData_Bathymetry::~HYDROData_Bathymetry()
106 //DEBTRACE("HYDROData_Bathymetry destructor start " << this);
108 // delete myQuadtree;
109 // Nodes_3D::iterator it = myListOfNodes.begin();
110 // for( ; it != myListOfNodes.end(); ++it)
112 // myListOfNodes.clear();
115 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
116 MapOfTreatedObjects& theTreatedObjects ) const
118 QStringList aResList = dumpObjectCreation( theTreatedObjects );
119 QString aBathymetryName = GetObjPyName();
121 aResList << QString( "%1.SetAltitudesInverted( %2 )" )
122 .arg( aBathymetryName ).arg( IsAltitudesInverted() );
124 TCollection_AsciiString aFilePath = GetFilePath();
125 aResList << QString( "if not(%1.ImportFromFile( \"%2\" )):" )
126 .arg( aBathymetryName ).arg( aFilePath.ToCString() );
127 aResList << QString( " raise ValueError('problem while loading bathymetry')" );
128 aResList << QString( "" );
129 aResList << QString( "%1.Update()" ).arg( aBathymetryName );
130 aResList << QString( "" );
135 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
137 RemoveAltitudePoints();
139 if ( thePoints.IsEmpty() )
143 Handle(TDataStd_RealArray) aCoordsArray =
144 TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
146 AltitudePoints::Iterator anIter( thePoints );
147 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
149 const AltitudePoint& aPoint = anIter.Value();
151 aCoordsArray->SetValue( i * 3, aPoint.X() );
152 aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
153 aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
159 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
161 AltitudePoints aPoints;
163 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
164 if ( aLabel.IsNull() )
167 Handle(TDataStd_RealArray) aCoordsArray;
168 if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
171 Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
172 for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
177 AltitudePoint aPoint;
178 aPoint.SetX( aCoordsArray->Value( i++ ) );
179 aPoint.SetY( aCoordsArray->Value( i++ ) );
180 aPoint.SetZ( aCoordsArray->Value( i++ ) );
182 if( IsConvertToGlobal )
183 aDoc->Transform( aPoint, false );
184 aPoints.Append( aPoint );
190 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
192 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
195 int labkey = myLab.Tag();
196 //int altkey = aLabel.Tag();
197 //DEBTRACE("GetQuadtreeNodes this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
198 // if (myQuadtree->isEmpty() )
199 if (myQuadtrees.find(labkey) == myQuadtrees.end())
201 DEBTRACE("GetQuadtreeNodes init " << this << " " << labkey);
202 HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
203 myQuadtrees[labkey] = aQuadtree;
204 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
208 Handle(TDataStd_RealArray) aCoordsArray;
209 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
212 Nodes_3D* aListOfNodes = new Nodes_3D();
215 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
220 double x = aCoordsArray->Value(i++);
221 double y = aCoordsArray->Value(i++);
222 double z = aCoordsArray->Value(i++);
223 gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
225 aListOfNodes->push_back(aPoint);
227 DEBTRACE(" GetQuadtreeNodes call setNodesAndCompute");
228 aQuadtree->setNodesAndCompute(aListOfNodes);
232 return myQuadtrees[labkey];
235 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
237 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
240 int labkey = myLab.Tag();
241 //int altkey = aLabel.Tag();
242 //DEBTRACE("GetVtkDelaunay2D this labkey altkey "<<this<<" "<<labkey<<" "<<altkey);
243 if (myDelaunay2D.find(labkey) == myDelaunay2D.end())
245 DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
247 TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
251 Handle(TDataStd_RealArray) aCoordsArray;
252 if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
256 vtkPoints *points = vtkPoints::New();
257 points->Allocate(aCoordsArray->Upper() +1);
258 for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
262 double x = aCoordsArray->Value(i++);
263 double y = aCoordsArray->Value(i++);
264 double z = aCoordsArray->Value(i++);
265 vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
266 //DEBTRACE(" " << index);
268 vtkPolyData* profile = vtkPolyData::New();
269 profile->SetPoints(points);
270 DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
272 hydroDelaunay2D* delaunay2D = hydroDelaunay2D::New();
273 delaunay2D->GetInputPortInformation(0)->Print(std::cerr);
274 DEBTRACE("set the input data");
275 delaunay2D->SetInputData(profile);
276 delaunay2D->GetInputPortInformation(0)->Print(std::cerr);
277 delaunay2D->GetOutputPortInformation(0)->Print(std::cerr);
279 delaunay2D->GetOutputInformation(0)->Print(std::cerr);
280 delaunay2D->GetExecutive()->GetOutputInformation(0)->Print(std::cerr);
282 vtkInformationVector** inv = delaunay2D->GetExecutive()->GetInputInformation();
283 vtkInformationVector* outv = delaunay2D->GetExecutive()->GetOutputInformation();
284 delaunay2D->GetExecutive()->GetInputInformation(0,0)->Print(std::cerr);
285 delaunay2D->GetExecutive()->GetOutputInformation(0)->Print(std::cerr);
286 vtkInformation* request = vtkInformation::New();
288 //delaunay2D->ProcessRequest(request, inv,outv);
290 delaunay2D->GetOutputPortInformation(0)->Print(std::cerr);
291 delaunay2D->GetExecutive()->GetInputInformation(0,0)->Print(std::cerr);
292 delaunay2D->GetExecutive()->GetOutputInformation(0)->Print(std::cerr);
295 delaunay2D->Update();
298 delaunay2D->UpdateDataObject();
300 vtkPolyData* data = delaunay2D->GetOutput();
303 myDelaunay2D[labkey] = data;
308 return myDelaunay2D[labkey];
312 void HYDROData_Bathymetry::RemoveAltitudePoints()
314 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
315 if ( !aLabel.IsNull() )
317 aLabel.ForgetAllAttributes();
322 void interpolateAltitudeForPoints( const gp_XY& thePoint,
323 const HYDROData_Bathymetry::AltitudePoint& theFirstPoint,
324 const HYDROData_Bathymetry::AltitudePoint& theSecPoint,
325 HYDROData_Bathymetry::AltitudePoint& theResPoint,
326 const bool& theIsVertical )
328 double aCoordX = thePoint.X();
329 double aCoordY = thePoint.Y();
333 aCoordX = theFirstPoint.X();
335 if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
337 // Recalculate X coordinate by equation of line from two points
338 aCoordX = ( ( ( thePoint.Y() - theFirstPoint.Y() ) * ( theSecPoint.X() - theFirstPoint.X() ) ) /
339 ( theSecPoint.Y() - theFirstPoint.Y() ) ) + theFirstPoint.X();
344 aCoordY = theFirstPoint.Y();
346 if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
348 // Recalculate y by equation of line from two points
349 aCoordY = ( ( ( thePoint.X() - theFirstPoint.X() ) * ( theSecPoint.Y() - theFirstPoint.Y() ) ) /
350 ( theSecPoint.X() - theFirstPoint.X() ) ) + theFirstPoint.Y();
354 theResPoint.SetX( aCoordX );
355 theResPoint.SetY( aCoordY );
357 // Calculate coefficient for interpolation
358 double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
359 Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
361 double aInterCoeff = 0;
363 aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
366 double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
367 Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
369 // Calculate interpolated value
370 double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
372 theResPoint.SetZ( aResVal );
375 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
378 int nbPts = triangle->GetNumberOfIds();
381 DEBTRACE("not a triangle ?");
385 double v[3][3]; // v[i][j] = j coordinate of node i
386 for (int i=0; i<3; i++)
388 s[i] = triangle->GetId(i);
389 delaunay2D->GetPoint(s[i],v[i]);
391 //DEBTRACE("triangle node id: " << s[0] << " " << s[1] << " " << s[2]);
392 //DEBTRACE("triangle node 0: " << v[0][0] << " " << v[0][1] << " " << v[0][2]);
393 //DEBTRACE("triangle node 1: " << v[1][0] << " " << v[1][1] << " " << v[1][2]);
394 //DEBTRACE("triangle node 2: " << v[2][0] << " " << v[2][1] << " " << v[2][2]);
396 // compute barycentric coordinates (https://en.wikipedia.org/wiki/Barycentric_coordinate_system)
397 // det = (y2-y3)(x1-x3)+(x3-x2)(y1-y3)
398 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]);
401 DEBTRACE("flat triangle ?");
405 // l0 = ((y2-y3)(x -x3)+(x3-x2)(y -y3))/det
406 double l0 = (v[1][1]-v[2][1])*(point.X()-v[2][0]) + (v[2][0]-v[1][0])*(point.Y()-v[2][1]);
409 // l1 = ((y3-y1)(x -x3)+(x1-x3)(y -y3))/det
410 double l1 = (v[2][1]-v[0][1])*(point.X()-v[2][0]) + (v[0][0]-v[2][0])*(point.Y()-v[2][1]);
413 double l2 = 1 -l0 -l1;
414 //DEBTRACE("l0, l1, l2: " << l0 << " " << l1 << " " << l2);
416 if ((l0>=0) and (l0<=1) and (l1>=0) and (l1<=1) and (l2>=0) and (l2<=1))
418 z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
424 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint) const
426 DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << ")");
427 double anInvalidAltitude = GetInvalidAltitude();
428 double aResAltitude = anInvalidAltitude;
430 // --- find the nearest point in the bathymetry cloud, with quadtree
432 HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
435 DEBTRACE(" no Quadtree");
439 std::map<double, const gpi_XYZ*> dist2nodes;
440 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
441 while (dist2nodes.size() == 0)
443 aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
444 DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
445 aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
447 std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
448 aResAltitude = it->second->Z();
449 int nodeIndex = it->second->getIndex();
450 DEBTRACE(" number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
452 // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
453 // interpolation is required.
454 // - get a Delaunay2D mesh on the bathymetry cloud,
455 // - get the triangle containing the point in the Delaunay2D mesh,
456 // - interpolate altitude
458 bool isBathyInterpolRequired = true;
459 if (isBathyInterpolRequired)
461 vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
462 vtkIdList* cells= vtkIdList::New();
464 vtkIdList* points= vtkIdList::New();
465 points->Allocate(64);
466 aDelaunay2D->GetPointCells(nodeIndex, cells);
467 vtkIdType nbCells = cells->GetNumberOfIds();
468 DEBTRACE(" triangles on nearest point: " << nbCells);
469 bool isInside = false;
470 for (int i=0; i<nbCells; i++)
472 aDelaunay2D->GetCellPoints(cells->GetId(i), points);
474 isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
478 DEBTRACE(" interpolated z: " << z);
482 if (!isInside) DEBTRACE(" point outside triangles, nearest z kept");
487 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
489 TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
492 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
494 TCollection_AsciiString aRes;
496 TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
497 if ( !aLabel.IsNull() )
499 Handle(TDataStd_AsciiString) anAsciiStr;
500 if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
501 aRes = anAsciiStr->Get();
507 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
508 const bool theIsUpdate )
510 bool anIsAltitudesInverted = IsAltitudesInverted();
511 if ( anIsAltitudesInverted == theIsInverted )
514 TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
521 // Update altitude points
522 AltitudePoints anAltitudePoints = GetAltitudePoints();
523 if ( anAltitudePoints.IsEmpty() )
526 AltitudePoints::Iterator anIter( anAltitudePoints );
527 for ( ; anIter.More(); anIter.Next() )
529 AltitudePoint& aPoint = anIter.ChangeValue();
530 aPoint.SetZ( aPoint.Z() * -1 );
533 SetAltitudePoints( anAltitudePoints );
536 bool HYDROData_Bathymetry::IsAltitudesInverted() const
540 TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
541 if ( !aLabel.IsNull() )
543 Handle(TDataStd_Integer) anIntVal;
544 if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
545 aRes = (bool)anIntVal->Get();
551 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
553 // Try to open the file
554 QFile aFile( theFileName.ToCString() );
555 if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
560 QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
562 AltitudePoints aPoints;
564 // Try to import the file
565 if ( aFileSuf == "xyz" )
566 aRes = importFromXYZFile( aFile, aPoints );
567 else if ( aFileSuf == "asc" )
568 aRes = importFromASCFile( aFile, aPoints );
574 // Convert from global to local CS
575 Handle_HYDROData_Document aDoc = HYDROData_Document::Document( myLab );
576 AltitudePoints::Iterator anIter( aPoints );
577 for ( ; anIter.More(); anIter.Next() )
579 AltitudePoint& aPoint = anIter.ChangeValue();
580 aDoc->Transform( aPoint, true );
585 // Update file path and altitude points of this Bathymetry
586 SetFilePath( theFileName );
587 SetAltitudePoints( aPoints );
590 return aRes && !aPoints.IsEmpty();
593 bool HYDROData_Bathymetry::importFromXYZFile( QFile& theFile,
594 AltitudePoints& thePoints ) const
596 if ( !theFile.isOpen() )
599 // Strings in file is written as:
600 // 1. X(float) Y(float) Z(float)
601 // 2. X(float) Y(float) Z(float)
609 bool anIsAltitudesInverted = IsAltitudesInverted();
610 while ( !theFile.atEnd() )
612 QString aLine = theFile.readLine().simplified();
613 if ( aLine.isEmpty() )
616 QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
617 if ( aValues.length() < 3 )
620 AltitudePoint aPoint;
622 QString anX = aValues.value( 0 );
623 QString anY = aValues.value( 1 );
624 QString aZ = aValues.value( 2 );
626 bool isXOk = false, isYOk = false, isZOk = false;
628 aPoint.SetX( anX.toDouble( &isXOk ) );
629 aPoint.SetY( anY.toDouble( &isYOk ) );
630 aPoint.SetZ( aZ.toDouble( &isZOk ) );
632 if ( !isXOk || !isYOk || !isZOk )
635 if ( HYDROData_Tool::IsNan( aPoint.X() ) || HYDROData_Tool::IsInf( aPoint.X() ) ||
636 HYDROData_Tool::IsNan( aPoint.Y() ) || HYDROData_Tool::IsInf( aPoint.Y() ) ||
637 HYDROData_Tool::IsNan( aPoint.Z() ) || HYDROData_Tool::IsInf( aPoint.Z() ) )
640 // Invert the z value if requested
641 if ( anIsAltitudesInverted )
642 aPoint.SetZ( -aPoint.Z() );
644 thePoints.Append( aPoint );
649 std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
650 aTimer.Show( stream );
656 bool HYDROData_Bathymetry::importFromASCFile( QFile& theFile,
657 AltitudePoints& thePoints ) const
659 if ( !theFile.isOpen() )
663 QStringList aStrList;
672 aLine = theFile.readLine().simplified();
673 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
674 if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
676 aNCols = aStrList[1].toInt();
678 aLine = theFile.readLine().simplified();
679 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
680 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
682 aNRows = aStrList[1].toInt();
684 aLine = theFile.readLine().simplified();
685 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
686 if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
688 anXllCorner = aStrList[1].toDouble();
690 aLine = theFile.readLine().simplified();
691 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
692 if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
694 anYllCorner = aStrList[1].toDouble();
696 aLine = theFile.readLine().simplified();
697 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
698 if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
700 aCellSize = aStrList[1].toDouble();
702 aLine = theFile.readLine().simplified();
703 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
704 if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
706 aNoDataValue = aStrList[1].toDouble();
708 bool anIsAltitudesInverted = IsAltitudesInverted();
712 while ( !theFile.atEnd() )
714 aLine = theFile.readLine().simplified();
715 aStrList = aLine.split( ' ', QString::SkipEmptyParts );
717 aStrLength = aStrList.length();
718 if ( aStrLength == 0 )
721 if ( aStrLength != aNRows )
724 for (int j = 0; j < aNCols; j++)
726 if (aStrList[j].toDouble() != aNoDataValue)
728 AltitudePoint aPoint;
729 aPoint.SetX(anXllCorner + aCellSize*(j + 0.5));
730 aPoint.SetY(anYllCorner + aCellSize*(aNRows - i + 0.5));
731 aPoint.SetZ(aStrList[j].toDouble());
733 if ( anIsAltitudesInverted )
734 aPoint.SetZ( -aPoint.Z() );
736 thePoints.Append(aPoint);
748 Handle_HYDROData_PolylineXY HYDROData_Bathymetry::CreateBoundaryPolyline() const
750 Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
751 Handle_HYDROData_PolylineXY aResult =
752 Handle_HYDROData_PolylineXY::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
754 if( aResult.IsNull() )
758 QString aPolylinePref = GetName() + "_Boundary";
759 QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
760 aResult->SetName( aPolylineName );
762 double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
764 AltitudePoints aPoints = GetAltitudePoints();
766 AltitudePoints::Iterator anIter( aPoints );
767 for ( ; anIter.More(); anIter.Next() )
769 const AltitudePoint& aPoint = anIter.Value();
771 double x = aPoint.X(), y = aPoint.Y();
772 if( isFirst || x<Xmin )
774 if( isFirst || x>Xmax )
776 if( isFirst || y<Ymin )
778 if( isFirst || y>Ymax )
783 aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
784 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
785 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
786 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
787 aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
789 aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
796 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
798 gp_XYZ aDelta( theDx, theDy, 0 );
799 AltitudePoints aPoints = GetAltitudePoints();
800 AltitudePoints::Iterator anIter( aPoints );
801 for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
803 AltitudePoint& aPoint = anIter.ChangeValue();
806 SetAltitudePoints( aPoints );