Subversion Repositories public iLand

Rev

Rev 797 | Rev 877 | Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1
 
671 werner 2
/********************************************************************************************
3
**    iLand - an individual based forest landscape and disturbance model
4
**    http://iland.boku.ac.at
5
**    Copyright (C) 2009-  Werner Rammer, Rupert Seidl
6
**
7
**    This program is free software: you can redistribute it and/or modify
8
**    it under the terms of the GNU General Public License as published by
9
**    the Free Software Foundation, either version 3 of the License, or
10
**    (at your option) any later version.
11
**
12
**    This program is distributed in the hope that it will be useful,
13
**    but WITHOUT ANY WARRANTY; without even the implied warranty of
14
**    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
15
**    GNU General Public License for more details.
16
**
17
**    You should have received a copy of the GNU General Public License
18
**    along with this program.  If not, see <http://www.gnu.org/licenses/>.
19
********************************************************************************************/
20
 
15 Werner 21
#ifndef GRID_H
22
#define GRID_H
23
 
22 Werner 24
#include <QtCore>
15 Werner 25
 
26
 
27
#include <stdexcept>
145 Werner 28
#include <limits>
150 iland 29
#include <cstring>
15 Werner 30
 
373 werner 31
#include "global.h"
32
 
247 werner 33
/** Grid class (template).
697 werner 34
@ingroup tools
74 Werner 35
Orientation
656 werner 36
The grid is oriented as typically coordinates on the northern hemisphere: higher y-values -> north, higher x-values-> east.
490 werner 37
The projection is reversed for drawing on screen (Viewport).
74 Werner 38
          N
490 werner 39
  (0/2) (1/2) (2/2)
656 werner 40
W (0/1) (1/1) (2/1)  E
74 Werner 41
  (0/0) (1/0) (2/0)
42
          S
43
*/
15 Werner 44
template <class T>
45
class Grid {
46
public:
47
 
48
    Grid();
49
    Grid(int cellsize, int sizex, int sizey) { mData=0; setup(cellsize, sizex, sizey); }
58 Werner 50
    Grid(const QRectF rect_metric, const float cellsize) { mData=0; setup(rect_metric,cellsize); }
33 Werner 51
    // copy ctor
52
    Grid(const Grid<T>& toCopy);
105 Werner 53
    ~Grid() { clear(); }
54
    void clear() { if (mData) delete[] mData; mData=0; }
15 Werner 55
 
18 Werner 56
    bool setup(const float cellsize, const int sizex, const int sizey);
22 Werner 57
    bool setup(const QRectF& rect, const double cellsize);
453 werner 58
    bool setup(Grid<T>& source) {     mData = 0;  mRect = source.mRect; return setup(source.cellsize(), source.sizeX(), source.sizeY()); }
75 Werner 59
    void initialize(const T& value) {for( T *p = begin();p!=end(); ++p) *p=value; }
150 iland 60
    void wipe(); ///< write 0-bytes with memcpy to the whole area
154 werner 61
    void wipe(const T value); ///< overwrite the whole area with "value" size of T must be the size of "int" ERRORNOUS!!!
764 werner 62
    /// copies the content of the source grid to this grid.
63
    /// no operation, if the grids are not of the same size.
64
    void copy(const Grid<T> source) { if (source.count()==count()) memcpy(mData, source.mData, count()*sizeof(T)); }
15 Werner 65
 
797 werner 66
    // get the number of cells in x and y direction
145 Werner 67
    int sizeX() const { return mSizeX; }
68
    int sizeY() const { return mSizeY; }
797 werner 69
    // get the size of the grid in metric coordinates (x and y direction)
145 Werner 70
    float metricSizeX() const { return mSizeX*mCellsize; }
71
    float metricSizeY() const { return mSizeY*mCellsize; }
797 werner 72
    /// get the metric rectangle of the grid
49 Werner 73
    QRectF metricRect() const { return mRect; }
797 werner 74
    /// get the rectangle of the grid in terms of indices
75
    QRect rectangle() const { return QRect(QPoint(0,0), QPoint(sizeX(), sizeY())); }
76
    /// get the length of one pixel of the grid
145 Werner 77
    float cellsize() const { return mCellsize; }
373 werner 78
    int count() const { return mCount; } ///< returns the number of elements of the grid
79
    bool isEmpty() const { return mData==NULL; } ///< returns false if the grid was not setup
32 Werner 80
    // operations
15 Werner 81
    // query
33 Werner 82
    /// access (const) with index variables. use int.
83
    inline const T& operator()(const int ix, const int iy) const { return constValueAtIndex(ix, iy); }
84
    /// access (const) using metric variables. use float.
85
    inline const T& operator()(const float x, const float y) const { return constValueAt(x, y); }
48 Werner 86
    inline const T& operator[] (const QPointF &p) const { return constValueAt(p); }
33 Werner 87
 
705 werner 88
    inline T& valueAtIndex(const QPoint& pos) {return valueAtIndex(pos.x(), pos.y());}  ///< value at position defined by a QPoint defining the two indices (x,y)
89
    T& valueAtIndex(const int ix, const int iy) { return mData[iy*mSizeX + ix];  } ///< const value at position defined by indices (x,y)
285 werner 90
    T& valueAtIndex(const int index) {return mData[index]; } ///< get a ref ot value at (one-dimensional) index 'index'.
33 Werner 91
 
705 werner 92
    /// value at position defined by a (integer) QPoint
93
    inline const T& constValueAtIndex(const QPoint& pos) const {return constValueAtIndex(pos.x(), pos.y()); }
94
    /// value at position defined by a pair of integer coordinates
95
    inline const T& constValueAtIndex(const int ix, const int iy) const { return mData[iy*mSizeX + ix];  }
96
    /// value at position defined by the index within the grid
549 werner 97
    const T& constValueAtIndex(const int index) const {return mData[index]; } ///< get a ref ot value at (one-dimensional) index 'index'.
33 Werner 98
 
99
    T& valueAt(const QPointF& posf); ///< value at position defined by metric coordinates (QPointF)
100
    const T& constValueAt(const QPointF& posf) const; ///< value at position defined by metric coordinates (QPointF)
101
 
102
    T& valueAt(const float x, const float y); ///< value at position defined by metric coordinates (x,y)
103
    const T& constValueAt(const float x, const float y) const; ///< value at position defined by metric coordinates (x,y)
104
 
105 Werner 105
    bool coordValid(const float x, const float y) const { return x>=mRect.left() && x<mRect.right()  && y>=mRect.top() && y<mRect.bottom(); }
49 Werner 106
    bool coordValid(const QPointF &pos) const { return coordValid(pos.x(), pos.y()); }
75 Werner 107
 
55 Werner 108
    QPoint indexAt(const QPointF& pos) const { return QPoint(int((pos.x()-mRect.left()) / mCellsize),  int((pos.y()-mRect.top())/mCellsize)); } ///< get index of value at position pos (metric)
538 werner 109
    /// get index (x/y) of the (linear) index 'index' (0..count-1)
110
    QPoint indexOf(const int index) const {return QPoint(index % mSizeX,  index / mSizeX); }
373 werner 111
    bool isIndexValid(const QPoint& pos) const { return (pos.x()>=0 && pos.x()<mSizeX && pos.y()>=0 && pos.y()<mSizeY); } ///< return true, if position is within the grid
112
    bool isIndexValid(const int x, const int y) const {return (x>=0 && x<mSizeX && y>=0 && y<mSizeY); } ///< return true, if index is within the grid
75 Werner 113
    /// force @param pos to contain valid indices with respect to this grid.
55 Werner 114
    void validate(QPoint &pos) const{ pos.setX( qMax(qMin(pos.x(), mSizeX-1), 0) );  pos.setY( qMax(qMin(pos.y(), mSizeY-1), 0) );} ///< ensure that "pos" is a valid key. if out of range, pos is set to minimum/maximum values.
105 Werner 115
    /// get the (metric) centerpoint of cell with index @p pos
549 werner 116
    QPointF cellCenterPoint(const QPoint &pos) const { return QPointF( (pos.x()+0.5)*mCellsize+mRect.left(), (pos.y()+0.5)*mCellsize + mRect.top());} ///< get metric coordinates of the cells center
105 Werner 117
    /// get the metric rectangle of the cell with index @pos
439 werner 118
    QRectF cellRect(const QPoint &pos) const { QRectF r( QPointF(mRect.left() + mCellsize*pos.x(), mRect.top() + pos.y()*mCellsize),
55 Werner 119
                                                   QSizeF(mCellsize, mCellsize)); return r; } ///< return coordinates of rect given by @param pos.
105 Werner 120
 
27 Werner 121
    inline  T* begin() const { return mData; } ///< get "iterator" pointer
37 Werner 122
    inline  T* end() const { return mEnd; } ///< get iterator end-pointer
717 werner 123
    inline QPoint indexOf(const T* element) const; ///< retrieve index (x/y) of the pointer element. returns -1/-1 if element is not valid.
27 Werner 124
    // special queries
33 Werner 125
    T max() const; ///< retrieve the maximum value of a grid
126
    T sum() const; ///< retrieve the sum of the grid
127
    T avg() const; ///< retrieve the average value of a grid
391 werner 128
    // modifying operations
129
    void add(const T& summand);
130
    void multiply(const T& factor);
33 Werner 131
    /// creates a grid with lower resolution and averaged cell values.
132
    /// @param factor factor by which grid size is reduced (e.g. 3 -> 3x3=9 pixels are averaged to 1 result pixel)
133
    /// @param offsetx, offsety: start averaging with an offset from 0/0 (e.g.: x=1, y=2, factor=3: -> 1/2-3/4 -> 0/0)
134
    /// @return Grid with size sizeX()/factor x sizeY()/factor
135
    Grid<T> averaged(const int factor, const int offsetx=0, const int offsety=0) const;
136
    /// normalized returns a normalized grid, in a way that the sum()  = @param targetvalue.
137
    /// if the grid is empty or the sum is 0, no modifications are performed.
138
    Grid<T> normalized(const T targetvalue) const;
373 werner 139
    T* ptr(int x, int y) { return &(mData[y*mSizeX + x]); } ///< get a pointer to the element denoted by "x" and "y"
140
    inline double distance(const QPoint &p1, const QPoint &p2); ///< distance (metric) between p1 and p2
141
    const QPoint randomPosition() const; ///< returns a (valid) random position within the grid
15 Werner 142
private:
77 Werner 143
 
15 Werner 144
    T* mData;
37 Werner 145
    T* mEnd; ///< pointer to 1 element behind the last
49 Werner 146
    QRectF mRect;
36 Werner 147
    float mCellsize; ///< size of a cell in meter
148
    int mSizeX; ///< count of cells in x-direction
149
    int mSizeY; ///< count of cells in y-direction
150
    int mCount; ///< total number of cells in the grid
15 Werner 151
};
152
 
