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