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};
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[
"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;
177 std::vector<short> m_class;
179 std::vector<double> m_noDataValues;
180 std::vector<double> m_threshold;
184 template<
class T1,
class T2>
void Filter2d::smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dim)
186 smooth(inputVector,outputVector,dim,dim);
189 template<
class T1,
class T2>
void Filter2d::smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dimX,
int 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;
197 filter(inputVector,outputVector);
202 outputVector.resize(inputVector.size());
203 int dimX=m_taps[0].size();
204 int dimY=m_taps.size();
206 std::vector<T2> outBuffer(inputVector[0].size());
211 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
212 inBuffer[indexJ]=inputVector[abs(j)];
216 for(
int y=0;y<inputVector.size();++y){
219 inBuffer.erase(inBuffer.begin());
221 if(y+dimY/2<inputVector.size()){
223 inBuffer.push_back(inputVector[y+dimY/2]);
226 int over=y+dimY/2-inputVector.nRows();
227 int index=(inBuffer.size()-1)-over;
229 assert(index<inBuffer.size());
230 inBuffer.push_back(inBuffer[index]);
233 for(
int x=0;x<inputVector.nCols();++x){
235 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
236 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
242 else if(x>=inputVector.nCols()-(dimX-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]);
253 outputVector[y]=outBuffer;
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)
259 const char* pszMessage;
260 void* pProgressArg=NULL;
261 GDALProgressFunc pfnProgress=GDALTermProgress;
263 pfnProgress(progress,pszMessage,pProgressArg);
265 double noDataValue=0;
266 if(m_noDataValues.size())
267 noDataValue=m_noDataValues[0];
273 outputVector.resize((inputVector.size()+down-1)/down);
275 std::vector<T2> outBuffer((inputVector[0].size()+down-1)/down);
279 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
280 inBuffer[indexJ]=inputVector[abs(j)];
283 for(
int y=0;y<inputVector.size();++y){
286 inBuffer.erase(inBuffer.begin());
288 if(y+dimY/2<inputVector.size())
289 inBuffer.push_back(inputVector[y+dimY/2]);
291 int over=y+dimY/2-inputVector.size();
292 int index=(inBuffer.size()-1)-over;
294 assert(index<inBuffer.size());
295 inBuffer.push_back(inBuffer[index]);
298 if((y+1+down/2)%down)
300 for(
int x=0;x<inputVector[0].size();++x){
301 if((x+1+down/2)%down)
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){
313 else if(indexI>=inputVector[0].size())
314 indexI=inputVector[0].size()-i;
317 else if(y+j>=inputVector.size())
318 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
322 for(
int imask=0;imask<m_noDataValues.size();++imask){
323 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
329 std::vector<short>::const_iterator vit=m_class.begin();
332 ++occurrence[inBuffer[indexJ][indexI]];
334 while(vit!=m_class.end()){
335 if(inBuffer[indexJ][indexI]==*(vit++))
336 ++occurrence[inBuffer[indexJ][indexI]];
339 windowBuffer.push_back(inBuffer[indexJ][indexI]);
343 switch(getFilterType(method)){
344 case(filter2d::median):
345 if(windowBuffer.empty())
346 outBuffer[x/down]=noDataValue;
348 outBuffer[x/down]=stat.median(windowBuffer);
350 case(filter2d::var):{
351 if(windowBuffer.empty())
352 outBuffer[x/down]=noDataValue;
354 outBuffer[x/down]=stat.var(windowBuffer);
357 case(filter2d::stdev):{
358 if(windowBuffer.empty())
359 outBuffer[x/down]=noDataValue;
361 outBuffer[x/down]=sqrt(stat.var(windowBuffer));
364 case(filter2d::mean):{
365 if(windowBuffer.empty())
366 outBuffer[x/down]=noDataValue;
368 outBuffer[x/down]=stat.mean(windowBuffer);
371 case(filter2d::min):{
372 if(windowBuffer.empty())
373 outBuffer[x/down]=noDataValue;
375 outBuffer[x/down]=stat.mymin(windowBuffer);
378 case(filter2d::ismin):{
379 if(windowBuffer.empty())
380 outBuffer[x/down]=noDataValue;
382 outBuffer[x/down]=(stat.mymin(windowBuffer)==windowBuffer[centre])? 1:0;
385 case(filter2d::minmax):{
388 if(windowBuffer.empty())
389 outBuffer[x/down]=noDataValue;
391 stat.minmax(windowBuffer,windowBuffer.begin(),windowBuffer.end(),min,max);
395 outBuffer[x/down]=windowBuffer[centre];
399 case(filter2d::max):{
400 if(windowBuffer.empty())
401 outBuffer[x/down]=noDataValue;
403 outBuffer[x/down]=stat.mymax(windowBuffer);
406 case(filter2d::ismax):{
407 if(windowBuffer.empty())
408 outBuffer[x/down]=noDataValue;
410 outBuffer[x/down]=(stat.mymax(windowBuffer)==windowBuffer[centre])? 1:0;
413 case(filter2d::order):{
414 if(windowBuffer.empty())
415 outBuffer[x/down]=noDataValue;
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);
426 case(filter2d::sum):{
427 outBuffer[x/down]=stat.sum(windowBuffer);
430 case(filter2d::percentile):{
431 assert(m_threshold.size());
432 outBuffer[x/down]=stat.percentile(windowBuffer,windowBuffer.begin(),windowBuffer.end(),m_threshold[0]);
435 case(filter2d::proportion):{
436 assert(m_threshold.size());
437 double sum=stat.sum(windowBuffer);
439 outBuffer[x/down]=windowBuffer[centre]/sum;
441 outBuffer[x/down]=noDataValue;
444 case(filter2d::homog):
445 if(occurrence.size()==1)
446 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
448 outBuffer[x/down]=noDataValue;
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)
454 else if(*wit!=inBuffer[(dimY-1)/2][x])
456 else if(*wit==inBuffer[(dimY-1)/2][x]){
457 outBuffer[x/down]=noDataValue;
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();
470 outBuffer[x/down]=noDataValue;
473 case(filter2d::countid):{
474 if(windowBuffer.size())
475 outBuffer[x/down]=occurrence.size();
477 outBuffer[x/down]=noDataValue;
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)
487 if(occurrence[inBuffer[(dimY-1)/2][x]]<maxit->second)
488 outBuffer[x/down]=maxit->first;
490 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
493 outBuffer[x/down]=noDataValue;
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];
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];
506 outBuffer[x/down]=noDataValue;
509 case(filter2d::scramble):{
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];
517 outBuffer[x/down]=windowBuffer[randomIndex];
520 outBuffer[x/down]=noDataValue;
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)){
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;
537 outBuffer[x/down]=MF;
543 outBuffer[x/down]=NF;
547 outBuffer[x/down]=inBuffer[indexJ][indexI];
554 progress=(1.0+y/down);
555 progress+=(outputVector.size());
556 progress/=outputVector.size();
557 pfnProgress(progress,pszMessage,pProgressArg);
559 outputVector[y/down]=outBuffer;
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;
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;
583 pfnProgress(progress,pszMessage,pProgressArg);
584 for(
int j=0;j<input.nRows();++j){
585 for(
int i=0;i<input.nCols();++i){
590 randomX=gsl_ran_gaussian(rangen,randomSigma);
591 randomY=gsl_ran_gaussian(rangen,randomSigma);
593 double readCol=i+offsetX+randomX;
594 double readRow=j+offsetY+randomY;
597 if(readRow>input.nRows()-1)
598 readRow=input.nRows()-1;
601 if(readCol>input.nCols()-1)
602 readCol=input.nCols()-1;
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);
614 assert(lowerRow<input.nRows());
616 assert(lowerCol<input.nCols());
618 assert(upperRow<input.nRows());
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;
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;
637 theValue=input[
static_cast<int>(readRow)][static_cast<int>(readCol)];
641 assert(j<output.nRows());
643 assert(i<output.nCols());
644 output[j][i]=theValue;
647 progress/=output.nRows();
648 pfnProgress(progress,pszMessage,pProgressArg);
650 gsl_rng_free(rangen);
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)
655 const char* pszMessage;
656 void* pProgressArg=NULL;
657 GDALProgressFunc pfnProgress=GDALTermProgress;
659 pfnProgress(progress,pszMessage,pProgressArg);
661 double noDataValue=0;
662 if(m_noDataValues.size())
663 noDataValue=m_noDataValues[0];
665 unsigned long int nchange=0;
671 output.resize(input.nRows(),input.nCols());
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];
680 for(
int y=0;y<input.nRows();++y){
683 inBuffer.erase(inBuffer.begin());
685 if(y+dimY/2<input.nRows()){
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];
692 int over=y+dimY/2-input.nRows();
693 int index=(inBuffer.size()-1)-over;
695 assert(index<inBuffer.size());
696 inBuffer.push_back(inBuffer[index]);
699 for(
int x=0;x<input.nCols();++x){
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]){
710 output[y][x]=currentValue;
712 output[y][x]=currentValue;
715 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
716 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
718 if(disc&&(d2>(dimX/2)*(dimY/2)))
724 else if(indexI>=input.nCols())
725 indexI=input.nCols()-i;
728 else if(y+j>=input.nRows())
729 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
732 if(inBuffer[indexJ][indexI]==noDataValue)
735 for(
int imask=0;imask<m_noDataValues.size();++imask){
736 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
743 for(
int iclass=0;iclass<m_class.size();++iclass){
744 if(inBuffer[indexJ][indexI]==m_class[iclass]){
750 statBuffer.push_back(binValue);
752 statBuffer.push_back(inBuffer[indexJ][indexI]);
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);
764 case(filter2d::erode):
765 if(output[y][x]>stat.mymin(statBuffer)+hThreshold){
766 output[y][x]=stat.mymin(statBuffer);
771 std::ostringstream ess;
772 ess <<
"Error: morphology method " << method <<
" not supported, choose " << filter2d::dilate <<
" (dilate) or " << filter2d::erode <<
" (erode)" << std::endl;
776 if(output[y][x]&&m_class.size())
777 output[y][x]=m_class[0];
784 output[y][x]=noDataValue;
788 progress/=output.nRows();
789 pfnProgress(progress,pszMessage,pProgressArg);
794 template<
class T>
unsigned long int Filter2d::dsm2dtm_nwse(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
796 const char* pszMessage;
797 void* pProgressArg=NULL;
798 GDALProgressFunc pfnProgress=GDALTermProgress;
800 pfnProgress(progress,pszMessage,pProgressArg);
803 double noDataValue=0;
804 if(m_noDataValues.size())
805 noDataValue=m_noDataValues[0];
807 unsigned long int nchange=0;
814 if(outputMask.size()!=inputDSM.nRows())
815 outputMask.resize(inputDSM.nRows());
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];
824 for(
int y=0;y<tmpDSM.nRows();++y){
827 inBuffer.erase(inBuffer.begin());
829 if(y+dimY/2<tmpDSM.nRows()){
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];
836 int over=y+dimY/2-tmpDSM.nRows();
837 int index=(inBuffer.size()-1)-over;
839 assert(index<inBuffer.size());
840 inBuffer.push_back(inBuffer[index]);
843 for(
int x=0;x<tmpDSM.nCols();++x){
844 double centerValue=inBuffer[(dimY-1)/2][x];
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){
853 else if(indexI>=tmpDSM.nCols())
854 indexI=tmpDSM.nCols()-i;
857 else if(y+j>=tmpDSM.nRows())
858 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
861 double difference=(centerValue-inBuffer[indexJ][indexI]);
863 neighbors.push_back(inBuffer[indexJ][indexI]);
864 if(difference>hThreshold)
875 sort(neighbors.begin(),neighbors.end());
876 assert(neighbors.size()>1);
877 inBuffer[(dimY-1)/2][x]=neighbors[1];
882 progress/=outputMask.nRows();
883 pfnProgress(progress,pszMessage,pProgressArg);
888 template<
class T>
unsigned long int Filter2d::dsm2dtm_nesw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
890 const char* pszMessage;
891 void* pProgressArg=NULL;
892 GDALProgressFunc pfnProgress=GDALTermProgress;
894 pfnProgress(progress,pszMessage,pProgressArg);
897 double noDataValue=0;
898 if(m_noDataValues.size())
899 noDataValue=m_noDataValues[0];
901 unsigned long int nchange=0;
908 if(outputMask.size()!=inputDSM.nRows())
909 outputMask.resize(inputDSM.nRows());
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];
918 for(
int y=0;y<tmpDSM.nRows();++y){
921 inBuffer.erase(inBuffer.begin());
923 if(y+dimY/2<tmpDSM.nRows()){
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];
930 int over=y+dimY/2-tmpDSM.nRows();
931 int index=(inBuffer.size()-1)-over;
933 assert(index<inBuffer.size());
934 inBuffer.push_back(inBuffer[index]);
937 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
938 double centerValue=inBuffer[(dimY-1)/2][x];
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){
947 else if(indexI>=tmpDSM.nCols())
948 indexI=tmpDSM.nCols()-i;
951 else if(y+j>=tmpDSM.nRows())
952 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
955 double difference=(centerValue-inBuffer[indexJ][indexI]);
957 neighbors.push_back(inBuffer[indexJ][indexI]);
958 if(difference>hThreshold)
969 sort(neighbors.begin(),neighbors.end());
970 assert(neighbors.size()>1);
971 inBuffer[(dimY-1)/2][x]=neighbors[1];
976 progress/=outputMask.nRows();
977 pfnProgress(progress,pszMessage,pProgressArg);
982 template<
class T>
unsigned long int Filter2d::dsm2dtm_senw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
984 const char* pszMessage;
985 void* pProgressArg=NULL;
986 GDALProgressFunc pfnProgress=GDALTermProgress;
988 pfnProgress(progress,pszMessage,pProgressArg);
991 double noDataValue=0;
992 if(m_noDataValues.size())
993 noDataValue=m_noDataValues[0];
995 unsigned long int nchange=0;
1002 if(outputMask.size()!=inputDSM.nRows())
1003 outputMask.resize(inputDSM.nRows());
1005 int indexJ=inputDSM.nRows()-1;
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];
1012 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1013 if(y<tmpDSM.nRows()-1){
1015 inBuffer.erase(inBuffer.end()-1);
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];
1024 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1027 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
1028 double centerValue=inBuffer[(dimY-1)/2][x];
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){
1037 else if(indexI>=tmpDSM.nCols())
1038 indexI=tmpDSM.nCols()-i;
1041 else if(y+j>=tmpDSM.nRows())
1042 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1044 indexJ=(dimY-1)/2+j;
1045 double difference=(centerValue-inBuffer[indexJ][indexI]);
1047 neighbors.push_back(inBuffer[indexJ][indexI]);
1048 if(difference>hThreshold)
1052 if(nmasked<=nlimit){
1059 sort(neighbors.begin(),neighbors.end());
1060 assert(neighbors.size()>1);
1061 inBuffer[(dimY-1)/2][x]=neighbors[1];
1066 progress/=outputMask.nRows();
1067 pfnProgress(progress,pszMessage,pProgressArg);
1072 template<
class T>
unsigned long int Filter2d::dsm2dtm_swne(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
1074 const char* pszMessage;
1075 void* pProgressArg=NULL;
1076 GDALProgressFunc pfnProgress=GDALTermProgress;
1078 pfnProgress(progress,pszMessage,pProgressArg);
1081 double noDataValue=0;
1082 if(m_noDataValues.size())
1083 noDataValue=m_noDataValues[0];
1085 unsigned long int nchange=0;
1092 if(outputMask.size()!=inputDSM.nRows())
1093 outputMask.resize(inputDSM.nRows());
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];
1102 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1103 if(y<tmpDSM.nRows()-1){
1105 inBuffer.erase(inBuffer.end()-1);
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];
1114 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1117 for(
int x=0;x<tmpDSM.nCols();++x){
1118 double centerValue=inBuffer[(dimY-1)/2][x];
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){
1127 else if(indexI>=tmpDSM.nCols())
1128 indexI=tmpDSM.nCols()-i;
1131 else if(y+j>=tmpDSM.nRows())
1132 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1134 indexJ=(dimY-1)/2+j;
1135 double difference=(centerValue-inBuffer[indexJ][indexI]);
1137 neighbors.push_back(inBuffer[indexJ][indexI]);
1138 if(difference>hThreshold)
1142 if(nmasked<=nlimit){
1149 sort(neighbors.begin(),neighbors.end());
1150 assert(neighbors.size()>1);
1151 inBuffer[(dimY-1)/2][x]=neighbors[1];
1156 progress/=outputMask.nRows();
1157 pfnProgress(progress,pszMessage,pProgressArg);
1162 template<
class T>
void Filter2d::shadowDsm(
const Vector2d<T>& input,
Vector2d<T>& output,
double sza,
double saa,
double pixelSize,
short shadowFlag)
1164 unsigned int ncols=input.nCols();
1166 output.resize(input.nRows(),ncols);
1173 const char* pszMessage;
1174 void* pProgressArg=NULL;
1175 GDALProgressFunc pfnProgress=GDALTermProgress;
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)));
1182 double theDir=DEG2RAD(saa)+PI/2.0;
1185 for(
int d=0;d<theDist;++d){
1186 indexI=x+d*cos(theDir);
1187 indexJ=y+d*sin(theDir);
1188 if(indexJ<0||indexJ>=input.size())
1190 if(indexI<0||indexI>=input[indexJ].size())
1192 if(input[indexJ][indexI]<currentValue-d*pixelSize/tan(DEG2RAD(sza))){
1193 output[indexJ][indexI]=shadowFlag;
1198 progress/=output.nRows();
1199 pfnProgress(progress,pszMessage,pProgressArg);
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;
1208 pfnProgress(progress,pszMessage,pProgressArg);
1210 int nRow=theBuffer.size();
1212 int nCol=theBuffer[0].size();
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];
1228 int nsize=theBuffer.size()*theBuffer[0].size();
1230 gsl_wavelet_workspace *work;
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];
1242 progress=(1.0+irow);
1243 progress/=theBuffer.nRows();
1244 pfnProgress(progress,pszMessage,pProgressArg);
1246 gsl_wavelet_free (w);
1247 gsl_wavelet_workspace_free (work);
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;
1255 pfnProgress(progress,pszMessage,pProgressArg);
1257 int nRow=theBuffer.size();
1259 int nCol=theBuffer[0].size();
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]);
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];
1276 int nsize=theBuffer.size()*theBuffer[0].size();
1278 gsl_wavelet_workspace *work;
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];
1290 progress=(1.0+irow);
1291 progress/=theBuffer.nRows();
1292 pfnProgress(progress,pszMessage,pProgressArg);
1294 gsl_wavelet_free (w);
1295 gsl_wavelet_workspace_free (work);
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;
1303 pfnProgress(progress,pszMessage,pProgressArg);
1305 int nRow=theBuffer.size();
1307 int nCol=theBuffer[0].size();
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];
1325 int nsize=theBuffer.size()*theBuffer[0].size();
1327 gsl_wavelet_workspace *work;
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]);
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++)
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];
1348 progress=(1.0+irow);
1349 progress/=theBuffer.nRows();
1350 pfnProgress(progress,pszMessage,pProgressArg);
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());
1358 gsl_wavelet_free (w);
1359 gsl_wavelet_workspace_free (work);