153
typedef Grid<float> FloatGrid;
154
 
643 werner 155
enum GridViewType { GridViewRainbow, GridViewRainbowReverse, GridViewGray, GridViewGrayReverse };
156
 
438 werner 157
/** @class GridRunner is a helper class to iterate over a rectangular fraction of a grid
158
*/
159
template <class T>
160
class GridRunner {
161
public:
650 werner 162
    // constructors with a QRectF (metric coordinates)
617 werner 163
    GridRunner(Grid<T> &target_grid, const QRectF &rectangle) {setup(&target_grid, rectangle);}
164
    GridRunner(const Grid<T> &target_grid, const QRectF &rectangle) {setup(&target_grid, rectangle);}
165
    GridRunner(Grid<T> *target_grid, const QRectF &rectangle) {setup(target_grid, rectangle);}
650 werner 166
    // constructors with a QRect (indices within the grid)
167
    GridRunner(Grid<T> &target_grid, const QRect &rectangle) {setup(&target_grid, rectangle);}
168
    GridRunner(const Grid<T> &target_grid, const QRect &rectangle) {setup(&target_grid, rectangle);}
169
    GridRunner(Grid<T> *target_grid, const QRect &rectangle) {setup(target_grid, rectangle);}
438 werner 170
    T* next(); ///< to to next element, return NULL if finished
662 werner 171
    T* current() const { return mCurrent; }
650 werner 172
    void reset() { mCurrent = mFirst-1; mCurrentCol = -1; }
173
    // helpers
174
    /// fill array with pointers to neighbors (north, east, west, south)
175
    /// or Null-pointers if out of range.
176
    /// the target array (rArray) is not checked and must be valid!
177
    void neighbors4(T** rArray);
178
    void neighbors8(T** rArray);
438 werner 179
private:
617 werner 180
    void setup(const Grid<T> *target_grid, const QRectF &rectangle);
650 werner 181
    void setup(const Grid<T> *target_grid, const QRect &rectangle);
182
    T* mFirst; // points to the first element of the grid
183
    T* mLast; // points to the last element of the grid
438 werner 184
    T* mCurrent;
185
    size_t mLineLength;
186
    size_t mCols;
187
    size_t mCurrentCol;
188
};
189
 
