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