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