646 werner 190
/** @class Vector3D is a simple 3d vector.
191
  QVector3D (from Qt) is in QtGui so we needed a replacement.
192
*/
193
class Vector3D
194
{
195
 public:
196
    Vector3D(): mX(0.), mY(0.), mZ(0.) {}
197
    Vector3D(const double x, const double y, const double z): mX(x), mY(y), mZ(z) {}
198
    double x() const { return mX; } ///< get x-coordinate
199
    double y() const { return mY; } ///< get y-coordinate
200
    double z() const { return mZ; } ///< get z-coordinate
201
    // set variables
202
    void setX(const double x) { mX=x; } ///< set value of the x-coordinate
203
    void setY(const double y) { mY=y; } ///< set value of the y-coordinate
204
    void setZ(const double z) { mZ=z; } ///< set value of the z-coordinate
205
private:
206
    double mX;
207
    double mY;
208
    double mZ;
209
};
438 werner 210
 
33 Werner 211
// copy constructor
212
template <class T>
213
Grid<T>::Grid(const Grid<T>& toCopy)
214
{
40 Werner 215
    mData = 0;
50 Werner 216
    mRect = toCopy.mRect;
33 Werner 217
    setup(toCopy.cellsize(), toCopy.sizeX(), toCopy.sizeY());
218
    const T* end = toCopy.end();
219
    T* ptr = begin();
220
    for (T* i= toCopy.begin(); i!=end; ++i, ++ptr)
221
       *ptr = *i;
222
}
22 Werner 223
 
