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