]> SALOME platform Git repositories - modules/hydro.git/blob - src/HYDROData/HYDROData_Bathymetry.cxx
Salome HOME
first part of the porting on OCCT 7.0
[modules/hydro.git] / src / HYDROData / HYDROData_Bathymetry.cxx
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.
6 //
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.
11 //
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
15 //
16 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
17 //
18
19 #include "HYDROData_Bathymetry.h"
20 #include "HYDROData_Document.h"
21 #include "HYDROData_Tool.h"
22 #include "HYDROData_PolylineXY.h"
23
24 #include <gp_XY.hxx>
25 #include <gp_XYZ.hxx>
26
27 #include <TDataStd_RealArray.hxx>
28 #include <TDataStd_AsciiString.hxx>
29 #include <TDataStd_Integer.hxx>
30
31 #include <QColor>
32 #include <QFile>
33 #include <QFileInfo>
34 #include <QPointF>
35 #include <QPolygonF>
36 #include <QStringList>
37
38 #include <vtkPoints.h>
39 #include <vtkDelaunay2D.h>
40 #include <vtkPolyData.h>
41 #include <vtkSmartPointer.h>
42 #include <vtkIdList.h>
43
44 #include <iostream>
45
46 #include <math.h>
47
48 // #define _TIMER
49 #ifdef _TIMER
50 #include <OSD_Timer.hxx>
51 #endif
52
53 #define _DEVDEBUG_
54 #include "HYDRO_trace.hxx"
55
56 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
57
58 //HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
59 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
60 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
61
62
63 HYDROData_Bathymetry::HYDROData_Bathymetry()
64 : HYDROData_IAltitudeObject()
65 {
66 }
67
68 HYDROData_Bathymetry::~HYDROData_Bathymetry()
69 {
70 }
71
72 QStringList HYDROData_Bathymetry::DumpToPython( const QString& thePyScriptPath,
73                                                 MapOfTreatedObjects& theTreatedObjects ) const
74 {
75   QStringList aResList = dumpObjectCreation( theTreatedObjects );
76   QString aBathymetryName = GetObjPyName();
77
78   aResList << QString( "%1.SetAltitudesInverted( %2 )" )
79               .arg( aBathymetryName ).arg( IsAltitudesInverted() );
80
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( "" );
88
89   return aResList;
90 }
91
92 void HYDROData_Bathymetry::SetAltitudePoints( const AltitudePoints& thePoints )
93 {
94   RemoveAltitudePoints();
95
96   if ( thePoints.IsEmpty() )
97     return;
98
99   // Save coordinates
100   Handle(TDataStd_RealArray) aCoordsArray = 
101     TDataStd_RealArray::Set( myLab.FindChild( DataTag_AltitudePoints ), 0, thePoints.Length() * 3 - 1 );
102
103   AltitudePoints::Iterator anIter( thePoints );
104   for ( int i = 0 ; anIter.More(); ++i, anIter.Next() )
105   {
106     const AltitudePoint& aPoint = anIter.Value();
107
108     aCoordsArray->SetValue( i * 3, aPoint.X() );
109     aCoordsArray->SetValue( i * 3 + 1, aPoint.Y() );
110     aCoordsArray->SetValue( i * 3 + 2, aPoint.Z() );
111   }
112
113   Changed( Geom_Z );
114 }
115
116 HYDROData_Bathymetry::AltitudePoints HYDROData_Bathymetry::GetAltitudePoints(bool IsConvertToGlobal) const
117 {
118   AltitudePoints aPoints;
119
120   TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
121   if ( aLabel.IsNull() )
122     return aPoints;
123
124   Handle(TDataStd_RealArray) aCoordsArray;
125   if ( !aLabel.FindAttribute( TDataStd_RealArray::GetID(), aCoordsArray ) )
126     return aPoints;
127
128   Handle(HYDROData_Document) aDoc = HYDROData_Document::Document( myLab );
129   for ( int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n; )
130   {
131     if ( i + 3 > n + 1 )
132       break;
133
134     AltitudePoint aPoint;
135     aPoint.SetX( aCoordsArray->Value( i++ ) );
136     aPoint.SetY( aCoordsArray->Value( i++ ) );
137     aPoint.SetZ( aCoordsArray->Value( i++ ) );
138
139     if( IsConvertToGlobal )
140       aDoc->Transform( aPoint, false );
141     aPoints.Append( aPoint );
142   }
143
144   return aPoints;
145 }
146
147 HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
148 {
149   TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
150   if (aLabel.IsNull())
151     return 0;
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())
157     {
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);
162       if (aLabel.IsNull())
163         return 0;
164
165       Handle(TDataStd_RealArray) aCoordsArray;
166       if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
167         return 0;
168
169       Nodes_3D* aListOfNodes = new Nodes_3D();
170
171       int index =0;
172       for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
173         {
174           if (i + 3 > n + 1)
175             break;
176
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);
181           index++;
182           aListOfNodes->push_back(aPoint);
183         }
184       DEBTRACE("  GetQuadtreeNodes call setNodesAndCompute");
185       aQuadtree->setNodesAndCompute(aListOfNodes);
186       return aQuadtree;
187     }
188   else
189     return myQuadtrees[labkey];
190 }
191
192 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
193 {
194   TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
195   if (aLabel.IsNull())
196     return 0;
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())
201     {
202       DEBTRACE("GetVtkDelaunay2D init " << this << " " << labkey);
203
204       TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
205       if (aLabel.IsNull())
206         return 0;
207       Handle(TDataStd_RealArray) aCoordsArray;
208       if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
209         return 0;
210
211       vtkPoints *points = vtkPoints::New();
212       points->Allocate(aCoordsArray->Upper() +1);
213       for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
214         {
215           if (i + 3 > n + 1)
216             break;
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);
222         }
223       vtkPolyData* profile = vtkPolyData::New();
224       profile->SetPoints(points);
225       DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
226
227       vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
228       delaunay2D->SetInputData(profile);
229       delaunay2D->Update();
230       vtkPolyData* data = delaunay2D->GetOutput();
231       data->BuildLinks();
232       myDelaunay2D[labkey] = data;
233       return data;
234     }
235   else
236     return myDelaunay2D[labkey];
237
238 }
239
240 void HYDROData_Bathymetry::RemoveAltitudePoints()
241 {
242   TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
243   if ( !aLabel.IsNull() )
244   {
245     aLabel.ForgetAllAttributes();
246     Changed( Geom_Z );
247   }
248 }
249
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 )
255 {
256   double aCoordX = thePoint.X();
257   double aCoordY = thePoint.Y();
258
259   if ( theIsVertical )
260   {
261     aCoordX = theFirstPoint.X();
262
263     if ( !ValuesEquals( theFirstPoint.X(), theSecPoint.X() ) )
264     {
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();
268     }
269   }
270   else
271   {
272     aCoordY = theFirstPoint.Y();
273
274     if ( !ValuesEquals( theFirstPoint.Y(), theSecPoint.Y() ) )
275     {
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();
279     }
280   }
281
282   theResPoint.SetX( aCoordX );
283   theResPoint.SetY( aCoordY );
284
285   // Calculate coefficient for interpolation
286   double aLength = Sqrt( Pow( theSecPoint.Y() - theFirstPoint.Y(), 2 ) +
287                          Pow( theSecPoint.X() - theFirstPoint.X(), 2 ) );
288
289   double aInterCoeff = 0;
290   if ( aLength != 0 )
291    aInterCoeff = ( theSecPoint.Z() - theFirstPoint.Z() ) / aLength;
292
293
294   double aNewLength = Sqrt( Pow( theResPoint.Y() - theFirstPoint.Y(), 2 ) +
295                             Pow( theResPoint.X() - theFirstPoint.X(), 2 ) );
296
297   // Calculate interpolated value
298   double aResVal = theFirstPoint.Z() + aInterCoeff * aNewLength;
299
300   theResPoint.SetZ( aResVal );
301 }
302
303 bool interpolZtriangle(const gp_XY& point, vtkPolyData* delaunay2D, vtkIdList* triangle, double& z)
304 {
305
306   int nbPts = triangle->GetNumberOfIds();
307   if (nbPts != 3)
308     {
309       DEBTRACE("not a triangle ?");
310       return false;
311     }
312   vtkIdType s[3];
313   double v[3][3]; // v[i][j] = j coordinate of node i
314   for (int i=0; i<3; i++)
315     {
316       s[i] = triangle->GetId(i);
317       delaunay2D->GetPoint(s[i],v[i]);
318     }
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]);
323
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]);
327   if (det == 0)
328     {
329       DEBTRACE("flat triangle ?");
330       return false;
331     }
332
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]);
335   l0 = l0/det;
336
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]);
339   l1 = l1/det;
340
341   double l2  = 1 -l0 -l1;
342   //DEBTRACE("l0, l1, l2: " << l0  << " "  << l1  << " "  << l2);
343
344   if ((l0>=0) and (l0<=1) and (l1>=0) and (l1<=1) and (l2>=0) and (l2<=1))
345     {
346       z = l0*v[0][2] + l1*v[1][2] + l2*v[2][2];
347       return true;
348     }
349   return false;
350 }
351
352 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
353 {
354   DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
355   double anInvalidAltitude = GetInvalidAltitude();
356   double aResAltitude = anInvalidAltitude;
357
358   // --- find the nearest point in the bathymetry cloud, with quadtree
359
360   HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
361   if (!aQuadtree)
362     {
363       DEBTRACE("  no Quadtree");
364       return aResAltitude;
365     }
366
367   std::map<double, const gpi_XYZ*> dist2nodes;
368   aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
369   while (dist2nodes.size() == 0)
370     {
371       aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
372       DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
373       aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
374     }
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);
379
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
385
386   bool isBathyInterpolRequired = false;
387   if (theMethod)
388     isBathyInterpolRequired =true;
389   if (isBathyInterpolRequired)
390     {
391       vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
392       vtkIdList* cells= vtkIdList::New();
393       cells->Allocate(64);
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++)
401         {
402           aDelaunay2D->GetCellPoints(cells->GetId(i), points);
403           double z = 0;
404           isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
405           if (isInside)
406             {
407               aResAltitude = z;
408               DEBTRACE("  interpolated z: " << z);
409               break;
410             }
411         }
412       if (!isInside) DEBTRACE("  point outside triangles, nearest z kept");
413     }
414   return aResAltitude;
415 }
416
417 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
418 {
419   TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
420 }
421
422 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
423 {
424   TCollection_AsciiString aRes;
425
426   TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
427   if ( !aLabel.IsNull() )
428   {
429     Handle(TDataStd_AsciiString) anAsciiStr;
430     if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
431       aRes = anAsciiStr->Get();
432   }
433
434   return aRes;
435 }
436
437 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
438                                                  const bool theIsUpdate )
439 {
440   bool anIsAltitudesInverted = IsAltitudesInverted();
441   if ( anIsAltitudesInverted == theIsInverted )
442     return;
443
444   TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
445
446   Changed( Geom_Z );
447
448   if ( !theIsUpdate )
449     return;
450
451   // Update altitude points
452   AltitudePoints anAltitudePoints = GetAltitudePoints();
453   if ( anAltitudePoints.IsEmpty() )
454     return;
455
456   AltitudePoints::Iterator anIter( anAltitudePoints );
457   for ( ; anIter.More(); anIter.Next() )
458   {
459     AltitudePoint& aPoint = anIter.ChangeValue();
460     aPoint.SetZ( aPoint.Z() * -1 );
461   }
462
463   SetAltitudePoints( anAltitudePoints );
464 }
465
466 bool HYDROData_Bathymetry::IsAltitudesInverted() const
467 {
468   bool aRes = false;
469
470   TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
471   if ( !aLabel.IsNull() )
472   {
473     Handle(TDataStd_Integer) anIntVal;
474     if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
475       aRes = (bool)anIntVal->Get();
476   }
477
478   return aRes;
479 }
480
481 bool HYDROData_Bathymetry::ImportFromFile( const TCollection_AsciiString& theFileName )
482 {
483   // Try to open the file
484   QFile aFile( theFileName.ToCString() );
485   if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
486     return false;
487
488   bool aRes = false;
489
490   QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
491
492   AltitudePoints aPoints;
493
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 );
499
500   // Close the file
501   aFile.close();
502   
503
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() )
508   {
509     AltitudePoint& aPoint = anIter.ChangeValue();
510     aDoc->Transform( aPoint, true );
511   }
512
513   if ( aRes )
514   {
515     // Update file path and altitude points of this Bathymetry
516     SetFilePath( theFileName );
517     SetAltitudePoints( aPoints );
518   }
519
520   return aRes && !aPoints.IsEmpty();
521 }
522
523 bool HYDROData_Bathymetry::importFromXYZFile( QFile&          theFile,
524                                               AltitudePoints& thePoints ) const
525 {
526   if ( !theFile.isOpen() )
527     return false;
528
529   // Strings in file is written as:
530   //  1. X(float) Y(float) Z(float)
531   //  2. X(float) Y(float) Z(float)
532   //  ...
533
534 #ifdef _TIMER
535   OSD_Timer aTimer;
536   aTimer.Start();
537 #endif
538
539   bool anIsAltitudesInverted = IsAltitudesInverted();
540   while ( !theFile.atEnd() )
541   {
542     QString aLine = theFile.readLine().simplified();
543     if ( aLine.isEmpty() )
544       continue;
545
546     QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
547     if ( aValues.length() < 3 )
548       return false;
549
550     AltitudePoint aPoint;
551     
552     QString anX = aValues.value( 0 );
553     QString anY = aValues.value( 1 );
554     QString aZ  = aValues.value( 2 );
555
556     bool isXOk = false, isYOk = false, isZOk = false;
557
558     aPoint.SetX( anX.toDouble( &isXOk ) );
559     aPoint.SetY( anY.toDouble( &isYOk ) );
560     aPoint.SetZ( aZ.toDouble( &isZOk ) );
561
562     if ( !isXOk || !isYOk || !isZOk )
563       return false;
564
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() ) )
568       return false;
569
570     // Invert the z value if requested
571     if ( anIsAltitudesInverted )
572       aPoint.SetZ( -aPoint.Z() );
573
574     thePoints.Append( aPoint );
575   }
576
577 #ifdef _TIMER
578   aTimer.Stop();
579   std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
580   aTimer.Show( stream );
581 #endif
582
583   return true;
584 }
585
586 bool HYDROData_Bathymetry::importFromASCFile( QFile&          theFile,
587                                               AltitudePoints& thePoints ) const
588 {
589   if ( !theFile.isOpen() )
590     return false;
591
592   QString aLine;
593   QStringList aStrList;
594
595   int aNCols;
596   int aNRows;
597   double anXllCorner; 
598   double anYllCorner; 
599   double aCellSize; 
600   double aNoDataValue;
601
602   aLine = theFile.readLine().simplified();
603   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
604   if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
605     return false;
606   aNCols = aStrList[1].toInt();
607
608   aLine = theFile.readLine().simplified();
609   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
610   if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
611     return false;
612   aNRows = aStrList[1].toInt();
613
614   aLine = theFile.readLine().simplified();
615   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
616   if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
617     return false;
618   anXllCorner = aStrList[1].toDouble();
619
620   aLine = theFile.readLine().simplified();
621   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
622   if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
623     return false;
624   anYllCorner = aStrList[1].toDouble();
625
626   aLine = theFile.readLine().simplified();
627   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
628   if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
629     return false;
630   aCellSize = aStrList[1].toDouble();
631
632   aLine = theFile.readLine().simplified();
633   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
634   if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
635     return false;
636   aNoDataValue = aStrList[1].toDouble();
637
638   bool anIsAltitudesInverted = IsAltitudesInverted();
639
640   int i = 0;
641   int aStrLength = 0;
642   while ( !theFile.atEnd() )
643   {
644     aLine = theFile.readLine().simplified();
645     aStrList = aLine.split( ' ', QString::SkipEmptyParts );
646
647     aStrLength =  aStrList.length();
648     if ( aStrLength == 0 )
649       continue;
650
651     if ( aStrLength != aNRows )
652       return false;
653
654     for (int j = 0; j < aNCols; j++)
655     {
656       if (aStrList[j].toDouble() != aNoDataValue)
657       {
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());
662
663         if ( anIsAltitudesInverted )
664          aPoint.SetZ( -aPoint.Z() );
665
666         thePoints.Append(aPoint);
667       }
668     }
669     i++;
670
671   }
672
673   return true;
674
675 }
676
677
678 Handle(HYDROData_PolylineXY) HYDROData_Bathymetry::CreateBoundaryPolyline() const
679 {
680   Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
681   Handle(HYDROData_PolylineXY) aResult = 
682     Handle(HYDROData_PolylineXY)::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
683
684   if( aResult.IsNull() )
685     return aResult;
686
687   //search free name
688   QString aPolylinePref = GetName() + "_Boundary";
689   QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
690   aResult->SetName( aPolylineName );
691
692   double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
693   bool isFirst = true;
694   AltitudePoints aPoints = GetAltitudePoints();
695
696   AltitudePoints::Iterator anIter( aPoints );
697   for ( ; anIter.More(); anIter.Next() )
698   {
699     const AltitudePoint& aPoint = anIter.Value();
700
701     double x = aPoint.X(), y = aPoint.Y();
702     if( isFirst || x<Xmin )
703       Xmin = x;
704     if( isFirst || x>Xmax )
705       Xmax = x;
706     if( isFirst || y<Ymin )
707       Ymin = y;
708     if( isFirst || y>Ymax )
709       Ymax = y;
710     isFirst = false;
711   }
712
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 ) );
718   
719   aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
720   
721   aResult->Update();
722
723   return aResult;
724 }
725
726 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
727 {
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() )
732   {
733     AltitudePoint& aPoint = anIter.ChangeValue();
734     aPoint += aDelta;
735   }
736   SetAltitudePoints( aPoints );
737 }
738