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