pktools  2.6.6
Processing Kernel for geospatial data
Filter2d.h
1 /**********************************************************************
2 Filter2d.h: class for filtering images
3 Copyright (C) 2008-2012 Pieter Kempeneers
4 
5 This file is part of pktools
6 
7 pktools 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 pktools 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 pktools. If not, see <http://www.gnu.org/licenses/>.
19 ***********************************************************************/
20 #ifndef _MYFILTER2D_H_
21 #define _MYFILTER2D_H_
22 
23 #ifndef PI
24 #define PI 3.1415926535897932384626433832795
25 #endif
26 
27 #ifndef DEG2RAD
28 #define DEG2RAD(DEG) (DEG/180.0*PI)
29 #endif
30 
31 #ifndef RAD2DEG
32 #define RAD2DEG(RAD) (RAD/PI*180)
33 #endif
34 
35 #ifdef WIN32
36 #include <process.h>
37 #define getpid _getpid
38 #endif
39 
40 #include <assert.h>
41 #include <math.h>
42 #include <limits>
43 #include <vector>
44 #include <string>
45 #include <map>
46 extern "C" {
47 #include <gsl/gsl_sort.h>
48 #include <gsl/gsl_wavelet.h>
49 #include <gsl/gsl_wavelet2d.h>
50 #include <gsl/gsl_rng.h>
51 #include <gsl/gsl_randist.h>
52 }
53 #include "base/Vector2d.h"
54 #include "Filter.h"
55 #include "imageclasses/ImgReaderGdal.h"
56 #include "imageclasses/ImgWriterGdal.h"
57 #include "algorithms/StatFactory.h"
58 
59 namespace filter2d
60 {
61  enum FILTER_TYPE { median=100, var=101 , min=102, max=103, sum=104, mean=105, minmax=106, dilate=107, erode=108, close=109, open=110, homog=111, sobelx=112, sobely=113, sobelxy=114, sobelyx=115, smooth=116, density=117, mode=118, mixed=119, threshold=120, ismin=121, ismax=122, heterog=123, order=124, stdev=125, mrf=126, dwt=127, dwti=128, dwt_cut=129, scramble=130, shift=131, linearfeature=132, smoothnodata=133, countid=134, dwt_cut_from=135, savgolay=136, percentile=137, proportion=138, nvalid=139};
62 
63  enum RESAMPLE { NEAR = 0, BILINEAR = 1, BICUBIC = 2 };//bicubic not supported yet...
64 
65 class Filter2d
66 {
67 public:
68  Filter2d(void);
69  Filter2d(const Vector2d<double> &taps);
70  virtual ~Filter2d(){};
71  static FILTER_TYPE getFilterType(const std::string filterType){
72  std::map<std::string, FILTER_TYPE> m_filterMap;
73  initMap(m_filterMap);
74  return m_filterMap[filterType];
75  };
76  static const RESAMPLE getResampleType(const std::string resampleType){
77  if(resampleType=="near") return(NEAR);
78  else if(resampleType=="bilinear") return(BILINEAR);
79  else{
80  std::string errorString="resampling type not supported: ";
81  errorString+=resampleType;
82  errorString+=" use near or bilinear";
83  throw(errorString);
84  }
85  };
86 
87  void setTaps(const Vector2d<double> &taps);
88  /* void setNoValue(double noValue=0){m_noValue=noValue;}; */
89  void pushClass(short theClass=1){m_class.push_back(theClass);};
90  int pushNoDataValue(double noDataValue=0);//{m_mask.push_back(theMask);};
91  void pushThreshold(double theThreshold){m_threshold.push_back(theThreshold);};
92  void setThresholds(const std::vector<double>& theThresholds){m_threshold=theThresholds;};
93  void setClasses(const std::vector<short>& theClasses){m_class=theClasses;};
94  void filter(const ImgReaderGdal& input, ImgWriterGdal& output, bool absolute=false, bool normalize=false, bool noData=false);
95  void smooth(const ImgReaderGdal& input, ImgWriterGdal& output,int dim);
96  void smooth(const ImgReaderGdal& input, ImgWriterGdal& output,int dimX, int dimY);
97  void smoothNoData(const ImgReaderGdal& input, ImgWriterGdal& output,int dim);
98  void smoothNoData(const ImgReaderGdal& input, ImgWriterGdal& output,int dimX, int dimY);
99  template<class T1, class T2> void filter(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector);
100  template<class T1, class T2> void smooth(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector,int dim);
101  template<class T1, class T2> void smooth(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector,int dimX, int dimY);
102  void dwtForward(const ImgReaderGdal& input, ImgWriterGdal& output, const std::string& wavelet_type, int family);
103  void dwtInverse(const ImgReaderGdal& input, ImgWriterGdal& output, const std::string& wavelet_type, int family);
104  void dwtCut(const ImgReaderGdal& input, ImgWriterGdal& output, const std::string& wavelet_type, int family, double cut, bool verbose=false);
105  template<class T> void dwtForward(Vector2d<T>& data, const std::string& wavelet_type, int family);
106  template<class T> void dwtInverse(Vector2d<T>& data, const std::string& wavelet_type, int family);
107  template<class T> void dwtCut(Vector2d<T>& data, const std::string& wavelet_type, int family, double cut);
108  void majorVoting(const std::string& inputFilename, const std::string& outputFilename,int dim=0,const std::vector<int> &prior=std::vector<int>());
109  /* void homogeneousSpatial(const std::string& inputFilename, const std::string& outputFilename, int dim, bool disc=false, int noValue=0); */
110  void doit(const ImgReaderGdal& input, ImgWriterGdal& output, const std::string& method, int dim, short down=1, bool disc=false);
111  void doit(const ImgReaderGdal& input, ImgWriterGdal& output, const std::string& method, int dimX, int dimY, short down=1, bool disc=false);
112  void mrf(const ImgReaderGdal& input, ImgWriterGdal& output, int dimX, int dimY, double beta, bool eightConnectivity=true, short down=1, bool verbose=false);
113  void mrf(const ImgReaderGdal& input, ImgWriterGdal& output, int dimX, int dimY, Vector2d<double> beta, bool eightConnectivity=true, short down=1, bool verbose=false);
114  template<class T1, class T2> void doit(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector, const std::string& method, int dimX, int dimY, short down=1, bool disc=false);
115  void median(const std::string& inputFilename, const std::string& outputFilename, int dim, bool disc=false);
116  void var(const std::string& inputFilename, const std::string& outputFilename, int dim, bool disc=false);
117  void morphology(const ImgReaderGdal& input, ImgWriterGdal& output, const std::string& method, int dimX, int dimY, const std::vector<double> &angle, bool disc=false);
118  template<class T> unsigned long int morphology(const Vector2d<T>& input, Vector2d<T>& output, const std::string& method, int dimX, int dimY, bool disc=false, double hThreshold=0);
119  template<class T> unsigned long int dsm2dtm_nwse(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim=3);
120  template<class T> unsigned long int dsm2dtm_nesw(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim=3);
121  template<class T> unsigned long int dsm2dtm_senw(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim=3);
122  template<class T> unsigned long int dsm2dtm_swne(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim=3);
123  template<class T> void shadowDsm(const Vector2d<T>& input, Vector2d<T>& output, double sza, double saa, double pixelSize, short shadowFlag=1);
124  void shadowDsm(const ImgReaderGdal& input, ImgWriterGdal& output, double sza, double saa, double pixelSize, short shadowFlag=1);
125  void dwt_texture(const std::string& inputFilename, const std::string& outputFilename,int dim, int scale, int down=1, int iband=0, bool verbose=false);
126  void shift(const ImgReaderGdal& input, ImgWriterGdal& output, double offsetX=0, double offsetY=0, double randomSigma=0, RESAMPLE resample=BILINEAR, bool verbose=false);
127  template<class T> void shift(const Vector2d<T>& input, Vector2d<T>& output, double offsetX=0, double offsetY=0, double randomSigma=0, RESAMPLE resample=NEAR, bool verbose=false);
128  void linearFeature(const Vector2d<float>& input, std::vector< Vector2d<float> >& output, float angle=361, float angleStep=1, float maxDistance=0, float eps=0, bool l1=true, bool a1=true, bool l2=true, bool a2=true, bool verbose=false);
129  void linearFeature(const ImgReaderGdal& input, ImgWriterGdal& output, float angle=361, float angleStep=1, float maxDistance=0, float eps=0, bool l1=true, bool a1=true, bool l2=true, bool a2=true, int band=0, bool verbose=false);
130 
131 private:
132  static void initMap(std::map<std::string, FILTER_TYPE>& m_filterMap){
133  //initialize selMap
134  m_filterMap["median"]=filter2d::median;
135  m_filterMap["nvalid"]=filter2d::nvalid;
136  m_filterMap["var"]=filter2d::var;
137  m_filterMap["min"]=filter2d::min;
138  m_filterMap["max"]=filter2d::max;
139  m_filterMap["sum"]=filter2d::sum;
140  m_filterMap["mean"]=filter2d::mean;
141  m_filterMap["minmax"]=filter2d::minmax;
142  m_filterMap["dilate"]=filter2d::dilate;
143  m_filterMap["erode"]=filter2d::erode;
144  m_filterMap["close"]=filter2d::close;
145  m_filterMap["open"]=filter2d::open;
146  m_filterMap["homog"]=filter2d::homog;
147  m_filterMap["sobelx"]=filter2d::sobelx;
148  m_filterMap["sobely"]=filter2d::sobely;
149  m_filterMap["sobelxy"]=filter2d::sobelxy;
150  m_filterMap["sobelyx"]=filter2d::sobelyx;
151  m_filterMap["smooth"]=filter2d::smooth;
152  m_filterMap["density"]=filter2d::density;
153  m_filterMap["mode"]=filter2d::mode;
154  m_filterMap["mixed"]=filter2d::mixed;
155  m_filterMap["smoothnodata"]=filter2d::smoothnodata;
156  m_filterMap["threshold"]=filter2d::threshold;
157  m_filterMap["ismin"]=filter2d::ismin;
158  m_filterMap["ismax"]=filter2d::ismax;
159  m_filterMap["heterog"]=filter2d::heterog;
160  m_filterMap["order"]=filter2d::order;
161  m_filterMap["stdev"]=filter2d::stdev;
162  m_filterMap["mrf"]=filter2d::mrf;
163  m_filterMap["dwt"]=filter2d::dwt;
164  m_filterMap["dwti"]=filter2d::dwti;
165  m_filterMap["dwt_cut"]=filter2d::dwt_cut;
166  m_filterMap["dwt_cut_from"]=filter2d::dwt_cut_from;
167  m_filterMap["scramble"]=filter2d::scramble;
168  m_filterMap["shift"]=filter2d::shift;
169  m_filterMap["linearfeature"]=filter2d::linearfeature;
170  m_filterMap["countid"]=filter2d::countid;
171  m_filterMap["savgolay"]=filter2d::savgolay;
172  m_filterMap["percentile"]=filter2d::percentile;
173  m_filterMap["proportion"]=filter2d::proportion;
174  }
175 
176  Vector2d<double> m_taps;
177  /* double m_noValue; */
178  std::vector<short> m_class;
179  /* std::vector<short> m_mask; */
180  std::vector<double> m_noDataValues;
181  std::vector<double> m_threshold;
182 };
183 
184 
185  template<class T1, class T2> void Filter2d::smooth(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector,int dim)
186  {
187  smooth(inputVector,outputVector,dim,dim);
188  }
189 
190  template<class T1, class T2> void Filter2d::smooth(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector,int dimX, int dimY)
191  {
192  m_taps.resize(dimY);
193  for(int j=0;j<dimY;++j){
194  m_taps[j].resize(dimX);
195  for(int i=0;i<dimX;++i)
196  m_taps[j][i]=1.0/dimX/dimY;
197  }
198  filter(inputVector,outputVector);
199  }
200 
201  template<class T1, class T2> void Filter2d::filter(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector)
202  {
203  outputVector.resize(inputVector.size());
204  int dimX=m_taps[0].size();//horizontal!!!
205  int dimY=m_taps.size();//vertical!!!
206  Vector2d<T1> inBuffer(dimY);
207  std::vector<T2> outBuffer(inputVector[0].size());
208  //initialize last half of inBuffer
209  int indexI=0;
210  int indexJ=0;
211  //initialize last half of inBuffer
212  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
213  inBuffer[indexJ]=inputVector[abs(j)];
214  ++indexJ;
215  }
216 
217  for(int y=0;y<inputVector.size();++y){
218  if(y){//inBuffer already initialized for y=0
219  //erase first line from inBuffer
220  inBuffer.erase(inBuffer.begin());
221  //read extra line and push back to inBuffer if not out of bounds
222  if(y+dimY/2<inputVector.size()){
223  //allocate buffer
224  inBuffer.push_back(inputVector[y+dimY/2]);
225  }
226  else{
227  int over=y+dimY/2-inputVector.nRows();
228  int index=(inBuffer.size()-1)-over;
229  assert(index>=0);
230  assert(index<inBuffer.size());
231  inBuffer.push_back(inBuffer[index]);
232  }
233  }
234  for(int x=0;x<inputVector.nCols();++x){
235  outBuffer[x]=0;
236  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
237  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
238  indexI=x+i;
239  indexJ=(dimY-1)/2+j;
240  //check if out of bounds
241  if(x<(dimX-1)/2)
242  indexI=x+abs(i);
243  else if(x>=inputVector.nCols()-(dimX-1)/2)
244  indexI=x-abs(i);
245  if(y<(dimY-1)/2)
246  indexJ=(dimY-1)/2+abs(j);
247  else if(y>=inputVector.nRows()-(dimY-1)/2)
248  indexJ=(dimY-1)/2-abs(j);
249  outBuffer[x]+=(m_taps[(dimY-1)/2+j][(dimX-1)/2+i]*inBuffer[indexJ][indexI]);
250  }
251  }
252  }
253  //copy outBuffer to outputVector
254  outputVector[y]=outBuffer;
255  }
256  }
257 
258 template<class T1, class T2> void Filter2d::doit(const Vector2d<T1>& inputVector, Vector2d<T2>& outputVector, const std::string& method, int dimX, int dimY, short down, bool disc)
259 {
260  const char* pszMessage;
261  void* pProgressArg=NULL;
262  GDALProgressFunc pfnProgress=GDALTermProgress;
263  double progress=0;
264  pfnProgress(progress,pszMessage,pProgressArg);
265 
266  double noDataValue=0;
267  if(m_noDataValues.size())
268  noDataValue=m_noDataValues[0];
269 
270  assert(dimX);
271  assert(dimY);
272 
274  outputVector.resize((inputVector.size()+down-1)/down);
275  Vector2d<T1> inBuffer(dimY);
276  std::vector<T2> outBuffer((inputVector[0].size()+down-1)/down);
277  int indexI=0;
278  int indexJ=0;
279  //initialize last half of inBuffer
280  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
281  inBuffer[indexJ]=inputVector[abs(j)];
282  ++indexJ;
283  }
284  for(int y=0;y<inputVector.size();++y){
285  if(y){//inBuffer already initialized for y=0
286  //erase first line from inBuffer
287  inBuffer.erase(inBuffer.begin());
288  //read extra line and push back to inBuffer if not out of bounds
289  if(y+dimY/2<inputVector.size())
290  inBuffer.push_back(inputVector[y+dimY/2]);
291  else{
292  int over=y+dimY/2-inputVector.size();
293  int index=(inBuffer.size()-1)-over;
294  assert(index>=0);
295  assert(index<inBuffer.size());
296  inBuffer.push_back(inBuffer[index]);
297  }
298  }
299  if((y+1+down/2)%down)
300  continue;
301  for(int x=0;x<inputVector[0].size();++x){
302  if((x+1+down/2)%down)
303  continue;
304  outBuffer[x/down]=0;
305  std::vector<double> windowBuffer;
306  std::map<int,int> occurrence;
307  int centre=dimX*(dimY-1)/2+(dimX-1)/2;
308  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
309  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
310  indexI=x+i;
311  //check if out of bounds
312  if(indexI<0)
313  indexI=-indexI;
314  else if(indexI>=inputVector[0].size())
315  indexI=inputVector[0].size()-i;
316  if(y+j<0)
317  indexJ=-j;
318  else if(y+j>=inputVector.size())
319  indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
320  else
321  indexJ=(dimY-1)/2+j;
322  bool masked=false;
323  for(int imask=0;imask<m_noDataValues.size();++imask){
324  if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
325  masked=true;
326  break;
327  }
328  }
329  if(!masked){
330  std::vector<short>::const_iterator vit=m_class.begin();
331  //todo: test if this works (only add occurrence if within defined classes)!
332  if(!m_class.size())
333  ++occurrence[inBuffer[indexJ][indexI]];
334  else{
335  while(vit!=m_class.end()){
336  if(inBuffer[indexJ][indexI]==*(vit++))
337  ++occurrence[inBuffer[indexJ][indexI]];
338  }
339  }
340  windowBuffer.push_back(inBuffer[indexJ][indexI]);
341  }
342  }
343  }
344  switch(getFilterType(method)){
345  case(filter2d::nvalid):
346  outBuffer[x/down]=stat.nvalid(windowBuffer);
347  break;
348  case(filter2d::median):
349  if(windowBuffer.empty())
350  outBuffer[x/down]=noDataValue;
351  else
352  outBuffer[x/down]=stat.median(windowBuffer);
353  break;
354  case(filter2d::var):{
355  if(windowBuffer.empty())
356  outBuffer[x/down]=noDataValue;
357  else
358  outBuffer[x/down]=stat.var(windowBuffer);
359  break;
360  }
361  case(filter2d::stdev):{
362  if(windowBuffer.empty())
363  outBuffer[x/down]=noDataValue;
364  else
365  outBuffer[x/down]=sqrt(stat.var(windowBuffer));
366  break;
367  }
368  case(filter2d::mean):{
369  if(windowBuffer.empty())
370  outBuffer[x/down]=noDataValue;
371  else
372  outBuffer[x/down]=stat.mean(windowBuffer);
373  break;
374  }
375  case(filter2d::min):{
376  if(windowBuffer.empty())
377  outBuffer[x/down]=noDataValue;
378  else
379  outBuffer[x/down]=stat.mymin(windowBuffer);
380  break;
381  }
382  case(filter2d::ismin):{
383  if(windowBuffer.empty())
384  outBuffer[x/down]=noDataValue;
385  else
386  outBuffer[x/down]=(stat.mymin(windowBuffer)==windowBuffer[centre])? 1:0;
387  break;
388  }
389  case(filter2d::minmax):{
390  double min=0;
391  double max=0;
392  if(windowBuffer.empty())
393  outBuffer[x/down]=noDataValue;
394  else{
395  stat.minmax(windowBuffer,windowBuffer.begin(),windowBuffer.end(),min,max);
396  if(min!=max)
397  outBuffer[x/down]=0;
398  else
399  outBuffer[x/down]=windowBuffer[centre];//centre pixels
400  }
401  break;
402  }
403  case(filter2d::max):{
404  if(windowBuffer.empty())
405  outBuffer[x/down]=noDataValue;
406  else
407  outBuffer[x/down]=stat.mymax(windowBuffer);
408  break;
409  }
410  case(filter2d::ismax):{
411  if(windowBuffer.empty())
412  outBuffer[x/down]=noDataValue;
413  else
414  outBuffer[x/down]=(stat.mymax(windowBuffer)==windowBuffer[centre])? 1:0;
415  break;
416  }
417  case(filter2d::order):{
418  if(windowBuffer.empty())
419  outBuffer[x/down]=noDataValue;
420  else{
421  double lbound=0;
422  double ubound=dimX*dimY;
423  double theMin=stat.mymin(windowBuffer);
424  double theMax=stat.mymax(windowBuffer);
425  double scale=(ubound-lbound)/(theMax-theMin);
426  outBuffer[x/down]=static_cast<short>(scale*(windowBuffer[centre]-theMin)+lbound);
427  }
428  break;
429  }
430  case(filter2d::sum):{
431  outBuffer[x/down]=stat.sum(windowBuffer);
432  break;
433  }
434  case(filter2d::percentile):{
435  assert(m_threshold.size());
436  outBuffer[x/down]=stat.percentile(windowBuffer,windowBuffer.begin(),windowBuffer.end(),m_threshold[0]);
437  break;
438  }
439  case(filter2d::proportion):{
440  assert(m_threshold.size());
441  double sum=stat.sum(windowBuffer);
442  if(sum)
443  outBuffer[x/down]=windowBuffer[centre]/sum;
444  else
445  outBuffer[x/down]=noDataValue;
446  break;
447  }
448  case(filter2d::homog):
449  if(occurrence.size()==1)//all values in window must be the same
450  outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
451  else//favorize original value in case of ties
452  outBuffer[x/down]=noDataValue;
453  break;
454  case(filter2d::heterog):{
455  for(std::vector<double>::const_iterator wit=windowBuffer.begin();wit!=windowBuffer.end();++wit){
456  if(wit==windowBuffer.begin()+windowBuffer.size()/2)
457  continue;
458  else if(*wit!=inBuffer[(dimY-1)/2][x])
459  outBuffer[x/down]=1;
460  else if(*wit==inBuffer[(dimY-1)/2][x]){//todo:wit mag niet central pixel zijn
461  outBuffer[x/down]=noDataValue;
462  break;
463  }
464  }
465  break;
466  }
467  case(filter2d::density):{
468  if(windowBuffer.size()){
469  std::vector<short>::const_iterator vit=m_class.begin();
470  while(vit!=m_class.end())
471  outBuffer[x/down]+=100.0*occurrence[*(vit++)]/windowBuffer.size();
472  }
473  else
474  outBuffer[x/down]=noDataValue;
475  break;
476  }
477  case(filter2d::countid):{
478  if(windowBuffer.size())
479  outBuffer[x/down]=occurrence.size();
480  else
481  outBuffer[x/down]=noDataValue;
482  break;
483  }
484  case(filter2d::mode):{
485  if(occurrence.size()){
486  std::map<int,int>::const_iterator maxit=occurrence.begin();
487  for(std::map<int,int>::const_iterator mit=occurrence.begin();mit!=occurrence.end();++mit){
488  if(mit->second>maxit->second)
489  maxit=mit;
490  }
491  if(occurrence[inBuffer[(dimY-1)/2][x]]<maxit->second)//
492  outBuffer[x/down]=maxit->first;
493  else//favorize original value in case of ties
494  outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
495  }
496  else
497  outBuffer[x/down]=noDataValue;
498  break;
499  }
500  case(filter2d::threshold):{
501  assert(m_class.size()==m_threshold.size());
502  if(windowBuffer.size()){
503  outBuffer[x/down]=inBuffer[(dimY-1)/2][x];//initialize with original value (in case thresholds not met)
504  for(int iclass=0;iclass<m_class.size();++iclass){
505  if(100.0*(occurrence[m_class[iclass]])/windowBuffer.size()>m_threshold[iclass])
506  outBuffer[x/down]=m_class[iclass];
507  }
508  }
509  else
510  outBuffer[x/down]=noDataValue;
511  break;
512  }
513  case(filter2d::scramble):{//could be done more efficiently window by window with random shuffling entire buffer and assigning entire buffer at once to output image...
514  if(windowBuffer.size()){
515  int randomIndex=std::rand()%windowBuffer.size();
516  if(randomIndex>=windowBuffer.size())
517  outBuffer[x/down]=windowBuffer.back();
518  else if(randomIndex<0)
519  outBuffer[x/down]=windowBuffer[0];
520  else
521  outBuffer[x/down]=windowBuffer[randomIndex];
522  }
523  else
524  outBuffer[x/down]=noDataValue;
525  break;
526  }
527  case(filter2d::mixed):{
528  enum MixType { BF=11, CF=12, MF=13, NF=20, W=30 };
529  double nBF=occurrence[BF];
530  double nCF=occurrence[CF];
531  double nMF=occurrence[MF];
532  double nNF=occurrence[NF];
533  double nW=occurrence[W];
534  if(windowBuffer.size()){
535  if((nBF+nCF+nMF)&&(nBF+nCF+nMF>=nNF+nW)){//forest
536  if(nBF/(nBF+nCF)>=0.75)
537  outBuffer[x/down]=BF;
538  else if(nCF/(nBF+nCF)>=0.75)
539  outBuffer[x/down]=CF;
540  else
541  outBuffer[x/down]=MF;
542  }
543  else{//non-forest
544  if(nW&&(nW>=nNF))
545  outBuffer[x/down]=W;
546  else
547  outBuffer[x/down]=NF;
548  }
549  }
550  else
551  outBuffer[x/down]=inBuffer[indexJ][indexI];
552  break;
553  }
554  default:
555  break;
556  }
557  }
558  progress=(1.0+y/down);
559  progress+=(outputVector.size());
560  progress/=outputVector.size();
561  pfnProgress(progress,pszMessage,pProgressArg);
562  //copy outBuffer to outputVector
563  outputVector[y/down]=outBuffer;
564  }
565 }
566 
567 // class Compare_mapValue{
568 // public:
569 // int operator() (const map<int,int>::value_type& v1, const map<int, int>::value_type& v2) const{
570 // return (v1.second)>(v2.second);
571 // }
572 // };
573 
574 template<class T> void Filter2d::shift(const Vector2d<T>& input, Vector2d<T>& output, double offsetX, double offsetY, double randomSigma, RESAMPLE resample, bool verbose){
575  output.resize(input.nRows(),input.nCols());
576  const gsl_rng_type *rangenType;
577  gsl_rng *rangen;
578  gsl_rng_env_setup();
579  rangenType=gsl_rng_default;
580  rangen=gsl_rng_alloc(rangenType);
581  long seed=time(NULL)*getpid();
582  gsl_rng_set(rangen,seed);
583  const char* pszMessage;
584  void* pProgressArg=NULL;
585  GDALProgressFunc pfnProgress=GDALTermProgress;
586  double progress=0;
587  pfnProgress(progress,pszMessage,pProgressArg);
588  for(int j=0;j<input.nRows();++j){
589  for(int i=0;i<input.nCols();++i){
590  T theValue=0;
591  double randomX=0;
592  double randomY=0;
593  if(randomSigma>0){
594  randomX=gsl_ran_gaussian(rangen,randomSigma);
595  randomY=gsl_ran_gaussian(rangen,randomSigma);
596  }
597  double readCol=i+offsetX+randomX;
598  double readRow=j+offsetY+randomY;
599  if(readRow<0)
600  readRow=0;
601  if(readRow>input.nRows()-1)
602  readRow=input.nRows()-1;
603  if(readCol<0)
604  readCol=0;
605  if(readCol>input.nCols()-1)
606  readCol=input.nCols()-1;
607  switch(resample){
608  case(BILINEAR):{
609  double lowerRow=readRow-0.5;
610  double upperRow=readRow+0.5;
611  lowerRow=static_cast<int>(lowerRow);
612  upperRow=static_cast<int>(upperRow);
613  double lowerCol=readCol-0.5;
614  double upperCol=readCol+0.5;
615  lowerCol=static_cast<int>(lowerCol);
616  upperCol=static_cast<int>(upperCol);
617  assert(lowerRow>=0);
618  assert(lowerRow<input.nRows());
619  assert(lowerCol>=0);
620  assert(lowerCol<input.nCols());
621  assert(upperRow>=0);
622  assert(upperRow<input.nRows());
623  assert(upperCol>=0);
624  if(upperCol>=input.nCols()){
625  std::cout << "upperCol: " << upperCol << std::endl;
626  std::cout << "readCol: " << readCol << std::endl;
627  std::cout << "readCol+0.5: " << readCol+0.5 << std::endl;
628  std::cout << "static_cast<int>(readCol+0.5): " << static_cast<int>(readCol+0.5) << std::endl;
629  }
630  assert(upperCol<input.nCols());
631  double c00=input[lowerRow][lowerCol];
632  double c11=input[upperRow][upperCol];
633  double c01=input[lowerRow][upperCol];
634  double c10=input[upperRow][lowerCol];
635  double a=(upperCol-readCol)*c00+(readCol-lowerCol)*c01;
636  double b=(upperCol-readCol)*c10+(readCol-lowerCol)*c11;
637  theValue=(upperRow-readRow)*a+(readRow-lowerRow)*b;
638  break;
639  }
640  default:
641  theValue=input[static_cast<int>(readRow)][static_cast<int>(readCol)];
642  break;
643  }
644  assert(j>=0);
645  assert(j<output.nRows());
646  assert(i>=0);
647  assert(i<output.nCols());
648  output[j][i]=theValue;
649  }
650  progress=(1.0+j);
651  progress/=output.nRows();
652  pfnProgress(progress,pszMessage,pProgressArg);
653  }
654  gsl_rng_free(rangen);
655 }
656 
657 template<class T> unsigned long int Filter2d::morphology(const Vector2d<T>& input, Vector2d<T>& output, const std::string& method, int dimX, int dimY, bool disc, double hThreshold)
658 {
659  const char* pszMessage;
660  void* pProgressArg=NULL;
661  GDALProgressFunc pfnProgress=GDALTermProgress;
662  double progress=0;
663  pfnProgress(progress,pszMessage,pProgressArg);
664 
665  double noDataValue=0;
666  if(m_noDataValues.size())
667  noDataValue=m_noDataValues[0];
668 
669  unsigned long int nchange=0;
670  assert(dimX);
671  assert(dimY);
673  Vector2d<T> inBuffer(dimY,input.nCols());
674  output.clear();
675  output.resize(input.nRows(),input.nCols());
676  int indexI=0;
677  int indexJ=0;
678  //initialize last half of inBuffer
679  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
680  for(int i=0;i<input.nCols();++i)
681  inBuffer[indexJ][i]=input[abs(j)][i];
682  ++indexJ;
683  }
684  for(int y=0;y<input.nRows();++y){
685  if(y){//inBuffer already initialized for y=0
686  //erase first line from inBuffer
687  inBuffer.erase(inBuffer.begin());
688  //read extra line and push back to inBuffer if not out of bounds
689  if(y+dimY/2<input.nRows()){
690  //allocate buffer
691  inBuffer.push_back(inBuffer.back());
692  for(int i=0;i<input.nCols();++i)
693  inBuffer[inBuffer.size()-1][i]=input[y+dimY/2][i];
694  }
695  else{
696  int over=y+dimY/2-input.nRows();
697  int index=(inBuffer.size()-1)-over;
698  assert(index>=0);
699  assert(index<inBuffer.size());
700  inBuffer.push_back(inBuffer[index]);
701  }
702  }
703  for(int x=0;x<input.nCols();++x){
704  output[y][x]=0;
705  double currentValue=inBuffer[(dimY-1)/2][x];
706  std::vector<double> statBuffer;
707  bool currentMasked=false;
708  for(int imask=0;imask<m_noDataValues.size();++imask){
709  if(currentValue==m_noDataValues[imask]){
710  currentMasked=true;
711  break;
712  }
713  }
714  output[y][x]=currentValue;//introduced due to hThreshold
715  if(currentMasked){
716  output[y][x]=currentValue;
717  }
718  else{
719  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
720  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
721  double d2=i*i+j*j;//square distance
722  if(disc&&(d2>(dimX/2)*(dimY/2)))
723  continue;
724  indexI=x+i;
725  //check if out of bounds
726  if(indexI<0)
727  indexI=-indexI;
728  else if(indexI>=input.nCols())
729  indexI=input.nCols()-i;
730  if(y+j<0)
731  indexJ=-j;
732  else if(y+j>=input.nRows())
733  indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
734  else
735  indexJ=(dimY-1)/2+j;
736  if(inBuffer[indexJ][indexI]==noDataValue)
737  continue;
738  bool masked=false;
739  for(int imask=0;imask<m_noDataValues.size();++imask){
740  if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
741  masked=true;
742  break;
743  }
744  }
745  if(!masked){
746  short binValue=0;
747  for(int iclass=0;iclass<m_class.size();++iclass){
748  if(inBuffer[indexJ][indexI]==m_class[iclass]){
749  binValue=1;
750  break;
751  }
752  }
753  if(m_class.size())
754  statBuffer.push_back(binValue);
755  else
756  statBuffer.push_back(inBuffer[indexJ][indexI]);
757  }
758  }
759  }
760  if(statBuffer.size()){
761  switch(getFilterType(method)){
762  case(filter2d::dilate):
763  if(output[y][x]<stat.mymax(statBuffer)-hThreshold){
764  output[y][x]=stat.mymax(statBuffer);
765  ++nchange;
766  }
767  break;
768  case(filter2d::erode):
769  if(output[y][x]>stat.mymin(statBuffer)+hThreshold){
770  output[y][x]=stat.mymin(statBuffer);
771  ++nchange;
772  }
773  break;
774  default:
775  std::ostringstream ess;
776  ess << "Error: morphology method " << method << " not supported, choose " << filter2d::dilate << " (dilate) or " << filter2d::erode << " (erode)" << std::endl;
777  throw(ess.str());
778  break;
779  }
780  if(output[y][x]&&m_class.size())
781  output[y][x]=m_class[0];
782  // else{
783  // assert(m_noDataValues.size());
784  // output[x]=m_noDataValues[0];
785  // }
786  }
787  else
788  output[y][x]=noDataValue;
789  }
790  }
791  progress=(1.0+y);
792  progress/=output.nRows();
793  pfnProgress(progress,pszMessage,pProgressArg);
794  }
795  return nchange;
796 }
797 
798  template<class T> unsigned long int Filter2d::dsm2dtm_nwse(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim)
799 {
800  const char* pszMessage;
801  void* pProgressArg=NULL;
802  GDALProgressFunc pfnProgress=GDALTermProgress;
803  double progress=0;
804  pfnProgress(progress,pszMessage,pProgressArg);
805 
806  Vector2d<T> tmpDSM(inputDSM);
807  double noDataValue=0;
808  if(m_noDataValues.size())
809  noDataValue=m_noDataValues[0];
810 
811  unsigned long int nchange=0;
812  int dimX=dim;
813  int dimY=dim;
814  assert(dimX);
815  assert(dimY);
817  Vector2d<T> inBuffer(dimY,inputDSM.nCols());
818  if(outputMask.size()!=inputDSM.nRows())
819  outputMask.resize(inputDSM.nRows());
820  int indexI=0;
821  int indexJ=0;
822  //initialize last half of inBuffer
823  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
824  for(int i=0;i<inputDSM.nCols();++i)
825  inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
826  ++indexJ;
827  }
828  for(int y=0;y<tmpDSM.nRows();++y){
829  if(y){//inBuffer already initialized for y=0
830  //erase first line from inBuffer
831  inBuffer.erase(inBuffer.begin());
832  //read extra line and push back to inBuffer if not out of bounds
833  if(y+dimY/2<tmpDSM.nRows()){
834  //allocate buffer
835  inBuffer.push_back(inBuffer.back());
836  for(int i=0;i<tmpDSM.nCols();++i)
837  inBuffer[inBuffer.size()-1][i]=tmpDSM[y+dimY/2][i];
838  }
839  else{
840  int over=y+dimY/2-tmpDSM.nRows();
841  int index=(inBuffer.size()-1)-over;
842  assert(index>=0);
843  assert(index<inBuffer.size());
844  inBuffer.push_back(inBuffer[index]);
845  }
846  }
847  for(int x=0;x<tmpDSM.nCols();++x){
848  double centerValue=inBuffer[(dimY-1)/2][x];
849  short nmasked=0;
850  std::vector<T> neighbors;
851  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
852  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
853  indexI=x+i;
854  //check if out of bounds
855  if(indexI<0)
856  indexI=-indexI;
857  else if(indexI>=tmpDSM.nCols())
858  indexI=tmpDSM.nCols()-i;
859  if(y+j<0)
860  indexJ=-j;
861  else if(y+j>=tmpDSM.nRows())
862  indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
863  else
864  indexJ=(dimY-1)/2+j;
865  double difference=(centerValue-inBuffer[indexJ][indexI]);
866  if(i||j)//skip centerValue
867  neighbors.push_back(inBuffer[indexJ][indexI]);
868  if(difference>hThreshold)
869  ++nmasked;
870  }
871  }
872  if(nmasked<=nlimit){
873  ++nchange;
874  //reset pixel in outputMask
875  outputMask[y][x]=0;
876  }
877  else{
878  //reset pixel height in tmpDSM
879  sort(neighbors.begin(),neighbors.end());
880  assert(neighbors.size()>1);
881  inBuffer[(dimY-1)/2][x]=neighbors[1];
882  /* inBuffer[(dimY-1)/2][x]=stat.mymin(neighbors); */
883  }
884  }
885  progress=(1.0+y);
886  progress/=outputMask.nRows();
887  pfnProgress(progress,pszMessage,pProgressArg);
888  }
889  return nchange;
890 }
891 
892  template<class T> unsigned long int Filter2d::dsm2dtm_nesw(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim)
893 {
894  const char* pszMessage;
895  void* pProgressArg=NULL;
896  GDALProgressFunc pfnProgress=GDALTermProgress;
897  double progress=0;
898  pfnProgress(progress,pszMessage,pProgressArg);
899 
900  Vector2d<T> tmpDSM(inputDSM);
901  double noDataValue=0;
902  if(m_noDataValues.size())
903  noDataValue=m_noDataValues[0];
904 
905  unsigned long int nchange=0;
906  int dimX=dim;
907  int dimY=dim;
908  assert(dimX);
909  assert(dimY);
911  Vector2d<T> inBuffer(dimY,inputDSM.nCols());
912  if(outputMask.size()!=inputDSM.nRows())
913  outputMask.resize(inputDSM.nRows());
914  int indexI=0;
915  int indexJ=0;
916  //initialize last half of inBuffer
917  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
918  for(int i=0;i<inputDSM.nCols();++i)
919  inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
920  ++indexJ;
921  }
922  for(int y=0;y<tmpDSM.nRows();++y){
923  if(y){//inBuffer already initialized for y=0
924  //erase first line from inBuffer
925  inBuffer.erase(inBuffer.begin());
926  //read extra line and push back to inBuffer if not out of bounds
927  if(y+dimY/2<tmpDSM.nRows()){
928  //allocate buffer
929  inBuffer.push_back(inBuffer.back());
930  for(int i=0;i<tmpDSM.nCols();++i)
931  inBuffer[inBuffer.size()-1][i]=tmpDSM[y+dimY/2][i];
932  }
933  else{
934  int over=y+dimY/2-tmpDSM.nRows();
935  int index=(inBuffer.size()-1)-over;
936  assert(index>=0);
937  assert(index<inBuffer.size());
938  inBuffer.push_back(inBuffer[index]);
939  }
940  }
941  for(int x=tmpDSM.nCols()-1;x>=0;--x){
942  double centerValue=inBuffer[(dimY-1)/2][x];
943  short nmasked=0;
944  std::vector<T> neighbors;
945  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
946  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
947  indexI=x+i;
948  //check if out of bounds
949  if(indexI<0)
950  indexI=-indexI;
951  else if(indexI>=tmpDSM.nCols())
952  indexI=tmpDSM.nCols()-i;
953  if(y+j<0)
954  indexJ=-j;
955  else if(y+j>=tmpDSM.nRows())
956  indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
957  else
958  indexJ=(dimY-1)/2+j;
959  double difference=(centerValue-inBuffer[indexJ][indexI]);
960  if(i||j)//skip centerValue
961  neighbors.push_back(inBuffer[indexJ][indexI]);
962  if(difference>hThreshold)
963  ++nmasked;
964  }
965  }
966  if(nmasked<=nlimit){
967  ++nchange;
968  //reset pixel in outputMask
969  outputMask[y][x]=0;
970  }
971  else{
972  //reset pixel height in tmpDSM
973  sort(neighbors.begin(),neighbors.end());
974  assert(neighbors.size()>1);
975  inBuffer[(dimY-1)/2][x]=neighbors[1];
976  /* inBuffer[(dimY-1)/2][x]=stat.mymin(neighbors); */
977  }
978  }
979  progress=(1.0+y);
980  progress/=outputMask.nRows();
981  pfnProgress(progress,pszMessage,pProgressArg);
982  }
983  return nchange;
984 }
985 
986  template<class T> unsigned long int Filter2d::dsm2dtm_senw(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim)
987 {
988  const char* pszMessage;
989  void* pProgressArg=NULL;
990  GDALProgressFunc pfnProgress=GDALTermProgress;
991  double progress=0;
992  pfnProgress(progress,pszMessage,pProgressArg);
993 
994  Vector2d<T> tmpDSM(inputDSM);
995  double noDataValue=0;
996  if(m_noDataValues.size())
997  noDataValue=m_noDataValues[0];
998 
999  unsigned long int nchange=0;
1000  int dimX=dim;
1001  int dimY=dim;
1002  assert(dimX);
1003  assert(dimY);
1005  Vector2d<T> inBuffer(dimY,inputDSM.nCols());
1006  if(outputMask.size()!=inputDSM.nRows())
1007  outputMask.resize(inputDSM.nRows());
1008  int indexI=0;
1009  int indexJ=inputDSM.nRows()-1;
1010  //initialize first half of inBuffer
1011  for(int j=inputDSM.nRows()-dimY/2;j<inputDSM.nRows();--j){
1012  for(int i=0;i<inputDSM.nCols();++i)
1013  inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
1014  ++indexJ;
1015  }
1016  for(int y=tmpDSM.nRows()-1;y>=0;--y){
1017  if(y<tmpDSM.nRows()-1){//inBuffer already initialized for y=tmpDSM.nRows()-1
1018  //erase last line from inBuffer
1019  inBuffer.erase(inBuffer.end()-1);
1020  //read extra line and insert to inBuffer if not out of bounds
1021  if(y-dimY/2>0){
1022  //allocate buffer
1023  inBuffer.insert(inBuffer.begin(),inBuffer.back());
1024  for(int i=0;i<tmpDSM.nCols();++i)
1025  inBuffer[0][i]=tmpDSM[y-dimY/2][i];
1026  }
1027  else{
1028  inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1029  }
1030  }
1031  for(int x=tmpDSM.nCols()-1;x>=0;--x){
1032  double centerValue=inBuffer[(dimY-1)/2][x];
1033  short nmasked=0;
1034  std::vector<T> neighbors;
1035  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
1036  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
1037  indexI=x+i;
1038  //check if out of bounds
1039  if(indexI<0)
1040  indexI=-indexI;
1041  else if(indexI>=tmpDSM.nCols())
1042  indexI=tmpDSM.nCols()-i;
1043  if(y+j<0)
1044  indexJ=-j;
1045  else if(y+j>=tmpDSM.nRows())
1046  indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1047  else
1048  indexJ=(dimY-1)/2+j;
1049  double difference=(centerValue-inBuffer[indexJ][indexI]);
1050  if(i||j)//skip centerValue
1051  neighbors.push_back(inBuffer[indexJ][indexI]);
1052  if(difference>hThreshold)
1053  ++nmasked;
1054  }
1055  }
1056  if(nmasked<=nlimit){
1057  ++nchange;
1058  //reset pixel in outputMask
1059  outputMask[y][x]=0;
1060  }
1061  else{
1062  //reset pixel height in tmpDSM
1063  sort(neighbors.begin(),neighbors.end());
1064  assert(neighbors.size()>1);
1065  inBuffer[(dimY-1)/2][x]=neighbors[1];
1066  /* inBuffer[(dimY-1)/2][x]=stat.mymin(neighbors); */
1067  }
1068  }
1069  progress=(1.0+y);
1070  progress/=outputMask.nRows();
1071  pfnProgress(progress,pszMessage,pProgressArg);
1072  }
1073  return nchange;
1074 }
1075 
1076  template<class T> unsigned long int Filter2d::dsm2dtm_swne(const Vector2d<T>& inputDSM, Vector2d<T>& outputMask, double hThreshold, int nlimit, int dim)
1077 {
1078  const char* pszMessage;
1079  void* pProgressArg=NULL;
1080  GDALProgressFunc pfnProgress=GDALTermProgress;
1081  double progress=0;
1082  pfnProgress(progress,pszMessage,pProgressArg);
1083 
1084  Vector2d<T> tmpDSM(inputDSM);
1085  double noDataValue=0;
1086  if(m_noDataValues.size())
1087  noDataValue=m_noDataValues[0];
1088 
1089  unsigned long int nchange=0;
1090  int dimX=dim;
1091  int dimY=dim;
1092  assert(dimX);
1093  assert(dimY);
1095  Vector2d<T> inBuffer(dimY,inputDSM.nCols());
1096  if(outputMask.size()!=inputDSM.nRows())
1097  outputMask.resize(inputDSM.nRows());
1098  int indexI=0;
1099  int indexJ=0;
1100  //initialize first half of inBuffer
1101  for(int j=inputDSM.nRows()-dimY/2;j<inputDSM.nRows();--j){
1102  for(int i=0;i<inputDSM.nCols();++i)
1103  inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
1104  ++indexJ;
1105  }
1106  for(int y=tmpDSM.nRows()-1;y>=0;--y){
1107  if(y<tmpDSM.nRows()-1){//inBuffer already initialized for y=0
1108  //erase last line from inBuffer
1109  inBuffer.erase(inBuffer.end()-1);
1110  //read extra line and insert to inBuffer if not out of bounds
1111  if(y-dimY/2>0){
1112  //allocate buffer
1113  inBuffer.insert(inBuffer.begin(),inBuffer.back());
1114  for(int i=0;i<tmpDSM.nCols();++i)
1115  inBuffer[0][i]=tmpDSM[y-dimY/2][i];
1116  }
1117  else{
1118  inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1119  }
1120  }
1121  for(int x=0;x<tmpDSM.nCols();++x){
1122  double centerValue=inBuffer[(dimY-1)/2][x];
1123  short nmasked=0;
1124  std::vector<T> neighbors;
1125  for(int j=-(dimY-1)/2;j<=dimY/2;++j){
1126  for(int i=-(dimX-1)/2;i<=dimX/2;++i){
1127  indexI=x+i;
1128  //check if out of bounds
1129  if(indexI<0)
1130  indexI=-indexI;
1131  else if(indexI>=tmpDSM.nCols())
1132  indexI=tmpDSM.nCols()-i;
1133  if(y+j<0)
1134  indexJ=-j;
1135  else if(y+j>=tmpDSM.nRows())
1136  indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1137  else
1138  indexJ=(dimY-1)/2+j;
1139  double difference=(centerValue-inBuffer[indexJ][indexI]);
1140  if(i||j)//skip centerValue
1141  neighbors.push_back(inBuffer[indexJ][indexI]);
1142  if(difference>hThreshold)
1143  ++nmasked;
1144  }
1145  }
1146  if(nmasked<=nlimit){
1147  ++nchange;
1148  //reset pixel in outputMask
1149  outputMask[y][x]=0;
1150  }
1151  else{
1152  //reset pixel height in tmpDSM
1153  sort(neighbors.begin(),neighbors.end());
1154  assert(neighbors.size()>1);
1155  inBuffer[(dimY-1)/2][x]=neighbors[1];
1156  /* inBuffer[(dimY-1)/2][x]=stat.mymin(neighbors); */
1157  }
1158  }
1159  progress=(1.0+y);
1160  progress/=outputMask.nRows();
1161  pfnProgress(progress,pszMessage,pProgressArg);
1162  }
1163  return nchange;
1164 }
1165 
1166  template<class T> void Filter2d::shadowDsm(const Vector2d<T>& input, Vector2d<T>& output, double sza, double saa, double pixelSize, short shadowFlag)
1167 {
1168  unsigned int ncols=input.nCols();
1169  output.clear();
1170  output.resize(input.nRows(),ncols);
1171  //do we need to initialize output?
1172  // for(int y=0;y<output.nRows();++y)
1173  // for(int x=0;x<output.nCols();++x)
1174  // output[y][x]=0;
1175  int indexI=0;
1176  int indexJ=0;
1177  const char* pszMessage;
1178  void* pProgressArg=NULL;
1179  GDALProgressFunc pfnProgress=GDALTermProgress;
1180  double progress=0;
1181  pfnProgress(progress,pszMessage,pProgressArg);
1182  for(int y=0;y<input.nRows();++y){
1183  for(int x=0;x<input.nCols();++x){
1184  double currentValue=input[y][x];
1185  int theDist=static_cast<int>(sqrt((currentValue*tan(DEG2RAD(sza))/pixelSize)*(currentValue*tan(DEG2RAD(sza))/pixelSize)));//in pixels
1186  double theDir=DEG2RAD(saa)+PI/2.0;
1187  if(theDir<0)
1188  theDir+=2*PI;
1189  for(int d=0;d<theDist;++d){//d in pixels
1190  indexI=x+d*cos(theDir);//in pixels
1191  indexJ=y+d*sin(theDir);//in pixels
1192  if(indexJ<0||indexJ>=input.size())
1193  continue;
1194  if(indexI<0||indexI>=input[indexJ].size())
1195  continue;
1196  if(input[indexJ][indexI]<currentValue-d*pixelSize/tan(DEG2RAD(sza))){//in m
1197  output[indexJ][indexI]=shadowFlag;
1198  }
1199  }
1200  }
1201  progress=(1.0+y);
1202  progress/=output.nRows();
1203  pfnProgress(progress,pszMessage,pProgressArg);
1204  }
1205 }
1206 
1207 template<class T> void Filter2d::dwtForward(Vector2d<T>& theBuffer, const std::string& wavelet_type, int family){
1208  const char* pszMessage;
1209  void* pProgressArg=NULL;
1210  GDALProgressFunc pfnProgress=GDALTermProgress;
1211  double progress=0;
1212  pfnProgress(progress,pszMessage,pProgressArg);
1213 
1214  int nRow=theBuffer.size();
1215  assert(nRow);
1216  int nCol=theBuffer[0].size();
1217  assert(nCol);
1218  //make sure data size if power of 2
1219  while(theBuffer.size()&(theBuffer.size()-1))
1220  theBuffer.push_back(theBuffer.back());
1221  for(int irow=0;irow<theBuffer.size();++irow)
1222  while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1223  theBuffer[irow].push_back(theBuffer[irow].back());
1224  std::vector<double> vdata(theBuffer.size()*theBuffer[0].size());
1225  double* data=&(vdata[0]);
1226  for(int irow=0;irow<theBuffer.size();++irow){
1227  for(int icol=0;icol<theBuffer[0].size();++icol){
1228  int index=irow*theBuffer[0].size()+icol;
1229  data[index]=theBuffer[irow][icol];
1230  }
1231  }
1232  int nsize=theBuffer.size()*theBuffer[0].size();
1233  gsl_wavelet *w;
1234  gsl_wavelet_workspace *work;
1235  assert(nsize);
1236  w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1237  work=gsl_wavelet_workspace_alloc(nsize);
1238  gsl_wavelet2d_nstransform_forward (w, data, theBuffer.size(), theBuffer.size(),theBuffer[0].size(), work);
1239  theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1240  for(int irow=0;irow<theBuffer.size();++irow){
1241  theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1242  for(int icol=0;icol<theBuffer[irow].size();++icol){
1243  int index=irow*theBuffer[irow].size()+icol;
1244  theBuffer[irow][icol]=data[index];
1245  }
1246  progress=(1.0+irow);
1247  progress/=theBuffer.nRows();
1248  pfnProgress(progress,pszMessage,pProgressArg);
1249  }
1250  gsl_wavelet_free (w);
1251  gsl_wavelet_workspace_free (work);
1252 }
1253 
1254 template<class T> void Filter2d::dwtInverse(Vector2d<T>& theBuffer, const std::string& wavelet_type, int family){
1255  const char* pszMessage;
1256  void* pProgressArg=NULL;
1257  GDALProgressFunc pfnProgress=GDALTermProgress;
1258  double progress=0;
1259  pfnProgress(progress,pszMessage,pProgressArg);
1260 
1261  int nRow=theBuffer.size();
1262  assert(nRow);
1263  int nCol=theBuffer[0].size();
1264  assert(nCol);
1265  //make sure data size if power of 2
1266  while(theBuffer.size()&(theBuffer.size()-1))
1267  theBuffer.push_back(theBuffer.back());
1268  for(int irow=0;irow<theBuffer.size();++irow)
1269  while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1270  theBuffer[irow].push_back(theBuffer[irow].back());
1271  std::vector<double> vdata(theBuffer.size()*theBuffer[0].size());
1272  double* data=&(vdata[0]);
1273  //double data[theBuffer.size()*theBuffer[0].size()];
1274  for(int irow=0;irow<theBuffer.size();++irow){
1275  for(int icol=0;icol<theBuffer[0].size();++icol){
1276  int index=irow*theBuffer[0].size()+icol;
1277  data[index]=theBuffer[irow][icol];
1278  }
1279  }
1280  int nsize=theBuffer.size()*theBuffer[0].size();
1281  gsl_wavelet *w;
1282  gsl_wavelet_workspace *work;
1283  assert(nsize);
1284  w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1285  work=gsl_wavelet_workspace_alloc(nsize);
1286  gsl_wavelet2d_nstransform_inverse (w, data, theBuffer.size(), theBuffer.size(),theBuffer[0].size(), work);
1287  theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1288  for(int irow=0;irow<theBuffer.size();++irow){
1289  theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1290  for(int icol=0;icol<theBuffer[irow].size();++icol){
1291  int index=irow*theBuffer[irow].size()+icol;
1292  theBuffer[irow][icol]=data[index];
1293  }
1294  progress=(1.0+irow);
1295  progress/=theBuffer.nRows();
1296  pfnProgress(progress,pszMessage,pProgressArg);
1297  }
1298  gsl_wavelet_free (w);
1299  gsl_wavelet_workspace_free (work);
1300 }
1301 
1302 template<class T> void Filter2d::dwtCut(Vector2d<T>& theBuffer, const std::string& wavelet_type, int family, double cut){
1303  const char* pszMessage;
1304  void* pProgressArg=NULL;
1305  GDALProgressFunc pfnProgress=GDALTermProgress;
1306  double progress=0;
1307  pfnProgress(progress,pszMessage,pProgressArg);
1308 
1309  int nRow=theBuffer.size();
1310  assert(nRow);
1311  int nCol=theBuffer[0].size();
1312  assert(nCol);
1313  //make sure data size if power of 2
1314  while(theBuffer.size()&(theBuffer.size()-1))
1315  theBuffer.push_back(theBuffer.back());
1316  for(int irow=0;irow<theBuffer.size();++irow)
1317  while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1318  theBuffer[irow].push_back(theBuffer[irow].back());
1319  double* data=new double[theBuffer.size()*theBuffer[0].size()];
1320  double* abscoeff=new double[theBuffer.size()*theBuffer[0].size()];
1321  size_t* p=new size_t[theBuffer.size()*theBuffer[0].size()];
1322  for(int irow=0;irow<theBuffer.size();++irow){
1323  for(int icol=0;icol<theBuffer[0].size();++icol){
1324  int index=irow*theBuffer[0].size()+icol;
1325  assert(index<theBuffer.size()*theBuffer[0].size());
1326  data[index]=theBuffer[irow][icol];
1327  }
1328  }
1329  int nsize=theBuffer.size()*theBuffer[0].size();
1330  gsl_wavelet *w;
1331  gsl_wavelet_workspace *work;
1332  assert(nsize);
1333  w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1334  work=gsl_wavelet_workspace_alloc(nsize);
1335  gsl_wavelet2d_nstransform_forward (w, data, theBuffer.size(), theBuffer[0].size(),theBuffer[0].size(), work);
1336  for(int irow=0;irow<theBuffer.size();++irow){
1337  for(int icol=0;icol<theBuffer[0].size();++icol){
1338  int index=irow*theBuffer[0].size()+icol;
1339  abscoeff[index]=fabs(data[index]);
1340  }
1341  }
1342  int nc=(100-cut)/100.0*nsize;
1343  gsl_sort_index(p,abscoeff,1,nsize);
1344  for(int i=0;(i+nc)<nsize;i++)
1345  data[p[i]]=0;
1346  gsl_wavelet2d_nstransform_inverse (w, data, theBuffer.size(), theBuffer[0].size(),theBuffer[0].size(), work);
1347  for(int irow=0;irow<theBuffer.size();++irow){
1348  for(int icol=0;icol<theBuffer[irow].size();++icol){
1349  int index=irow*theBuffer[irow].size()+icol;
1350  theBuffer[irow][icol]=data[index];
1351  }
1352  progress=(1.0+irow);
1353  progress/=theBuffer.nRows();
1354  pfnProgress(progress,pszMessage,pProgressArg);
1355  }
1356  theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1357  for(int irow=0;irow<theBuffer.size();++irow)
1358  theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1359  delete[] data;
1360  delete[] abscoeff;
1361  delete[] p;
1362  gsl_wavelet_free (w);
1363  gsl_wavelet_workspace_free (work);
1364 
1365 }
1366 
1367 }
1368 
1369 #endif /* _MYFILTER_H_ */