33 Werner 224
// normalize function
32 Werner 225
template <class T>
33 Werner 226
Grid<T> Grid<T>::normalized(const T targetvalue) const
32 Werner 227
{
33 Werner 228
    Grid<T> target(*this);
229
    T total = sum();
230
    T multiplier;
231
    if (total)
232
        multiplier = targetvalue / total;
233
    else
234
        return target;
235
    for (T* p=target.begin();p!=target.end();++p)
236
        *p *= multiplier;
40 Werner 237
    return target;
33 Werner 238
}
239
 
240
 
241
template <class T>
242
Grid<T> Grid<T>::averaged(const int factor, const int offsetx, const int offsety) const
243
{
32 Werner 244
    Grid<T> target;
245
    target.setup(cellsize()*factor, sizeX()/factor, sizeY()/factor);
246
    int x,y;
247
    T sum=0;
248
    target.initialize(sum);
249
    // sum over array of 2x2, 3x3, 4x4, ...
250
    for (x=offsetx;x<mSizeX;x++)
251
        for (y=offsety;y<mSizeY;y++) {
252
            target.valueAtIndex((x-offsetx)/factor, (y-offsety)/factor) += constValueAtIndex(x,y);
253
        }
254
    // divide
255
    double fsquare = factor*factor;
256
    for (T* p=target.begin();p!=target.end();++p)
257
        *p /= fsquare;
258
    return target;
259
}
22 Werner 260
 
