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