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