36 Werner 261
 
27 Werner 262
template <class T>
33 Werner 263
T&  Grid<T>::valueAt(const float x, const float y)
264
{
265
    return valueAtIndex( indexAt(QPointF(x,y)) );
266
}
36 Werner 267
 
33 Werner 268
template <class T>
269
const T&  Grid<T>::constValueAt(const float x, const float y) const
270
{
271
    return constValueAtIndex( indexAt(QPointF(x,y)) );
272
}
36 Werner 273
 
33 Werner 274
template <class T>
22 Werner 275
T&  Grid<T>::valueAt(const QPointF& posf)
276
{
277
    return valueAtIndex( indexAt(posf) );
278
}
36 Werner 279
 
33 Werner 280
template <class T>
281
const T&  Grid<T>::constValueAt(const QPointF& posf) const
282
{
283
    return constValueAtIndex( indexAt(posf) );
284
}
22 Werner 285
 
286
template <class T>
15 Werner 287
Grid<T>::Grid()
288
{
37 Werner 289
    mData = 0; mCellsize=0.f;
290
    mEnd = 0;
15 Werner 291
}
292
 
293
template <class T>
18 Werner 294
bool Grid<T>::setup(const float cellsize, const int sizex, const int sizey)
15 Werner 295
{
37 Werner 296
    mSizeX=sizex; mSizeY=sizey; mCellsize=cellsize;
50 Werner 297
    if (mRect.isNull()) // only set rect if not set before
298
        mRect.setCoords(0., 0., cellsize*sizex, cellsize*sizey);
15 Werner 299
    mCount = mSizeX*mSizeY;
37 Werner 300
    if (mData) {
301
         delete[] mData; mData=NULL;
302
     }
15 Werner 303
   if (mCount>0)
37 Werner 304
        mData = new T[mCount];
305
   mEnd = &(mData[mCount]);
15 Werner 306
   return true;
307
}
308
 
309
template <class T>
22 Werner 310
bool Grid<T>::setup(const QRectF& rect, const double cellsize)
15 Werner 311
{
49 Werner 312
    mRect = rect;
22 Werner 313
    int dx = int(rect.width()/cellsize);
49 Werner 314
    if (mRect.left()+cellsize*dx<rect.right())
22 Werner 315
        dx++;
316
    int dy = int(rect.height()/cellsize);
49 Werner 317
    if (mRect.top()+cellsize*dy<rect.bottom())
22 Werner 318
        dy++;
319
    return setup(cellsize, dx, dy);
15 Werner 320
}
321
 
261 werner 322
/** retrieve from the index from an element reversely from a pointer to that element.
323
    The internal memory layout is (for dimx=6, dimy=3):
324
 
325
6  7  8  9  10 11
326
12 13 14 15 16 17
327
Note: north and south are reversed, thus the item with index 0 is located in the south-western edge of the grid! */
487 werner 328
template <class T> inline
717 werner 329
QPoint Grid<T>::indexOf(const T* element) const
25 Werner 330
{
487 werner 331
//    QPoint result(-1,-1);
25 Werner 332
    if (element==NULL || element<mData || element>=end())
487 werner 333
        return QPoint(-1, -1);
25 Werner 334
    int idx = element - mData;
487 werner 335
    return QPoint(idx % mSizeX,  idx / mSizeX);
336
//    result.setX( idx % mSizeX);
337
//    result.setY( idx / mSizeX);
338
//    return result;
25 Werner 339
}
22 Werner 340
 
27 Werner 341
template <class T>
342
T  Grid<T>::max() const
343
{
143 Werner 344
    T maxv = -std::numeric_limits<T>::max();
27 Werner 345
    T* p;
346
    T* pend = end();
347
    for (p=begin(); p!=pend;++p)
348
       maxv = std::max(maxv, *p);
349
    return maxv;
350
}
351
 
