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