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};
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;
176 std::vector<short> m_class;
178 std::vector<double> m_noDataValues;
179 std::vector<double> m_threshold;
183 template<
class T1,
class T2>
void Filter2d::smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dim)
185 smooth(inputVector,outputVector,dim,dim);
188 template<
class T1,
class T2>
void Filter2d::smooth(
const Vector2d<T1>& inputVector,
Vector2d<T2>& outputVector,
int dimX,
int dimY)
191 for(
int j=0;j<dimY;++j){
192 m_taps[j].resize(dimX);
193 for(
int i=0;i<dimX;++i)
194 m_taps[j][i]=1.0/dimX/dimY;
196 filter(inputVector,outputVector);
201 outputVector.resize(inputVector.size());
202 int dimX=m_taps[0].size();
203 int dimY=m_taps.size();
205 std::vector<T2> outBuffer(inputVector[0].size());
210 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
211 inBuffer[indexJ]=inputVector[abs(j)];
215 for(
int y=0;y<inputVector.size();++y){
218 inBuffer.erase(inBuffer.begin());
220 if(y+dimY/2<inputVector.size()){
222 inBuffer.push_back(inputVector[y+dimY/2]);
225 int over=y+dimY/2-inputVector.nRows();
226 int index=(inBuffer.size()-1)-over;
228 assert(index<inBuffer.size());
229 inBuffer.push_back(inBuffer[index]);
232 for(
int x=0;x<inputVector.nCols();++x){
234 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
235 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
241 else if(x>=inputVector.nCols()-(dimX-1)/2)
244 indexJ=(dimY-1)/2+abs(j);
245 else if(y>=inputVector.nRows()-(dimY-1)/2)
246 indexJ=(dimY-1)/2-abs(j);
247 outBuffer[x]+=(m_taps[(dimY-1)/2+j][(dimX-1)/2+i]*inBuffer[indexJ][indexI]);
252 outputVector[y]=outBuffer;
256 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 const char* pszMessage;
259 void* pProgressArg=NULL;
260 GDALProgressFunc pfnProgress=GDALTermProgress;
262 pfnProgress(progress,pszMessage,pProgressArg);
264 double noDataValue=0;
265 if(m_noDataValues.size())
266 noDataValue=m_noDataValues[0];
272 outputVector.resize((inputVector.size()+down-1)/down);
274 std::vector<T2> outBuffer((inputVector[0].size()+down-1)/down);
278 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
279 inBuffer[indexJ]=inputVector[abs(j)];
282 for(
int y=0;y<inputVector.size();++y){
285 inBuffer.erase(inBuffer.begin());
287 if(y+dimY/2<inputVector.size())
288 inBuffer.push_back(inputVector[y+dimY/2]);
290 int over=y+dimY/2-inputVector.size();
291 int index=(inBuffer.size()-1)-over;
293 assert(index<inBuffer.size());
294 inBuffer.push_back(inBuffer[index]);
297 if((y+1+down/2)%down)
299 for(
int x=0;x<inputVector[0].size();++x){
300 if((x+1+down/2)%down)
303 std::vector<double> windowBuffer;
304 std::map<int,int> occurrence;
305 int centre=dimX*(dimY-1)/2+(dimX-1)/2;
306 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
307 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
312 else if(indexI>=inputVector[0].size())
313 indexI=inputVector[0].size()-i;
316 else if(y+j>=inputVector.size())
317 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
321 for(
int imask=0;imask<m_noDataValues.size();++imask){
322 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
328 std::vector<short>::const_iterator vit=m_class.begin();
331 ++occurrence[inBuffer[indexJ][indexI]];
333 while(vit!=m_class.end()){
334 if(inBuffer[indexJ][indexI]==*(vit++))
335 ++occurrence[inBuffer[indexJ][indexI]];
338 windowBuffer.push_back(inBuffer[indexJ][indexI]);
342 switch(getFilterType(method)){
343 case(filter2d::median):
344 if(windowBuffer.empty())
345 outBuffer[x/down]=noDataValue;
347 outBuffer[x/down]=stat.median(windowBuffer);
349 case(filter2d::var):{
350 if(windowBuffer.empty())
351 outBuffer[x/down]=noDataValue;
353 outBuffer[x/down]=stat.var(windowBuffer);
356 case(filter2d::stdev):{
357 if(windowBuffer.empty())
358 outBuffer[x/down]=noDataValue;
360 outBuffer[x/down]=sqrt(stat.var(windowBuffer));
363 case(filter2d::mean):{
364 if(windowBuffer.empty())
365 outBuffer[x/down]=noDataValue;
367 outBuffer[x/down]=stat.mean(windowBuffer);
370 case(filter2d::min):{
371 if(windowBuffer.empty())
372 outBuffer[x/down]=noDataValue;
374 outBuffer[x/down]=stat.mymin(windowBuffer);
377 case(filter2d::ismin):{
378 if(windowBuffer.empty())
379 outBuffer[x/down]=noDataValue;
381 outBuffer[x/down]=(stat.mymin(windowBuffer)==windowBuffer[centre])? 1:0;
384 case(filter2d::minmax):{
387 if(windowBuffer.empty())
388 outBuffer[x/down]=noDataValue;
390 stat.minmax(windowBuffer,windowBuffer.begin(),windowBuffer.end(),min,max);
394 outBuffer[x/down]=windowBuffer[centre];
398 case(filter2d::max):{
399 if(windowBuffer.empty())
400 outBuffer[x/down]=noDataValue;
402 outBuffer[x/down]=stat.mymax(windowBuffer);
405 case(filter2d::ismax):{
406 if(windowBuffer.empty())
407 outBuffer[x/down]=noDataValue;
409 outBuffer[x/down]=(stat.mymax(windowBuffer)==windowBuffer[centre])? 1:0;
412 case(filter2d::order):{
413 if(windowBuffer.empty())
414 outBuffer[x/down]=noDataValue;
417 double ubound=dimX*dimY;
418 double theMin=stat.mymin(windowBuffer);
419 double theMax=stat.mymax(windowBuffer);
420 double scale=(ubound-lbound)/(theMax-theMin);
421 outBuffer[x/down]=
static_cast<short>(scale*(windowBuffer[centre]-theMin)+lbound);
425 case(filter2d::sum):{
426 outBuffer[x/down]=stat.sum(windowBuffer);
429 case(filter2d::percentile):{
430 assert(m_threshold.size());
431 outBuffer[x/down]=stat.percentile(windowBuffer,windowBuffer.begin(),windowBuffer.end(),m_threshold[0]);
434 case(filter2d::homog):
435 if(occurrence.size()==1)
436 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
438 outBuffer[x/down]=noDataValue;
440 case(filter2d::heterog):{
441 for(std::vector<double>::const_iterator wit=windowBuffer.begin();wit!=windowBuffer.end();++wit){
442 if(wit==windowBuffer.begin()+windowBuffer.size()/2)
444 else if(*wit!=inBuffer[(dimY-1)/2][x])
446 else if(*wit==inBuffer[(dimY-1)/2][x]){
447 outBuffer[x/down]=noDataValue;
453 case(filter2d::density):{
454 if(windowBuffer.size()){
455 std::vector<short>::const_iterator vit=m_class.begin();
456 while(vit!=m_class.end())
457 outBuffer[x/down]+=100.0*occurrence[*(vit++)]/windowBuffer.size();
460 outBuffer[x/down]=noDataValue;
463 case(filter2d::countid):{
464 if(windowBuffer.size())
465 outBuffer[x/down]=occurrence.size();
467 outBuffer[x/down]=noDataValue;
470 case(filter2d::mode):{
471 if(occurrence.size()){
472 std::map<int,int>::const_iterator maxit=occurrence.begin();
473 for(std::map<int,int>::const_iterator mit=occurrence.begin();mit!=occurrence.end();++mit){
474 if(mit->second>maxit->second)
477 if(occurrence[inBuffer[(dimY-1)/2][x]]<maxit->second)
478 outBuffer[x/down]=maxit->first;
480 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
483 outBuffer[x/down]=noDataValue;
486 case(filter2d::threshold):{
487 assert(m_class.size()==m_threshold.size());
488 if(windowBuffer.size()){
489 outBuffer[x/down]=inBuffer[(dimY-1)/2][x];
490 for(
int iclass=0;iclass<m_class.size();++iclass){
491 if(100.0*(occurrence[m_class[iclass]])/windowBuffer.size()>m_threshold[iclass])
492 outBuffer[x/down]=m_class[iclass];
496 outBuffer[x/down]=noDataValue;
499 case(filter2d::scramble):{
500 if(windowBuffer.size()){
501 int randomIndex=std::rand()%windowBuffer.size();
502 if(randomIndex>=windowBuffer.size())
503 outBuffer[x/down]=windowBuffer.back();
504 else if(randomIndex<0)
505 outBuffer[x/down]=windowBuffer[0];
507 outBuffer[x/down]=windowBuffer[randomIndex];
510 outBuffer[x/down]=noDataValue;
513 case(filter2d::mixed):{
514 enum MixType { BF=11, CF=12, MF=13, NF=20, W=30 };
515 double nBF=occurrence[BF];
516 double nCF=occurrence[CF];
517 double nMF=occurrence[MF];
518 double nNF=occurrence[NF];
519 double nW=occurrence[W];
520 if(windowBuffer.size()){
521 if((nBF+nCF+nMF)&&(nBF+nCF+nMF>=nNF+nW)){
522 if(nBF/(nBF+nCF)>=0.75)
523 outBuffer[x/down]=BF;
524 else if(nCF/(nBF+nCF)>=0.75)
525 outBuffer[x/down]=CF;
527 outBuffer[x/down]=MF;
533 outBuffer[x/down]=NF;
537 outBuffer[x/down]=inBuffer[indexJ][indexI];
544 progress=(1.0+y/down);
545 progress+=(outputVector.size());
546 progress/=outputVector.size();
547 pfnProgress(progress,pszMessage,pProgressArg);
549 outputVector[y/down]=outBuffer;
560 template<
class T>
void Filter2d::shift(
const Vector2d<T>& input,
Vector2d<T>& output,
double offsetX,
double offsetY,
double randomSigma, RESAMPLE resample,
bool verbose){
561 output.resize(input.nRows(),input.nCols());
562 const gsl_rng_type *rangenType;
565 rangenType=gsl_rng_default;
566 rangen=gsl_rng_alloc(rangenType);
567 long seed=time(NULL)*getpid();
568 gsl_rng_set(rangen,seed);
569 const char* pszMessage;
570 void* pProgressArg=NULL;
571 GDALProgressFunc pfnProgress=GDALTermProgress;
573 pfnProgress(progress,pszMessage,pProgressArg);
574 for(
int j=0;j<input.nRows();++j){
575 for(
int i=0;i<input.nCols();++i){
580 randomX=gsl_ran_gaussian(rangen,randomSigma);
581 randomY=gsl_ran_gaussian(rangen,randomSigma);
583 double readCol=i+offsetX+randomX;
584 double readRow=j+offsetY+randomY;
587 if(readRow>input.nRows()-1)
588 readRow=input.nRows()-1;
591 if(readCol>input.nCols()-1)
592 readCol=input.nCols()-1;
595 double lowerRow=readRow-0.5;
596 double upperRow=readRow+0.5;
597 lowerRow=
static_cast<int>(lowerRow);
598 upperRow=
static_cast<int>(upperRow);
599 double lowerCol=readCol-0.5;
600 double upperCol=readCol+0.5;
601 lowerCol=
static_cast<int>(lowerCol);
602 upperCol=
static_cast<int>(upperCol);
604 assert(lowerRow<input.nRows());
606 assert(lowerCol<input.nCols());
608 assert(upperRow<input.nRows());
610 if(upperCol>=input.nCols()){
611 std::cout <<
"upperCol: " << upperCol << std::endl;
612 std::cout <<
"readCol: " << readCol << std::endl;
613 std::cout <<
"readCol+0.5: " << readCol+0.5 << std::endl;
614 std::cout <<
"static_cast<int>(readCol+0.5): " <<
static_cast<int>(readCol+0.5) << std::endl;
616 assert(upperCol<input.nCols());
617 double c00=input[lowerRow][lowerCol];
618 double c11=input[upperRow][upperCol];
619 double c01=input[lowerRow][upperCol];
620 double c10=input[upperRow][lowerCol];
621 double a=(upperCol-readCol)*c00+(readCol-lowerCol)*c01;
622 double b=(upperCol-readCol)*c10+(readCol-lowerCol)*c11;
623 theValue=(upperRow-readRow)*a+(readRow-lowerRow)*b;
627 theValue=input[
static_cast<int>(readRow)][static_cast<int>(readCol)];
631 assert(j<output.nRows());
633 assert(i<output.nCols());
634 output[j][i]=theValue;
637 progress/=output.nRows();
638 pfnProgress(progress,pszMessage,pProgressArg);
640 gsl_rng_free(rangen);
643 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)
645 const char* pszMessage;
646 void* pProgressArg=NULL;
647 GDALProgressFunc pfnProgress=GDALTermProgress;
649 pfnProgress(progress,pszMessage,pProgressArg);
651 double noDataValue=0;
652 if(m_noDataValues.size())
653 noDataValue=m_noDataValues[0];
655 unsigned long int nchange=0;
661 output.resize(input.nRows(),input.nCols());
665 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
666 for(
int i=0;i<input.nCols();++i)
667 inBuffer[indexJ][i]=input[abs(j)][i];
670 for(
int y=0;y<input.nRows();++y){
673 inBuffer.erase(inBuffer.begin());
675 if(y+dimY/2<input.nRows()){
677 inBuffer.push_back(inBuffer.back());
678 for(
int i=0;i<input.nCols();++i)
679 inBuffer[inBuffer.size()-1][i]=input[y+dimY/2][i];
682 int over=y+dimY/2-input.nRows();
683 int index=(inBuffer.size()-1)-over;
685 assert(index<inBuffer.size());
686 inBuffer.push_back(inBuffer[index]);
689 for(
int x=0;x<input.nCols();++x){
691 double currentValue=inBuffer[(dimY-1)/2][x];
692 std::vector<double> statBuffer;
693 bool currentMasked=
false;
694 for(
int imask=0;imask<m_noDataValues.size();++imask){
695 if(currentValue==m_noDataValues[imask]){
700 output[y][x]=currentValue;
702 output[y][x]=currentValue;
705 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
706 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
708 if(disc&&(d2>(dimX/2)*(dimY/2)))
714 else if(indexI>=input.nCols())
715 indexI=input.nCols()-i;
718 else if(y+j>=input.nRows())
719 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
722 if(inBuffer[indexJ][indexI]==noDataValue)
725 for(
int imask=0;imask<m_noDataValues.size();++imask){
726 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
733 for(
int iclass=0;iclass<m_class.size();++iclass){
734 if(inBuffer[indexJ][indexI]==m_class[iclass]){
740 statBuffer.push_back(binValue);
742 statBuffer.push_back(inBuffer[indexJ][indexI]);
746 if(statBuffer.size()){
747 switch(getFilterType(method)){
748 case(filter2d::dilate):
749 if(output[y][x]<stat.mymax(statBuffer)-hThreshold){
750 output[y][x]=stat.mymax(statBuffer);
754 case(filter2d::erode):
755 if(output[y][x]>stat.mymin(statBuffer)+hThreshold){
756 output[y][x]=stat.mymin(statBuffer);
761 std::ostringstream ess;
762 ess <<
"Error: morphology method " << method <<
" not supported, choose " << filter2d::dilate <<
" (dilate) or " << filter2d::erode <<
" (erode)" << std::endl;
766 if(output[y][x]&&m_class.size())
767 output[y][x]=m_class[0];
774 output[y][x]=noDataValue;
778 progress/=output.nRows();
779 pfnProgress(progress,pszMessage,pProgressArg);
784 template<
class T>
unsigned long int Filter2d::dsm2dtm_nwse(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
786 const char* pszMessage;
787 void* pProgressArg=NULL;
788 GDALProgressFunc pfnProgress=GDALTermProgress;
790 pfnProgress(progress,pszMessage,pProgressArg);
793 double noDataValue=0;
794 if(m_noDataValues.size())
795 noDataValue=m_noDataValues[0];
797 unsigned long int nchange=0;
804 if(outputMask.size()!=inputDSM.nRows())
805 outputMask.resize(inputDSM.nRows());
809 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
810 for(
int i=0;i<inputDSM.nCols();++i)
811 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
814 for(
int y=0;y<tmpDSM.nRows();++y){
817 inBuffer.erase(inBuffer.begin());
819 if(y+dimY/2<tmpDSM.nRows()){
821 inBuffer.push_back(inBuffer.back());
822 for(
int i=0;i<tmpDSM.nCols();++i)
823 inBuffer[inBuffer.size()-1][i]=tmpDSM[y+dimY/2][i];
826 int over=y+dimY/2-tmpDSM.nRows();
827 int index=(inBuffer.size()-1)-over;
829 assert(index<inBuffer.size());
830 inBuffer.push_back(inBuffer[index]);
833 for(
int x=0;x<tmpDSM.nCols();++x){
834 double centerValue=inBuffer[(dimY-1)/2][x];
836 std::vector<T> neighbors;
837 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
838 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
843 else if(indexI>=tmpDSM.nCols())
844 indexI=tmpDSM.nCols()-i;
847 else if(y+j>=tmpDSM.nRows())
848 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
851 double difference=(centerValue-inBuffer[indexJ][indexI]);
853 neighbors.push_back(inBuffer[indexJ][indexI]);
854 if(difference>hThreshold)
865 sort(neighbors.begin(),neighbors.end());
866 assert(neighbors.size()>1);
867 inBuffer[(dimY-1)/2][x]=neighbors[1];
872 progress/=outputMask.nRows();
873 pfnProgress(progress,pszMessage,pProgressArg);
878 template<
class T>
unsigned long int Filter2d::dsm2dtm_nesw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
880 const char* pszMessage;
881 void* pProgressArg=NULL;
882 GDALProgressFunc pfnProgress=GDALTermProgress;
884 pfnProgress(progress,pszMessage,pProgressArg);
887 double noDataValue=0;
888 if(m_noDataValues.size())
889 noDataValue=m_noDataValues[0];
891 unsigned long int nchange=0;
898 if(outputMask.size()!=inputDSM.nRows())
899 outputMask.resize(inputDSM.nRows());
903 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
904 for(
int i=0;i<inputDSM.nCols();++i)
905 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
908 for(
int y=0;y<tmpDSM.nRows();++y){
911 inBuffer.erase(inBuffer.begin());
913 if(y+dimY/2<tmpDSM.nRows()){
915 inBuffer.push_back(inBuffer.back());
916 for(
int i=0;i<tmpDSM.nCols();++i)
917 inBuffer[inBuffer.size()-1][i]=tmpDSM[y+dimY/2][i];
920 int over=y+dimY/2-tmpDSM.nRows();
921 int index=(inBuffer.size()-1)-over;
923 assert(index<inBuffer.size());
924 inBuffer.push_back(inBuffer[index]);
927 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
928 double centerValue=inBuffer[(dimY-1)/2][x];
930 std::vector<T> neighbors;
931 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
932 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
937 else if(indexI>=tmpDSM.nCols())
938 indexI=tmpDSM.nCols()-i;
941 else if(y+j>=tmpDSM.nRows())
942 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
945 double difference=(centerValue-inBuffer[indexJ][indexI]);
947 neighbors.push_back(inBuffer[indexJ][indexI]);
948 if(difference>hThreshold)
959 sort(neighbors.begin(),neighbors.end());
960 assert(neighbors.size()>1);
961 inBuffer[(dimY-1)/2][x]=neighbors[1];
966 progress/=outputMask.nRows();
967 pfnProgress(progress,pszMessage,pProgressArg);
972 template<
class T>
unsigned long int Filter2d::dsm2dtm_senw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
974 const char* pszMessage;
975 void* pProgressArg=NULL;
976 GDALProgressFunc pfnProgress=GDALTermProgress;
978 pfnProgress(progress,pszMessage,pProgressArg);
981 double noDataValue=0;
982 if(m_noDataValues.size())
983 noDataValue=m_noDataValues[0];
985 unsigned long int nchange=0;
992 if(outputMask.size()!=inputDSM.nRows())
993 outputMask.resize(inputDSM.nRows());
995 int indexJ=inputDSM.nRows()-1;
997 for(
int j=inputDSM.nRows()-dimY/2;j<inputDSM.nRows();--j){
998 for(
int i=0;i<inputDSM.nCols();++i)
999 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
1002 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1003 if(y<tmpDSM.nRows()-1){
1005 inBuffer.erase(inBuffer.end()-1);
1009 inBuffer.insert(inBuffer.begin(),inBuffer.back());
1010 for(
int i=0;i<tmpDSM.nCols();++i)
1011 inBuffer[0][i]=tmpDSM[y-dimY/2][i];
1014 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1017 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
1018 double centerValue=inBuffer[(dimY-1)/2][x];
1020 std::vector<T> neighbors;
1021 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
1022 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
1027 else if(indexI>=tmpDSM.nCols())
1028 indexI=tmpDSM.nCols()-i;
1031 else if(y+j>=tmpDSM.nRows())
1032 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1034 indexJ=(dimY-1)/2+j;
1035 double difference=(centerValue-inBuffer[indexJ][indexI]);
1037 neighbors.push_back(inBuffer[indexJ][indexI]);
1038 if(difference>hThreshold)
1042 if(nmasked<=nlimit){
1049 sort(neighbors.begin(),neighbors.end());
1050 assert(neighbors.size()>1);
1051 inBuffer[(dimY-1)/2][x]=neighbors[1];
1056 progress/=outputMask.nRows();
1057 pfnProgress(progress,pszMessage,pProgressArg);
1062 template<
class T>
unsigned long int Filter2d::dsm2dtm_swne(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
1064 const char* pszMessage;
1065 void* pProgressArg=NULL;
1066 GDALProgressFunc pfnProgress=GDALTermProgress;
1068 pfnProgress(progress,pszMessage,pProgressArg);
1071 double noDataValue=0;
1072 if(m_noDataValues.size())
1073 noDataValue=m_noDataValues[0];
1075 unsigned long int nchange=0;
1082 if(outputMask.size()!=inputDSM.nRows())
1083 outputMask.resize(inputDSM.nRows());
1087 for(
int j=inputDSM.nRows()-dimY/2;j<inputDSM.nRows();--j){
1088 for(
int i=0;i<inputDSM.nCols();++i)
1089 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
1092 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1093 if(y<tmpDSM.nRows()-1){
1095 inBuffer.erase(inBuffer.end()-1);
1099 inBuffer.insert(inBuffer.begin(),inBuffer.back());
1100 for(
int i=0;i<tmpDSM.nCols();++i)
1101 inBuffer[0][i]=tmpDSM[y-dimY/2][i];
1104 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1107 for(
int x=0;x<tmpDSM.nCols();++x){
1108 double centerValue=inBuffer[(dimY-1)/2][x];
1110 std::vector<T> neighbors;
1111 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
1112 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
1117 else if(indexI>=tmpDSM.nCols())
1118 indexI=tmpDSM.nCols()-i;
1121 else if(y+j>=tmpDSM.nRows())
1122 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1124 indexJ=(dimY-1)/2+j;
1125 double difference=(centerValue-inBuffer[indexJ][indexI]);
1127 neighbors.push_back(inBuffer[indexJ][indexI]);
1128 if(difference>hThreshold)
1132 if(nmasked<=nlimit){
1139 sort(neighbors.begin(),neighbors.end());
1140 assert(neighbors.size()>1);
1141 inBuffer[(dimY-1)/2][x]=neighbors[1];
1146 progress/=outputMask.nRows();
1147 pfnProgress(progress,pszMessage,pProgressArg);
1152 template<
class T>
void Filter2d::shadowDsm(
const Vector2d<T>& input,
Vector2d<T>& output,
double sza,
double saa,
double pixelSize,
short shadowFlag)
1154 unsigned int ncols=input.nCols();
1156 output.resize(input.nRows(),ncols);
1163 const char* pszMessage;
1164 void* pProgressArg=NULL;
1165 GDALProgressFunc pfnProgress=GDALTermProgress;
1167 pfnProgress(progress,pszMessage,pProgressArg);
1168 for(
int y=0;y<input.nRows();++y){
1169 for(
int x=0;x<input.nCols();++x){
1170 double currentValue=input[y][x];
1171 int theDist=
static_cast<int>(sqrt((currentValue*tan(DEG2RAD(sza))/pixelSize)*(currentValue*tan(DEG2RAD(sza))/pixelSize)));
1172 double theDir=DEG2RAD(saa)+PI/2.0;
1175 for(
int d=0;d<theDist;++d){
1176 indexI=x+d*cos(theDir);
1177 indexJ=y+d*sin(theDir);
1178 if(indexJ<0||indexJ>=input.size())
1180 if(indexI<0||indexI>=input[indexJ].size())
1182 if(input[indexJ][indexI]<currentValue-d*pixelSize/tan(DEG2RAD(sza))){
1183 output[indexJ][indexI]=shadowFlag;
1188 progress/=output.nRows();
1189 pfnProgress(progress,pszMessage,pProgressArg);
1193 template<
class T>
void Filter2d::dwtForward(
Vector2d<T>& theBuffer,
const std::string& wavelet_type,
int family){
1194 const char* pszMessage;
1195 void* pProgressArg=NULL;
1196 GDALProgressFunc pfnProgress=GDALTermProgress;
1198 pfnProgress(progress,pszMessage,pProgressArg);
1200 int nRow=theBuffer.size();
1202 int nCol=theBuffer[0].size();
1205 while(theBuffer.size()&(theBuffer.size()-1))
1206 theBuffer.push_back(theBuffer.back());
1207 for(
int irow=0;irow<theBuffer.size();++irow)
1208 while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1209 theBuffer[irow].push_back(theBuffer[irow].back());
1210 std::vector<double> vdata(theBuffer.size()*theBuffer[0].size());
1211 double* data=&(vdata[0]);
1212 for(
int irow=0;irow<theBuffer.size();++irow){
1213 for(
int icol=0;icol<theBuffer[0].size();++icol){
1214 int index=irow*theBuffer[0].size()+icol;
1215 data[index]=theBuffer[irow][icol];
1218 int nsize=theBuffer.size()*theBuffer[0].size();
1220 gsl_wavelet_workspace *work;
1222 w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1223 work=gsl_wavelet_workspace_alloc(nsize);
1224 gsl_wavelet2d_nstransform_forward (w, data, theBuffer.size(), theBuffer.size(),theBuffer[0].size(), work);
1225 theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1226 for(
int irow=0;irow<theBuffer.size();++irow){
1227 theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1228 for(
int icol=0;icol<theBuffer[irow].size();++icol){
1229 int index=irow*theBuffer[irow].size()+icol;
1230 theBuffer[irow][icol]=data[index];
1232 progress=(1.0+irow);
1233 progress/=theBuffer.nRows();
1234 pfnProgress(progress,pszMessage,pProgressArg);
1236 gsl_wavelet_free (w);
1237 gsl_wavelet_workspace_free (work);
1240 template<
class T>
void Filter2d::dwtInverse(
Vector2d<T>& theBuffer,
const std::string& wavelet_type,
int family){
1241 const char* pszMessage;
1242 void* pProgressArg=NULL;
1243 GDALProgressFunc pfnProgress=GDALTermProgress;
1245 pfnProgress(progress,pszMessage,pProgressArg);
1247 int nRow=theBuffer.size();
1249 int nCol=theBuffer[0].size();
1252 while(theBuffer.size()&(theBuffer.size()-1))
1253 theBuffer.push_back(theBuffer.back());
1254 for(
int irow=0;irow<theBuffer.size();++irow)
1255 while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1256 theBuffer[irow].push_back(theBuffer[irow].back());
1257 std::vector<double> vdata(theBuffer.size()*theBuffer[0].size());
1258 double* data=&(vdata[0]);
1260 for(
int irow=0;irow<theBuffer.size();++irow){
1261 for(
int icol=0;icol<theBuffer[0].size();++icol){
1262 int index=irow*theBuffer[0].size()+icol;
1263 data[index]=theBuffer[irow][icol];
1266 int nsize=theBuffer.size()*theBuffer[0].size();
1268 gsl_wavelet_workspace *work;
1270 w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1271 work=gsl_wavelet_workspace_alloc(nsize);
1272 gsl_wavelet2d_nstransform_inverse (w, data, theBuffer.size(), theBuffer.size(),theBuffer[0].size(), work);
1273 theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1274 for(
int irow=0;irow<theBuffer.size();++irow){
1275 theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1276 for(
int icol=0;icol<theBuffer[irow].size();++icol){
1277 int index=irow*theBuffer[irow].size()+icol;
1278 theBuffer[irow][icol]=data[index];
1280 progress=(1.0+irow);
1281 progress/=theBuffer.nRows();
1282 pfnProgress(progress,pszMessage,pProgressArg);
1284 gsl_wavelet_free (w);
1285 gsl_wavelet_workspace_free (work);
1288 template<
class T>
void Filter2d::dwtCut(
Vector2d<T>& theBuffer,
const std::string& wavelet_type,
int family,
double cut){
1289 const char* pszMessage;
1290 void* pProgressArg=NULL;
1291 GDALProgressFunc pfnProgress=GDALTermProgress;
1293 pfnProgress(progress,pszMessage,pProgressArg);
1295 int nRow=theBuffer.size();
1297 int nCol=theBuffer[0].size();
1300 while(theBuffer.size()&(theBuffer.size()-1))
1301 theBuffer.push_back(theBuffer.back());
1302 for(
int irow=0;irow<theBuffer.size();++irow)
1303 while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1304 theBuffer[irow].push_back(theBuffer[irow].back());
1305 double* data=
new double[theBuffer.size()*theBuffer[0].size()];
1306 double* abscoeff=
new double[theBuffer.size()*theBuffer[0].size()];
1307 size_t* p=
new size_t[theBuffer.size()*theBuffer[0].size()];
1308 for(
int irow=0;irow<theBuffer.size();++irow){
1309 for(
int icol=0;icol<theBuffer[0].size();++icol){
1310 int index=irow*theBuffer[0].size()+icol;
1311 assert(index<theBuffer.size()*theBuffer[0].size());
1312 data[index]=theBuffer[irow][icol];
1315 int nsize=theBuffer.size()*theBuffer[0].size();
1317 gsl_wavelet_workspace *work;
1319 w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1320 work=gsl_wavelet_workspace_alloc(nsize);
1321 gsl_wavelet2d_nstransform_forward (w, data, theBuffer.size(), theBuffer[0].size(),theBuffer[0].size(), work);
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 abscoeff[index]=fabs(data[index]);
1328 int nc=(100-cut)/100.0*nsize;
1329 gsl_sort_index(p,abscoeff,1,nsize);
1330 for(
int i=0;(i+nc)<nsize;i++)
1332 gsl_wavelet2d_nstransform_inverse (w, data, theBuffer.size(), theBuffer[0].size(),theBuffer[0].size(), work);
1333 for(
int irow=0;irow<theBuffer.size();++irow){
1334 for(
int icol=0;icol<theBuffer[irow].size();++icol){
1335 int index=irow*theBuffer[irow].size()+icol;
1336 theBuffer[irow][icol]=data[index];
1338 progress=(1.0+irow);
1339 progress/=theBuffer.nRows();
1340 pfnProgress(progress,pszMessage,pProgressArg);
1342 theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1343 for(
int irow=0;irow<theBuffer.size();++irow)
1344 theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1348 gsl_wavelet_free (w);
1349 gsl_wavelet_workspace_free (work);