33 Werner 352
template <class T>
353
T  Grid<T>::sum() const
354
{
355
    T* pend = end();
356
    T total = 0;
357
    for (T *p=begin(); p!=pend;++p)
358
       total += *p;
359
    return total;
360
}
361
 
362
template <class T>
363
T  Grid<T>::avg() const
364
{
365
    if (count())
366
        return sum() / T(count());
367
    else return 0;
368
}
369
 
150 iland 370
template <class T>
391 werner 371
void Grid<T>::add(const T& summand)
372
{
373
    T* pend = end();
374
    for (T *p=begin(); p!=pend;*p+=summand,++p)
375
       ;
376
}
377
 
378
template <class T>
379
void Grid<T>::multiply(const T& factor)
380
{
381
    T* pend = end();
382
    for (T *p=begin(); p!=pend;*p*=factor,++p)
383
       ;
384
}
385
 
386
 
387
 
388
template <class T>
150 iland 389
void  Grid<T>::wipe()
390
{
391
    memset(mData, 0, mCount*sizeof(T));
392
}
393
template <class T>
394
void  Grid<T>::wipe(const T value)
395
{
154 werner 396
    /* this does not work properly !!! */
153 werner 397
    if (sizeof(T)==sizeof(int)) {
398
        float temp = value;
399
        float *pf = &temp;
400
 
401
        memset(mData, *((int*)pf), mCount*sizeof(T));
402
    } else
150 iland 403
        initialize(value);
404
}
405
 
373 werner 406
template <class T>
407
double Grid<T>::distance(const QPoint &p1, const QPoint &p2)
408
{
409
    QPointF fp1=cellCenterPoint(p1);
410
    QPointF fp2=cellCenterPoint(p2);
411
    double distance = sqrt( (fp1.x()-fp2.x())*(fp1.x()-fp2.x()) + (fp1.y()-fp2.y())*(fp1.y()-fp2.y()));
412
    return distance;
413
}
414
 
415
template <class T>
416
const QPoint Grid<T>::randomPosition() const
417
{
418
    return QPoint(irandom(0,mSizeX-1), irandom(0, mSizeY-1));
419
}
438 werner 420
 
373 werner 421
////////////////////////////////////////////////////////////
438 werner 422
// grid runner
423
////////////////////////////////////////////////////////////
424
template <class T>
650 werner 425
void GridRunner<T>::setup(const Grid<T> *target_grid, const QRect &rectangle)
438 werner 426
{
650 werner 427
    QPoint upper_left = rectangle.topLeft();
651 werner 428
    // due to the strange behavior of QRect::bottom() and right():
650 werner 429
    QPoint lower_right = rectangle.bottomRight();
617 werner 430
    mCurrent = const_cast<Grid<T> *>(target_grid)->ptr(upper_left.x(), upper_left.y());
650 werner 431
    mFirst = mCurrent;
585 werner 432
    mCurrent--; // point to first element -1
617 werner 433
    mLast = const_cast<Grid<T> *>(target_grid)->ptr(lower_right.x()-1, lower_right.y()-1);
438 werner 434
    mCols = lower_right.x() - upper_left.x(); //
617 werner 435
    mLineLength =  target_grid->sizeX() - mCols;
585 werner 436
    mCurrentCol = -1;
437
//    qDebug() << "GridRunner: rectangle:" << rectangle
438
//             << "upper_left:" << target_grid.cellCenterPoint(target_grid.indexOf(mCurrent))
439
//             << "lower_right:" << target_grid.cellCenterPoint(target_grid.indexOf(mLast));
438 werner 440
}
441
 
442
template <class T>
650 werner 443
void GridRunner<T>::setup(const Grid<T> *target_grid, const QRectF &rectangle_metric)
444
{
445
    QRect rect(target_grid->indexAt(rectangle_metric.topLeft()),
446
               target_grid->indexAt(rectangle_metric.bottomRight()) );
447
    setup (target_grid, rect);
448
}
449
 
