20 #ifndef _MYFILTER2D_H_
21 #define _MYFILTER2D_H_
24 #define PI 3.1415926535897932384626433832795
28 #define DEG2RAD(DEG) (DEG/180.0*PI)
32 #define RAD2DEG(RAD) (RAD/PI*180)
37 #define getpid _getpid
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>
53 #include "base/Vector2d.h"
55 #include "imageclasses/ImgReaderGdal.h"
56 #include "imageclasses/ImgWriterGdal.h"
57 #include "algorithms/StatFactory.h"
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};
63 enum RESAMPLE { NEAR = 0, BILINEAR = 1, BICUBIC = 2 };
71 static FILTER_TYPE getFilterType(
const std::string filterType){
72 std::map<std::string, FILTER_TYPE> m_filterMap;
74 return m_filterMap[filterType];
76 static const RESAMPLE getResampleType(
const std::string resampleType){
77 if(resampleType==
"near")
return(NEAR);
78 else if(resampleType==
"bilinear")
return(BILINEAR);
80 std::string errorString=
"resampling type not supported: ";
81 errorString+=resampleType;
82 errorString+=
" use near or bilinear";
89 void pushClass(
short theClass=1){m_class.push_back(theClass);};
90 int pushNoDataValue(
double noDataValue=0);
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;};
101 template<
class T1,
class T2>
void smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dimX,
int dimY);
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>());
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);
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);
132 static void initMap(std::map<std::string, FILTER_TYPE>& m_filterMap){
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;
178 std::vector<short> m_class;
180 std::vector<double> m_noDataValues;
181 std::vector<double> m_threshold;
185 template<
class T1,
class T2>
void Filter2d::smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dim)
187 smooth(inputVector,outputVector,dim,dim);
190 template<
class T1,
class T2>
void Filter2d::smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dimX,
int 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;
198 filter(inputVector,outputVector);
203 outputVector.resize(inputVector.size());
204 int dimX=m_taps[0].size();
205 int dimY=m_taps.size();
207 std::vector<T2> outBuffer(inputVector[0].size());
212 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
213 inBuffer[indexJ]=inputVector[abs(j)];
217 for(
int y=0;y<inputVector.size();++y){
220 inBuffer.erase(inBuffer.begin());
222 if(y+dimY/2<inputVector.size()){
224 inBuffer.push_back(inputVector[y+dimY/2]);
227 int over=y+dimY/2-inputVector.nRows();
228 int index=(inBuffer.size()-1)-over;
230 assert(index<inBuffer.size());
231 inBuffer.push_back(inBuffer[index]);
234 for(
int x=0;x<inputVector.nCols();++x){
236 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
237 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
243 else if(x>=inputVector.nCols()-(dimX-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]);
254 outputVector[y]=outBuffer;
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)
260 const char* pszMessage;
261 void* pProgressArg=NULL;
262 GDALProgressFunc pfnProgress=GDALTermProgress;
264 pfnProgress(progress,pszMessage,pProgressArg);
266 double noDataValue=0;
267 if(m_noDataValues.size())
268 noDataValue=m_noDataValues[0];
274 outputVector.resize((inputVector.size()+down-1)/down);
276 std::vector<T2> outBuffer((inputVector[0].size()+down-1)/down);
280 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
281 inBuffer[indexJ]=inputVector[abs(j)];
284 for(
int y=0;y<inputVector.size();++y){
287 inBuffer.erase(inBuffer.begin());
289 if(y+dimY/2<inputVector.size())
290 inBuffer.push_back(inputVector[y+dimY/2]);
292 int over=y+dimY/2-inputVector.size();
293 int index=(inBuffer.size()-1)-over;
295 assert(index<inBuffer.size());
296 inBuffer.push_back(inBuffer[index]);
299 if((y+1+down/2)%down)
301 for(
int x=0;x<inputVector[0].size();++x){
302 if((x+1+down/2)%down)
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){
314 else if(indexI>=inputVector[0].size())
315 indexI=inputVector[0].size()-i;
318 else if(y+j>=inputVector.size())
319 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
323 for(
int imask=0;imask<m_noDataValues.size();++imask){
324 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
330 std::vector<short>::const_iterator vit=m_class.begin();
333 ++occurrence[inBuffer[indexJ][indexI]];
335 while(vit!=m_class.end()){
336 if(inBuffer[indexJ][indexI]==*(vit++))
337 ++occurrence[inBuffer[indexJ][indexI]];
340 windowBuffer.push_back(inBuffer[indexJ][indexI]);
344 switch(getFilterType(method)){
345 case(filter2d::nvalid):
346 outBuffer[x/down]=stat.nvalid(windowBuffer);
348 case(filter2d::median):
349 if(windowBuffer.empty())
350 outBuffer[x/down]=noDataValue;
352 outBuffer[x/down]=stat.median(windowBuffer);
354 case(filter2d::var):{
355 if(windowBuffer.empty())
356 outBuffer[x/down]=noDataValue;
358 outBuffer[x/down]=stat.var(windowBuffer);
361 case(filter2d::stdev):{
362 if(windowBuffer.empty())
363 outBuffer[x/down]=noDataValue;
365 outBuffer[x/down]=sqrt(stat.var(windowBuffer));
368 case(filter2d::mean):{
369 if(windowBuffer.empty())
370 outBuffer[x/down]=noDataValue;
372 outBuffer[x/down]=stat.mean(windowBuffer);
375 case(filter2d::min):{
376 if(windowBuffer.empty())
377 outBuffer[x/down]=noDataValue;
379 outBuffer[x/down]=stat.mymin(windowBuffer);
382 case(filter2d::ismin):{
383 if(windowBuffer.empty())
384 outBuffer[x/down]=noDataValue;
386 outBuffer[x/down]=(stat.mymin(windowBuffer)==windowBuffer[centre])? 1:0;
389 case(filter2d::minmax):{
392 if(windowBuffer.empty())
393 outBuffer[x/down]=noDataValue;
395 stat.minmax(windowBuffer,windowBuffer.begin(),windowBuffer.end(),min,max);
399 outBuffer[x/down]=windowBuffer[centre];
403 case(filter2d::max):{
404 if(windowBuffer.empty())
405 outBuffer[x/down]=noDataValue;
407 outBuffer[x/down]=stat.mymax(windowBuffer);
410 case(filter2d::ismax):{
411 if(windowBuffer.empty())
412 outBuffer[x/down]=noDataValue;
414 outBuffer[x/down]=(stat.mymax(windowBuffer)==windowBuffer[centre])? 1:0;
417 case(filter2d::order):{
418 if(windowBuffer.empty())
419 outBuffer[x/down]=noDataValue;
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);
430 case(filter2d::sum):{
431 outBuffer[x/down]=stat.sum(windowBuffer);
434 case(filter2d::percentile):{
435 assert(m_threshold.size());
436 outBuffer[x/down]=stat.percentile(windowBuffer,windowBuffer.begin(),windowBuffer.end(),m_threshold[0]);
439 case(filter2d::proportion):{
440 assert(m_threshold.size());
441 double sum=stat.sum(windowBuffer);
443 outBuffer[x/down]=windowBuffer[centre]/sum;
445 outBuffer[x/down]=noDataValue;
448 case(filter2d::homog):
449 if(occurrence.size()==1)
450 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
452 outBuffer[x/down]=noDataValue;
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)
458 else if(*wit!=inBuffer[(dimY-1)/2][x])
460 else if(*wit==inBuffer[(dimY-1)/2][x]){
461 outBuffer[x/down]=noDataValue;
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();
474 outBuffer[x/down]=noDataValue;
477 case(filter2d::countid):{
478 if(windowBuffer.size())
479 outBuffer[x/down]=occurrence.size();
481 outBuffer[x/down]=noDataValue;
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)
491 if(occurrence[inBuffer[(dimY-1)/2][x]]<maxit->second)
492 outBuffer[x/down]=maxit->first;
494 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
497 outBuffer[x/down]=noDataValue;
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];
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];
510 outBuffer[x/down]=noDataValue;
513 case(filter2d::scramble):{
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];
521 outBuffer[x/down]=windowBuffer[randomIndex];
524 outBuffer[x/down]=noDataValue;
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)){
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;
541 outBuffer[x/down]=MF;
547 outBuffer[x/down]=NF;
551 outBuffer[x/down]=inBuffer[indexJ][indexI];
558 progress=(1.0+y/down);
559 progress+=(outputVector.size());
560 progress/=outputVector.size();
561 pfnProgress(progress,pszMessage,pProgressArg);
563 outputVector[y/down]=outBuffer;
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;
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;
587 pfnProgress(progress,pszMessage,pProgressArg);
588 for(
int j=0;j<input.nRows();++j){
589 for(
int i=0;i<input.nCols();++i){
594 randomX=gsl_ran_gaussian(rangen,randomSigma);
595 randomY=gsl_ran_gaussian(rangen,randomSigma);
597 double readCol=i+offsetX+randomX;
598 double readRow=j+offsetY+randomY;
601 if(readRow>input.nRows()-1)
602 readRow=input.nRows()-1;
605 if(readCol>input.nCols()-1)
606 readCol=input.nCols()-1;
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);
618 assert(lowerRow<input.nRows());
620 assert(lowerCol<input.nCols());
622 assert(upperRow<input.nRows());
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;
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;
641 theValue=input[
static_cast<int>(readRow)][static_cast<int>(readCol)];
645 assert(j<output.nRows());
647 assert(i<output.nCols());
648 output[j][i]=theValue;
651 progress/=output.nRows();
652 pfnProgress(progress,pszMessage,pProgressArg);
654 gsl_rng_free(rangen);
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)
659 const char* pszMessage;
660 void* pProgressArg=NULL;
661 GDALProgressFunc pfnProgress=GDALTermProgress;
663 pfnProgress(progress,pszMessage,pProgressArg);
665 double noDataValue=0;
666 if(m_noDataValues.size())
667 noDataValue=m_noDataValues[0];
669 unsigned long int nchange=0;
675 output.resize(input.nRows(),input.nCols());
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];
684 for(
int y=0;y<input.nRows();++y){
687 inBuffer.erase(inBuffer.begin());
689 if(y+dimY/2<input.nRows()){
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];
696 int over=y+dimY/2-input.nRows();
697 int index=(inBuffer.size()-1)-over;
699 assert(index<inBuffer.size());
700 inBuffer.push_back(inBuffer[index]);
703 for(
int x=0;x<input.nCols();++x){
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]){
714 output[y][x]=currentValue;
716 output[y][x]=currentValue;
719 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
720 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
722 if(disc&&(d2>(dimX/2)*(dimY/2)))
728 else if(indexI>=input.nCols())
729 indexI=input.nCols()-i;
732 else if(y+j>=input.nRows())
733 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
736 if(inBuffer[indexJ][indexI]==noDataValue)
739 for(
int imask=0;imask<m_noDataValues.size();++imask){
740 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
747 for(
int iclass=0;iclass<m_class.size();++iclass){
748 if(inBuffer[indexJ][indexI]==m_class[iclass]){
754 statBuffer.push_back(binValue);
756 statBuffer.push_back(inBuffer[indexJ][indexI]);
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);
768 case(filter2d::erode):
769 if(output[y][x]>stat.mymin(statBuffer)+hThreshold){
770 output[y][x]=stat.mymin(statBuffer);
775 std::ostringstream ess;
776 ess <<
"Error: morphology method " << method <<
" not supported, choose " << filter2d::dilate <<
" (dilate) or " << filter2d::erode <<
" (erode)" << std::endl;
780 if(output[y][x]&&m_class.size())
781 output[y][x]=m_class[0];
788 output[y][x]=noDataValue;
792 progress/=output.nRows();
793 pfnProgress(progress,pszMessage,pProgressArg);
798 template<
class T>
unsigned long int Filter2d::dsm2dtm_nwse(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
800 const char* pszMessage;
801 void* pProgressArg=NULL;
802 GDALProgressFunc pfnProgress=GDALTermProgress;
804 pfnProgress(progress,pszMessage,pProgressArg);
807 double noDataValue=0;
808 if(m_noDataValues.size())
809 noDataValue=m_noDataValues[0];
811 unsigned long int nchange=0;
818 if(outputMask.size()!=inputDSM.nRows())
819 outputMask.resize(inputDSM.nRows());
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];
828 for(
int y=0;y<tmpDSM.nRows();++y){
831 inBuffer.erase(inBuffer.begin());
833 if(y+dimY/2<tmpDSM.nRows()){
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];
840 int over=y+dimY/2-tmpDSM.nRows();
841 int index=(inBuffer.size()-1)-over;
843 assert(index<inBuffer.size());
844 inBuffer.push_back(inBuffer[index]);
847 for(
int x=0;x<tmpDSM.nCols();++x){
848 double centerValue=inBuffer[(dimY-1)/2][x];
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){
857 else if(indexI>=tmpDSM.nCols())
858 indexI=tmpDSM.nCols()-i;
861 else if(y+j>=tmpDSM.nRows())
862 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
865 double difference=(centerValue-inBuffer[indexJ][indexI]);
867 neighbors.push_back(inBuffer[indexJ][indexI]);
868 if(difference>hThreshold)
879 sort(neighbors.begin(),neighbors.end());
880 assert(neighbors.size()>1);
881 inBuffer[(dimY-1)/2][x]=neighbors[1];
886 progress/=outputMask.nRows();
887 pfnProgress(progress,pszMessage,pProgressArg);
892 template<
class T>
unsigned long int Filter2d::dsm2dtm_nesw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
894 const char* pszMessage;
895 void* pProgressArg=NULL;
896 GDALProgressFunc pfnProgress=GDALTermProgress;
898 pfnProgress(progress,pszMessage,pProgressArg);
901 double noDataValue=0;
902 if(m_noDataValues.size())
903 noDataValue=m_noDataValues[0];
905 unsigned long int nchange=0;
912 if(outputMask.size()!=inputDSM.nRows())
913 outputMask.resize(inputDSM.nRows());
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];
922 for(
int y=0;y<tmpDSM.nRows();++y){
925 inBuffer.erase(inBuffer.begin());
927 if(y+dimY/2<tmpDSM.nRows()){
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];
934 int over=y+dimY/2-tmpDSM.nRows();
935 int index=(inBuffer.size()-1)-over;
937 assert(index<inBuffer.size());
938 inBuffer.push_back(inBuffer[index]);
941 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
942 double centerValue=inBuffer[(dimY-1)/2][x];
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){
951 else if(indexI>=tmpDSM.nCols())
952 indexI=tmpDSM.nCols()-i;
955 else if(y+j>=tmpDSM.nRows())
956 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
959 double difference=(centerValue-inBuffer[indexJ][indexI]);
961 neighbors.push_back(inBuffer[indexJ][indexI]);
962 if(difference>hThreshold)
973 sort(neighbors.begin(),neighbors.end());
974 assert(neighbors.size()>1);
975 inBuffer[(dimY-1)/2][x]=neighbors[1];
980 progress/=outputMask.nRows();
981 pfnProgress(progress,pszMessage,pProgressArg);
986 template<
class T>
unsigned long int Filter2d::dsm2dtm_senw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
988 const char* pszMessage;
989 void* pProgressArg=NULL;
990 GDALProgressFunc pfnProgress=GDALTermProgress;
992 pfnProgress(progress,pszMessage,pProgressArg);
995 double noDataValue=0;
996 if(m_noDataValues.size())
997 noDataValue=m_noDataValues[0];
999 unsigned long int nchange=0;
1006 if(outputMask.size()!=inputDSM.nRows())
1007 outputMask.resize(inputDSM.nRows());
1009 int indexJ=inputDSM.nRows()-1;
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];
1016 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1017 if(y<tmpDSM.nRows()-1){
1019 inBuffer.erase(inBuffer.end()-1);
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];
1028 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1031 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
1032 double centerValue=inBuffer[(dimY-1)/2][x];
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){
1041 else if(indexI>=tmpDSM.nCols())
1042 indexI=tmpDSM.nCols()-i;
1045 else if(y+j>=tmpDSM.nRows())
1046 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1048 indexJ=(dimY-1)/2+j;
1049 double difference=(centerValue-inBuffer[indexJ][indexI]);
1051 neighbors.push_back(inBuffer[indexJ][indexI]);
1052 if(difference>hThreshold)
1056 if(nmasked<=nlimit){
1063 sort(neighbors.begin(),neighbors.end());
1064 assert(neighbors.size()>1);
1065 inBuffer[(dimY-1)/2][x]=neighbors[1];
1070 progress/=outputMask.nRows();
1071 pfnProgress(progress,pszMessage,pProgressArg);
1076 template<
class T>
unsigned long int Filter2d::dsm2dtm_swne(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
1078 const char* pszMessage;
1079 void* pProgressArg=NULL;
1080 GDALProgressFunc pfnProgress=GDALTermProgress;
1082 pfnProgress(progress,pszMessage,pProgressArg);
1085 double noDataValue=0;
1086 if(m_noDataValues.size())
1087 noDataValue=m_noDataValues[0];
1089 unsigned long int nchange=0;
1096 if(outputMask.size()!=inputDSM.nRows())
1097 outputMask.resize(inputDSM.nRows());
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];
1106 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1107 if(y<tmpDSM.nRows()-1){
1109 inBuffer.erase(inBuffer.end()-1);
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];
1118 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1121 for(
int x=0;x<tmpDSM.nCols();++x){
1122 double centerValue=inBuffer[(dimY-1)/2][x];
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){
1131 else if(indexI>=tmpDSM.nCols())
1132 indexI=tmpDSM.nCols()-i;
1135 else if(y+j>=tmpDSM.nRows())
1136 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1138 indexJ=(dimY-1)/2+j;
1139 double difference=(centerValue-inBuffer[indexJ][indexI]);
1141 neighbors.push_back(inBuffer[indexJ][indexI]);
1142 if(difference>hThreshold)
1146 if(nmasked<=nlimit){
1153 sort(neighbors.begin(),neighbors.end());
1154 assert(neighbors.size()>1);
1155 inBuffer[(dimY-1)/2][x]=neighbors[1];
1160 progress/=outputMask.nRows();
1161 pfnProgress(progress,pszMessage,pProgressArg);
1166 template<
class T>
void Filter2d::shadowDsm(
const Vector2d<T>& input,
Vector2d<T>& output,
double sza,
double saa,
double pixelSize,
short shadowFlag)
1168 unsigned int ncols=input.nCols();
1170 output.resize(input.nRows(),ncols);
1177 const char* pszMessage;
1178 void* pProgressArg=NULL;
1179 GDALProgressFunc pfnProgress=GDALTermProgress;
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)));
1186 double theDir=DEG2RAD(saa)+PI/2.0;
1189 for(
int d=0;d<theDist;++d){
1190 indexI=x+d*cos(theDir);
1191 indexJ=y+d*sin(theDir);
1192 if(indexJ<0||indexJ>=input.size())
1194 if(indexI<0||indexI>=input[indexJ].size())
1196 if(input[indexJ][indexI]<currentValue-d*pixelSize/tan(DEG2RAD(sza))){
1197 output[indexJ][indexI]=shadowFlag;
1202 progress/=output.nRows();
1203 pfnProgress(progress,pszMessage,pProgressArg);
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;
1212 pfnProgress(progress,pszMessage,pProgressArg);
1214 int nRow=theBuffer.size();
1216 int nCol=theBuffer[0].size();
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];
1232 int nsize=theBuffer.size()*theBuffer[0].size();
1234 gsl_wavelet_workspace *work;
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];
1246 progress=(1.0+irow);
1247 progress/=theBuffer.nRows();
1248 pfnProgress(progress,pszMessage,pProgressArg);
1250 gsl_wavelet_free (w);
1251 gsl_wavelet_workspace_free (work);
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;
1259 pfnProgress(progress,pszMessage,pProgressArg);
1261 int nRow=theBuffer.size();
1263 int nCol=theBuffer[0].size();
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]);
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];
1280 int nsize=theBuffer.size()*theBuffer[0].size();
1282 gsl_wavelet_workspace *work;
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];
1294 progress=(1.0+irow);
1295 progress/=theBuffer.nRows();
1296 pfnProgress(progress,pszMessage,pProgressArg);
1298 gsl_wavelet_free (w);
1299 gsl_wavelet_workspace_free (work);
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;
1307 pfnProgress(progress,pszMessage,pProgressArg);
1309 int nRow=theBuffer.size();
1311 int nCol=theBuffer[0].size();
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];
1329 int nsize=theBuffer.size()*theBuffer[0].size();
1331 gsl_wavelet_workspace *work;
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]);
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++)
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];
1352 progress=(1.0+irow);
1353 progress/=theBuffer.nRows();
1354 pfnProgress(progress,pszMessage,pProgressArg);
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());
1362 gsl_wavelet_free (w);
1363 gsl_wavelet_workspace_free (work);