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