450
template <class T>
438 werner 451
T* GridRunner<T>::next()
452
{
453
    if (mCurrent>mLast)
454
        return NULL;
455
    mCurrent++;
456
    mCurrentCol++;
585 werner 457
 
438 werner 458
    if (mCurrentCol >= mCols) {
459
        mCurrent += mLineLength; // skip to next line
460
        mCurrentCol = 0;
461
    }
585 werner 462
    if (mCurrent>mLast)
463
        return NULL;
464
    else
465
        return mCurrent;
438 werner 466
}
467
 
650 werner 468
template <class T>
656 werner 469
/// get pointers the the 4-neighborhood
470
/// north, east, south, west
803 werner 471
/// 0-pointers are returned for edge pixels.
650 werner 472
void GridRunner<T>::neighbors4(T** rArray)
473
{
474
    // north:
651 werner 475
    rArray[0] = mCurrent + mCols + mLineLength > mLast?0: mCurrent + mCols + mLineLength;
650 werner 476
    // south:
651 werner 477
    rArray[3] = mCurrent - (mCols + mLineLength) < mFirst?0: mCurrent -  (mCols + mLineLength);
650 werner 478
    // east / west
656 werner 479
    rArray[1] = mCurrentCol<mCols? mCurrent + 1 : 0;
480
    rArray[2] = mCurrentCol>0? mCurrent-1 : 0;
650 werner 481
}
482
 
483
/// get pointers to the 8-neighbor-hood
484
/// north/east/west/south/NE/NW/SE/SW
803 werner 485
/// 0-pointers are returned for edge pixels.
650 werner 486
template <class T>
487
void GridRunner<T>::neighbors8(T** rArray)
488
{
489
    neighbors4(rArray);
490
    // north-east
656 werner 491
    rArray[4] = rArray[0] && rArray[1]? rArray[0]+1: 0;
650 werner 492
    // north-west
656 werner 493
    rArray[5] = rArray[0] && rArray[2]? rArray[0]-1: 0;
650 werner 494
    // south-east
656 werner 495
    rArray[6] = rArray[3] && rArray[1]? rArray[3]+1: 0;
650 werner 496
    // south-west
656 werner 497
    rArray[7] = rArray[3] && rArray[2]? rArray[3]-1: 0;
650 werner 498
 
499
}
500
 
438 werner 501
////////////////////////////////////////////////////////////
36 Werner 502
// global functions
373 werner 503
////////////////////////////////////////////////////////////
36 Werner 504
 
505
/// dumps a FloatGrid to a String.
46 Werner 506
/// rows will be y-lines, columns x-values. (see grid.cpp)
599 werner 507
QString gridToString(const FloatGrid &grid, const QChar sep=QChar(';'), const int newline_after=-1);
36 Werner 508
 
509
/// creates and return a QImage from Grid-Data.
510
/// @param black_white true: max_value = white, min_value = black, false: color-mode: uses a HSV-color model from blue (min_value) to red (max_value), default: color mode (false)
511
/// @param min_value, max_value min/max bounds for color calcuations. values outside bounds are limited to these values. defaults: min=0, max=1
512
/// @param reverse if true, color ramps are inversed (to: min_value = white (black and white mode) or red (color mode). default = false.
513
/// @return a QImage with the Grids size of pixels. Pixel coordinates relate to the index values of the grid.
514
QImage gridToImage(const FloatGrid &grid,
515
                   bool black_white=false,
516
                   double min_value=0., double max_value=1.,
517
                   bool reverse=false);
518
 
556 werner 519
 
285 werner 520
/** load into 'rGrid' the content of the image pointed at by 'fileName'.
521
    Pixels are converted to grey-scale and then transformend to a value ranging from 0..1 (black..white).
522
  */
523
bool loadGridFromImage(const QString &fileName, FloatGrid &rGrid);
524
 
