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