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