46 Werner 525
/// template version for non-float grids (see also version for FloatGrid)
599 werner 526
/// @param sep string separator
527
/// @param newline_after if <>-1 a newline is added after every 'newline_after' data values
36 Werner 528
template <class T>
599 werner 529
        QString gridToString(const Grid<T> &grid, const QChar sep=QChar(';'), const int newline_after=-1)
36 Werner 530
{
531
    QString res;
532
    QTextStream ts(&res);
533
 
599 werner 534
    int newl_counter = newline_after;
708 werner 535
    for (int y=grid.sizeY()-1;y>=0;--y){
46 Werner 536
        for (int x=0;x<grid.sizeX();x++){
599 werner 537
            ts << grid.constValueAtIndex(x,y) << sep;
538
            if (--newl_counter==0) {
539
                ts << "\r\n";
540
                newl_counter = newline_after;
541
            }
36 Werner 542
        }
543
        ts << "\r\n";
544
    }
545
 
546
    return res;
547
}
46 Werner 548
 
599 werner 549
/// template version for non-float grids (see also version for FloatGrid)
550
/// @param valueFunction pointer to a function with the signature: QString func(const T&) : this should return a QString
551
/// @param sep string separator
552
/// @param newline_after if <>-1 a newline is added after every 'newline_after' data values
553
template <class T>
554
        QString gridToString(const Grid<T> &grid, QString (*valueFunction)(const T& value), const QChar sep=QChar(';'), const int newline_after=-1 )
708 werner 555
        {
556
            QString res;
557
            QTextStream ts(&res);
599 werner 558
 
708 werner 559
            int newl_counter = newline_after;
560
            for (int y=grid.sizeY()-1;y>=0;--y){
561
                for (int x=0;x<grid.sizeX();x++){
562
                    ts << (*valueFunction)(grid.constValueAtIndex(x,y)) << sep;
599 werner 563
 
708 werner 564
                    if (--newl_counter==0) {
565
                        ts << "\r\n";
566
                        newl_counter = newline_after;
567
                    }
568
                }
599 werner 569
                ts << "\r\n";
570
            }
708 werner 571
 
572
            return res;
599 werner 573
        }
646 werner 574
void modelToWorld(const Vector3D &From, Vector3D &To);
599 werner 575
 
576
template <class T>
577
    QString gridToESRIRaster(const Grid<T> &grid, QString (*valueFunction)(const T& value) )
578
{
646 werner 579
        Vector3D model(grid.metricRect().left(), grid.metricRect().top(), 0.);
580
        Vector3D world;
599 werner 581
        modelToWorld(model, world);
607 werner 582
        QString result = QString("ncols %1\r\nnrows %2\r\nxllcorner %3\r\nyllcorner %4\r\ncellsize %5\r\nNODATA_value %6\r\n")
599 werner 583
                                .arg(grid.sizeX())
584
                                .arg(grid.sizeY())
600 werner 585
                                .arg(world.x(),0,'f').arg(world.y(),0,'f')
599 werner 586
                                .arg(grid.cellsize()).arg(-9999);
694 werner 587
        QString line =  gridToString(grid, valueFunction, QChar(' ')); // for special grids
588
        return result + line;
589
}
599 werner 590
 
591
    template <class T>
592
        QString gridToESRIRaster(const Grid<T> &grid )
593
{
646 werner 594
            Vector3D model(grid.metricRect().left(), grid.metricRect().top(), 0.);
595
            Vector3D world;
599 werner 596
            modelToWorld(model, world);
683 werner 597
            QString result = QString("ncols %1\r\nnrows %2\r\nxllcorner %3\r\nyllcorner %4\r\ncellsize %5\r\nNODATA_value %6\r\n")
599 werner 598
                    .arg(grid.sizeX())
599
                    .arg(grid.sizeY())
680 werner 600
                    .arg(world.x(),0,'f').arg(world.y(),0,'f')
599 werner 601
                    .arg(grid.cellsize()).arg(-9999);
694 werner 602
            QString line = gridToString(grid, QChar(' ')); // for normal grids (e.g. float)
603
            return result + line;
599 werner 604
}
605
 
15 Werner 606
#endif // GRID_H