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