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