Salome HOME
Update for occt 7.5.0
[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   Message_ProgressRange aZIRange = HYDROData_Tool::GetZIProgress();
256   if ( aZIRange.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();
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
469 NCollection_Sequence<double> HYDROData_Bathymetry::GetAltitudesForPoints( const NCollection_Sequence<gp_XY>& thePoints, int theMethod) const
470 {
471   DEBTRACE("HYDROData_Bathymetry::GetAltitudesForPoints " << GetName().toStdString());
472   NCollection_Sequence<double> aResSeq;
473   for ( int i = 1, n = thePoints.Length(); i <= n; ++i )
474   {
475     const gp_XY& thePnt = thePoints.Value(i);
476     double anAltitude = GetAltitudeForPoint( thePnt, theMethod );
477     aResSeq.Append( anAltitude );
478   }
479   return aResSeq;
480 }
481
482 double HYDROData_Bathymetry::GetAltitudeForPoint(const gp_XY& thePoint, int theMethod) const
483 {
484   DEBTRACE("GetAltitudeForPoint p(" << thePoint.X() << ", " << thePoint.Y() << "), interpolation method: " << theMethod);
485   double anInvalidAltitude = GetInvalidAltitude();
486   double aResAltitude = anInvalidAltitude;
487
488   // --- find the nearest point in the bathymetry cloud, with quadtree
489   Message_ProgressRange aZIRange = HYDROData_Tool::GetZIProgress();
490   HYDROData_QuadtreeNode* aQuadtree = GetQuadtreeNodes();
491   if (!aQuadtree || (aZIRange.UserBreak()))
492     {
493       DEBTRACE("  no Quadtree");
494       return aResAltitude;
495     }
496
497   std::map<double, const gpi_XYZ*> dist2nodes;
498   aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
499   while (dist2nodes.size() == 0)
500     {
501       aQuadtree->setPrecision(aQuadtree->getPrecision() *2);
502       DEBTRACE("adjust precision to: " << aQuadtree->getPrecision());
503       aQuadtree->NodesAround(thePoint, dist2nodes, aQuadtree->getPrecision());
504     }
505   std::map<double, const gpi_XYZ*>::const_iterator it = dist2nodes.begin();
506   aResAltitude = it->second->Z();
507   int nodeIndex = it->second->getIndex();
508   DEBTRACE("  number of points found: " << dist2nodes.size() << " nearest z: " << aResAltitude << " point index: " << nodeIndex);
509
510   // --- for coarse bathymetry clouds (when the TELEMAC mesh is more refined than the bathymetry cloud)
511   //     interpolation is required.
512   //     - get a Delaunay2D mesh on the bathymetry cloud,
513   //     - get the triangle containing the point in the Delaunay2D mesh,
514   //     - interpolate altitude
515
516   bool isBathyInterpolRequired = false;
517   if (theMethod)
518     isBathyInterpolRequired =true;
519
520 #ifndef LIGHT_MODE
521   if (isBathyInterpolRequired)
522     {
523       vtkPolyData* aDelaunay2D = GetVtkDelaunay2D();
524       vtkIdList* cells= vtkIdList::New();
525       cells->Allocate(64);
526       vtkIdList* points= vtkIdList::New();
527       points->Allocate(64);
528       aDelaunay2D->GetPointCells(nodeIndex, cells);
529       vtkIdType nbCells = cells->GetNumberOfIds();
530       DEBTRACE("  triangles on nearest point: " << nbCells);
531       bool isInside = false;
532       for (int i=0; i<nbCells; i++)
533         {
534           aDelaunay2D->GetCellPoints(cells->GetId(i), points);
535           double z = 0;
536           isInside = interpolZtriangle(thePoint, aDelaunay2D, points, z);
537           if (isInside)
538             {
539               aResAltitude = z;
540               DEBTRACE("  interpolated z: " << z);
541               break;
542             }
543         }
544       if (!isInside)
545       {
546           DEBTRACE("  point outside triangles, nearest z kept");
547       }
548     }
549   #endif
550   return aResAltitude;
551 }
552
553 void HYDROData_Bathymetry::SetFilePath( const TCollection_AsciiString& theFilePath )
554 {
555   Handle(TDataStd_AsciiString) anAttr = TDataStd_AsciiString::Set( myLab.FindChild( DataTag_FilePath ), theFilePath );
556   anAttr->SetID(TDataStd_AsciiString::GetID());
557 }
558
559 void HYDROData_Bathymetry::SetFilePaths( const QStringList& theFilePaths )
560 {
561   int i = 1;
562   Handle_TDataStd_ExtStringArray TExtStrArr = TDataStd_ExtStringArray::Set( myLab.FindChild( DataTag_FilePaths ), 1, theFilePaths.size() );
563   TExtStrArr->SetID(TDataStd_ExtStringArray::GetID());
564   foreach (QString filepath, theFilePaths)
565   {
566     std::string sstr = filepath.toStdString();
567     const char* Val = sstr.c_str();
568     TExtStrArr->SetValue(i, TCollection_ExtendedString(Val));
569     i++;
570   }
571 }
572
573 TCollection_AsciiString HYDROData_Bathymetry::GetFilePath() const
574 {
575   TCollection_AsciiString aRes;
576
577   TDF_Label aLabel = myLab.FindChild( DataTag_FilePath, false );
578   if ( !aLabel.IsNull() )
579   {
580     Handle(TDataStd_AsciiString) anAsciiStr;
581     if ( aLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
582       aRes = anAsciiStr->Get();
583   }
584   else
585   {
586     aLabel = myLab.FindChild( DataTag_FilePaths, false );
587     if ( !aLabel.IsNull() )
588     {
589       Handle(TDataStd_ExtStringArray) anExtStrArr;
590       if ( aLabel.FindAttribute( TDataStd_ExtStringArray::GetID(), anExtStrArr ) )
591         aRes = anExtStrArr->Value(1); //try take the first; convert extstring to asciistring
592     }
593   }
594
595   return aRes;
596 }
597
598 QStringList HYDROData_Bathymetry::GetFilePaths() const
599 {
600   QStringList aResL;
601
602   TDF_Label aLabel = myLab.FindChild( DataTag_FilePaths, false );
603   if ( !aLabel.IsNull() )
604   {
605     Handle(TDataStd_ExtStringArray) anExtStrArr;
606     if ( aLabel.FindAttribute( TDataStd_ExtStringArray::GetID(), anExtStrArr ) )
607     {
608       for (int i = anExtStrArr->Lower(); i <= anExtStrArr->Upper(); i++ )
609       {
610         Standard_ExtString str = anExtStrArr->Value(i).ToExtString();
611         TCollection_AsciiString aText (str);
612         aResL << QString(aText.ToCString());
613       }
614     }
615   }
616   else //backward compatibility
617   {
618     TDF_Label anOldLabel = myLab.FindChild( DataTag_FilePath, false );
619     if ( !anOldLabel.IsNull() )
620     {
621       Handle(TDataStd_AsciiString) anAsciiStr;
622       if ( anOldLabel.FindAttribute( TDataStd_AsciiString::GetID(), anAsciiStr ) )
623         aResL << QString(anAsciiStr->Get().ToCString());
624     }
625   }
626
627   return aResL;
628 }
629
630 void HYDROData_Bathymetry::SetAltitudesInverted( const bool theIsInverted,
631                                                  const bool theIsUpdate )
632 {
633   bool anIsAltitudesInverted = IsAltitudesInverted();
634   if ( anIsAltitudesInverted == theIsInverted )
635     return;
636
637   Handle(TDataStd_Integer) anAttr = TDataStd_Integer::Set( myLab.FindChild( DataTag_AltitudesInverted ), (Standard_Integer)theIsInverted );
638   anAttr->SetID(TDataStd_Integer::GetID());
639   Changed( Geom_Z );
640
641   if ( !theIsUpdate )
642     return;
643
644   // Update altitude points
645   HYDROData_Bathymetry::AltitudePoints anAltitudePoints = GetAltitudePoints();
646   if ( anAltitudePoints.empty() )
647     return;
648
649   HYDROData_Bathymetry::AltitudePoints::iterator anIter = anAltitudePoints.begin(), aLast = anAltitudePoints.end();
650   for ( ; anIter!=aLast; ++anIter )
651   {
652     HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
653     aPoint.Z *= -1;
654   }
655
656   SetAltitudePoints( anAltitudePoints );
657 }
658
659 bool HYDROData_Bathymetry::IsAltitudesInverted() const
660 {
661   bool aRes = false;
662
663   TDF_Label aLabel = myLab.FindChild( DataTag_AltitudesInverted, false );
664   if ( !aLabel.IsNull() )
665   {
666     Handle(TDataStd_Integer) anIntVal;
667     if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anIntVal ) )
668       aRes = (bool)anIntVal->Get();
669   }
670
671   return aRes;
672 }
673
674 bool HYDROData_Bathymetry::ImportFromFile( const QString& theFileName )
675 {
676   return ImportFromFiles(QStringList(theFileName));
677 }
678
679 bool HYDROData_Bathymetry::ImportFromFiles( const QStringList& theFileNames )
680 {
681   AltitudePoints AllPoints;
682   bool Stat = false;
683
684   foreach (QString theFileName, theFileNames)
685   {
686     // Try to open the file
687     QFile aFile( theFileName );
688     if ( !aFile.exists() || !aFile.open( QIODevice::ReadOnly ) )
689       continue;
690
691     QString aFileSuf = QFileInfo( aFile ).suffix().toLower();
692
693     HYDROData_Bathymetry::AltitudePoints aPoints;
694
695     // Try to import the file
696     if ( aFileSuf == "xyz" )
697       Stat = importFromXYZFile( aFile, aPoints );
698     else if ( aFileSuf == "asc" )
699       Stat = importFromASCFile( aFile, aPoints );
700
701     if (!Stat)
702       continue; //ignore this points
703
704     // Close the file
705     aFile.close();
706
707     AllPoints.insert(AllPoints.end(), aPoints.begin(), aPoints.end());
708   }
709
710   // Convert from global to local CS
711   Handle(HYDROData_Document) aDoc = HYDROData_Document::Document();
712   HYDROData_Bathymetry::AltitudePoints::iterator anIter = AllPoints.begin(), aLast = AllPoints.end();
713   for ( ; anIter!=aLast; ++anIter )
714   {
715     HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
716     aDoc->Transform( aPoint.X, aPoint.Y, aPoint.Z, true );
717   }
718
719   if ( Stat )
720   {
721     // Update file path and altitude points of this Bathymetry
722     SetFilePaths (theFileNames );
723     SetAltitudePoints( AllPoints );
724   }
725
726   return Stat && !AllPoints.empty();
727 }
728
729 bool HYDROData_Bathymetry::importFromXYZFile( QFile&          theFile,
730                                               HYDROData_Bathymetry::AltitudePoints& thePoints ) const
731 {
732   if ( !theFile.isOpen() )
733     return false;
734
735   // Strings in file is written as:
736   //  1. X(float) Y(float) Z(float)
737   //  2. X(float) Y(float) Z(float)
738   //  ...
739
740 #ifdef _TIMER
741   OSD_Timer aTimer;
742   aTimer.Start();
743 #endif
744
745   bool anIsAltitudesInverted = IsAltitudesInverted();
746   while ( !theFile.atEnd() )
747   {
748     std::string aLine = theFile.readLine().simplified().toStdString();
749     if ( aLine.empty() )
750       continue;
751
752     HYDROData_Bathymetry::AltitudePoint aPoint;
753     if( sscanf( aLine.c_str(), "%lf %lf %lf", &aPoint.X, &aPoint.Y, &aPoint.Z )!=3 )
754       return false;
755
756     /*QStringList aValues = aLine.split( ' ', QString::SkipEmptyParts );
757     if ( aValues.length() < 3 )
758       return false;
759
760     QString anX = aValues.value( 0 );
761     QString anY = aValues.value( 1 );
762     QString aZ  = aValues.value( 2 );
763
764     bool isXOk = false, isYOk = false, isZOk = false;
765
766     aPoint.X = anX.toDouble( &isXOk );
767     aPoint.Y = anY.toDouble( &isYOk );
768     aPoint.Z = aZ.toDouble( &isZOk );
769
770     if ( !isXOk || !isYOk || !isZOk )
771       return false;*/
772
773     if ( HYDROData_Tool::IsNan( aPoint.X ) || HYDROData_Tool::IsInf( aPoint.X ) ||
774          HYDROData_Tool::IsNan( aPoint.Y ) || HYDROData_Tool::IsInf( aPoint.Y ) ||
775          HYDROData_Tool::IsNan( aPoint.Z ) || HYDROData_Tool::IsInf( aPoint.Z ) )
776       return false;
777
778     // Invert the z value if requested
779     if ( anIsAltitudesInverted )
780       aPoint.Z = -aPoint.Z;
781
782     if( thePoints.size()>=thePoints.capacity() )
783       thePoints.reserve( thePoints.size()+BLOCK_SIZE );
784
785     thePoints.push_back( aPoint );
786   }
787
788 #ifdef _TIMER
789   aTimer.Stop();
790   std::ofstream stream( "W:/HYDRO/WORK/log.txt", std::ofstream::out );
791   aTimer.Show( stream );
792 #endif
793
794   return true;
795 }
796
797 bool HYDROData_Bathymetry::importFromASCFile( QFile&          theFile,
798                                               HYDROData_Bathymetry::AltitudePoints& thePoints ) const
799 {
800   if ( !theFile.isOpen() )
801     return false;
802
803   QString aLine;
804   QStringList aStrList;
805
806   int aNCols;
807   int aNRows;
808   double anXllCorner;
809   double anYllCorner;
810   double aCellSize;
811   double aNoDataValue;
812
813   aLine = theFile.readLine().simplified();
814   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
815   if ( aStrList.length() != 2 && aStrList[0].toLower() != "ncols" )
816     return false;
817   aNCols = aStrList[1].toInt();
818
819   aLine = theFile.readLine().simplified();
820   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
821   if ( aStrList.length() != 2 && aStrList[0].toLower() != "nrows" )
822     return false;
823   aNRows = aStrList[1].toInt();
824
825   aLine = theFile.readLine().simplified();
826   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
827   if ( aStrList.length() != 2 && aStrList[0].toLower() != "xllcorner" )
828     return false;
829   anXllCorner = aStrList[1].toDouble();
830
831   aLine = theFile.readLine().simplified();
832   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
833   if ( aStrList.length() != 2 && aStrList[0].toLower() != "yllcorner" )
834     return false;
835   anYllCorner = aStrList[1].toDouble();
836
837   aLine = theFile.readLine().simplified();
838   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
839   if ( aStrList.length() != 2 && aStrList[0].toLower() != "cellsize" )
840     return false;
841   aCellSize = aStrList[1].toDouble();
842
843   aLine = theFile.readLine().simplified();
844   aStrList = aLine.split( ' ', QString::SkipEmptyParts );
845   if ( aStrList.length() != 2 && aStrList[0].toLower() != "nodata_value" )
846     return false;
847   aNoDataValue = aStrList[1].toDouble();
848
849   bool anIsAltitudesInverted = IsAltitudesInverted();
850
851   int i = 0;
852   int aStrLength = 0;
853   while ( !theFile.atEnd() )
854   {
855     aLine = theFile.readLine().simplified();
856     aStrList = aLine.split( ' ', QString::SkipEmptyParts );
857
858     aStrLength =  aStrList.length();
859     if ( aStrLength == 0 )
860       continue;
861
862     if ( aStrLength != aNRows )
863       return false;
864
865     for (int j = 0; j < aNCols; j++)
866     {
867       if (aStrList[j].toDouble() != aNoDataValue)
868       {
869         HYDROData_Bathymetry::AltitudePoint aPoint;
870         aPoint.X = anXllCorner + aCellSize*(j + 0.5);
871         aPoint.Y = anYllCorner + aCellSize*(aNRows - i + 0.5);
872         aPoint.Z = aStrList[j].toDouble();
873
874         if ( anIsAltitudesInverted )
875          aPoint.Z = -aPoint.Z;
876
877         if( thePoints.size()>=thePoints.capacity() )
878           thePoints.reserve( thePoints.size()+BLOCK_SIZE );
879         thePoints.push_back(aPoint);
880       }
881     }
882     i++;
883   }
884
885   return true;
886 }
887
888
889 Handle(HYDROData_PolylineXY) HYDROData_Bathymetry::CreateBoundaryPolyline() const
890 {
891   Handle(HYDROData_Document) aDocument = HYDROData_Document::Document();
892   Handle(HYDROData_PolylineXY) aResult =
893     Handle(HYDROData_PolylineXY)::DownCast( aDocument->CreateObject( KIND_POLYLINEXY ) );
894
895   if( aResult.IsNull() )
896     return aResult;
897
898   //search free name
899   QString aPolylinePref = GetName() + "_Boundary";
900   QString aPolylineName = HYDROData_Tool::GenerateObjectName( aDocument, aPolylinePref );
901   aResult->SetName( aPolylineName );
902
903   double Xmin = 0.0, Xmax = 0.0, Ymin = 0.0, Ymax = 0.0;
904   bool isFirst = true;
905   HYDROData_Bathymetry::AltitudePoints aPoints = GetAltitudePoints();
906
907   HYDROData_Bathymetry::AltitudePoints::const_iterator anIter = aPoints.begin(), aLast = aPoints.end();
908   for ( ; anIter!=aLast; ++anIter )
909   {
910     const HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
911
912     double x = aPoint.X, y = aPoint.Y;
913     if( isFirst || x<Xmin )
914       Xmin = x;
915     if( isFirst || x>Xmax )
916       Xmax = x;
917     if( isFirst || y<Ymin )
918       Ymin = y;
919     if( isFirst || y>Ymax )
920       Ymax = y;
921     isFirst = false;
922   }
923
924   aResult->AddSection( "bound", HYDROData_IPolyline::SECTION_POLYLINE, true );
925   aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymin ) );
926   aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmin, Ymax ) );
927   aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymax ) );
928   aResult->AddPoint( 0, HYDROData_IPolyline::Point( Xmax, Ymin ) );
929
930   aResult->SetWireColor( HYDROData_PolylineXY::DefaultWireColor() );
931
932   aResult->Update();
933
934   return aResult;
935 }
936
937 void HYDROData_Bathymetry::UpdateLocalCS( double theDx, double theDy )
938 {
939   gp_XYZ aDelta( theDx, theDy, 0 );
940   HYDROData_Bathymetry::AltitudePoints aPoints = GetAltitudePoints();
941   HYDROData_Bathymetry::AltitudePoints::iterator anIter = aPoints.begin(), aLast = aPoints.end();
942   for ( int i = 0; anIter!=aLast; ++i, ++anIter )
943   {
944     HYDROData_Bathymetry::AltitudePoint& aPoint = *anIter;
945     aPoint.X += aDelta.X();
946     aPoint.Y += aDelta.Y();
947     aPoint.Z += aDelta.Z();
948   }
949   SetAltitudePoints( aPoints );
950 }
951