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