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 outBuffer[x/down]=windowBuffer[randomIndex];
505 outBuffer[x/down]=noDataValue;
508 case(filter2d::mixed):{
509 enum MixType { BF=11, CF=12, MF=13, NF=20, W=30 };
510 double nBF=occurrence[BF];
511 double nCF=occurrence[CF];
512 double nMF=occurrence[MF];
513 double nNF=occurrence[NF];
514 double nW=occurrence[W];
515 if(windowBuffer.size()){
516 if((nBF+nCF+nMF)&&(nBF+nCF+nMF>=nNF+nW)){
517 if(nBF/(nBF+nCF)>=0.75)
518 outBuffer[x/down]=BF;
519 else if(nCF/(nBF+nCF)>=0.75)
520 outBuffer[x/down]=CF;
522 outBuffer[x/down]=MF;
528 outBuffer[x/down]=NF;
532 outBuffer[x/down]=inBuffer[indexJ][indexI];
539 progress=(1.0+y/down);
540 progress+=(outputVector.size());
541 progress/=outputVector.size();
542 pfnProgress(progress,pszMessage,pProgressArg);
544 outputVector[y/down]=outBuffer;
555 template<
class T>
void Filter2d::shift(
const Vector2d<T>& input,
Vector2d<T>& output,
double offsetX,
double offsetY,
double randomSigma, RESAMPLE resample,
bool verbose){
556 output.resize(input.nRows(),input.nCols());
557 const gsl_rng_type *rangenType;
560 rangenType=gsl_rng_default;
561 rangen=gsl_rng_alloc(rangenType);
562 long seed=time(NULL)*getpid();
563 gsl_rng_set(rangen,seed);
564 const char* pszMessage;
565 void* pProgressArg=NULL;
566 GDALProgressFunc pfnProgress=GDALTermProgress;
568 pfnProgress(progress,pszMessage,pProgressArg);
569 for(
int j=0;j<input.nRows();++j){
570 for(
int i=0;i<input.nCols();++i){
575 randomX=gsl_ran_gaussian(rangen,randomSigma);
576 randomY=gsl_ran_gaussian(rangen,randomSigma);
578 double readCol=i+offsetX+randomX;
579 double readRow=j+offsetY+randomY;
582 if(readRow>input.nRows()-1)
583 readRow=input.nRows()-1;
586 if(readCol>input.nCols()-1)
587 readCol=input.nCols()-1;
590 double lowerRow=readRow-0.5;
591 double upperRow=readRow+0.5;
592 lowerRow=
static_cast<int>(lowerRow);
593 upperRow=
static_cast<int>(upperRow);
594 double lowerCol=readCol-0.5;
595 double upperCol=readCol+0.5;
596 lowerCol=
static_cast<int>(lowerCol);
597 upperCol=
static_cast<int>(upperCol);
599 assert(lowerRow<input.nRows());
601 assert(lowerCol<input.nCols());
603 assert(upperRow<input.nRows());
605 if(upperCol>=input.nCols()){
606 std::cout <<
"upperCol: " << upperCol << std::endl;
607 std::cout <<
"readCol: " << readCol << std::endl;
608 std::cout <<
"readCol+0.5: " << readCol+0.5 << std::endl;
609 std::cout <<
"static_cast<int>(readCol+0.5): " <<
static_cast<int>(readCol+0.5) << std::endl;
611 assert(upperCol<input.nCols());
612 double c00=input[lowerRow][lowerCol];
613 double c11=input[upperRow][upperCol];
614 double c01=input[lowerRow][upperCol];
615 double c10=input[upperRow][lowerCol];
616 double a=(upperCol-readCol)*c00+(readCol-lowerCol)*c01;
617 double b=(upperCol-readCol)*c10+(readCol-lowerCol)*c11;
618 theValue=(upperRow-readRow)*a+(readRow-lowerRow)*b;
622 theValue=input[
static_cast<int>(readRow)][static_cast<int>(readCol)];
626 assert(j<output.nRows());
628 assert(i<output.nCols());
629 output[j][i]=theValue;
632 progress/=output.nRows();
633 pfnProgress(progress,pszMessage,pProgressArg);
635 gsl_rng_free(rangen);
638 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)
640 const char* pszMessage;
641 void* pProgressArg=NULL;
642 GDALProgressFunc pfnProgress=GDALTermProgress;
644 pfnProgress(progress,pszMessage,pProgressArg);
646 double noDataValue=0;
647 if(m_noDataValues.size())
648 noDataValue=m_noDataValues[0];
650 unsigned long int nchange=0;
656 output.resize(input.nRows(),input.nCols());
660 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
661 for(
int i=0;i<input.nCols();++i)
662 inBuffer[indexJ][i]=input[abs(j)][i];
665 for(
int y=0;y<input.nRows();++y){
668 inBuffer.erase(inBuffer.begin());
670 if(y+dimY/2<input.nRows()){
672 inBuffer.push_back(inBuffer.back());
673 for(
int i=0;i<input.nCols();++i)
674 inBuffer[inBuffer.size()-1][i]=input[y+dimY/2][i];
677 int over=y+dimY/2-input.nRows();
678 int index=(inBuffer.size()-1)-over;
680 assert(index<inBuffer.size());
681 inBuffer.push_back(inBuffer[index]);
684 for(
int x=0;x<input.nCols();++x){
686 double currentValue=inBuffer[(dimY-1)/2][x];
687 std::vector<double> statBuffer;
688 bool currentMasked=
false;
689 for(
int imask=0;imask<m_noDataValues.size();++imask){
690 if(currentValue==m_noDataValues[imask]){
695 output[y][x]=currentValue;
697 output[y][x]=currentValue;
700 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
701 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
703 if(disc&&(d2>(dimX/2)*(dimY/2)))
709 else if(indexI>=input.nCols())
710 indexI=input.nCols()-i;
713 else if(y+j>=input.nRows())
714 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
717 if(inBuffer[indexJ][indexI]==noDataValue)
720 for(
int imask=0;imask<m_noDataValues.size();++imask){
721 if(inBuffer[indexJ][indexI]==m_noDataValues[imask]){
728 for(
int iclass=0;iclass<m_class.size();++iclass){
729 if(inBuffer[indexJ][indexI]==m_class[iclass]){
735 statBuffer.push_back(binValue);
737 statBuffer.push_back(inBuffer[indexJ][indexI]);
741 if(statBuffer.size()){
742 switch(getFilterType(method)){
743 case(filter2d::dilate):
744 if(output[y][x]<stat.mymax(statBuffer)-hThreshold){
745 output[y][x]=stat.mymax(statBuffer);
749 case(filter2d::erode):
750 if(output[y][x]>stat.mymin(statBuffer)+hThreshold){
751 output[y][x]=stat.mymin(statBuffer);
756 std::ostringstream ess;
757 ess <<
"Error: morphology method " << method <<
" not supported, choose " << filter2d::dilate <<
" (dilate) or " << filter2d::erode <<
" (erode)" << std::endl;
761 if(output[y][x]&&m_class.size())
762 output[y][x]=m_class[0];
769 output[y][x]=noDataValue;
773 progress/=output.nRows();
774 pfnProgress(progress,pszMessage,pProgressArg);
779 template<
class T>
unsigned long int Filter2d::dsm2dtm_nwse(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
781 const char* pszMessage;
782 void* pProgressArg=NULL;
783 GDALProgressFunc pfnProgress=GDALTermProgress;
785 pfnProgress(progress,pszMessage,pProgressArg);
788 double noDataValue=0;
789 if(m_noDataValues.size())
790 noDataValue=m_noDataValues[0];
792 unsigned long int nchange=0;
799 if(outputMask.size()!=inputDSM.nRows())
800 outputMask.resize(inputDSM.nRows());
804 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
805 for(
int i=0;i<inputDSM.nCols();++i)
806 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
809 for(
int y=0;y<tmpDSM.nRows();++y){
812 inBuffer.erase(inBuffer.begin());
814 if(y+dimY/2<tmpDSM.nRows()){
816 inBuffer.push_back(inBuffer.back());
817 for(
int i=0;i<tmpDSM.nCols();++i)
818 inBuffer[inBuffer.size()-1][i]=tmpDSM[y+dimY/2][i];
821 int over=y+dimY/2-tmpDSM.nRows();
822 int index=(inBuffer.size()-1)-over;
824 assert(index<inBuffer.size());
825 inBuffer.push_back(inBuffer[index]);
828 for(
int x=0;x<tmpDSM.nCols();++x){
829 double centerValue=inBuffer[(dimY-1)/2][x];
831 std::vector<T> neighbors;
832 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
833 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
838 else if(indexI>=tmpDSM.nCols())
839 indexI=tmpDSM.nCols()-i;
842 else if(y+j>=tmpDSM.nRows())
843 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
846 double difference=(centerValue-inBuffer[indexJ][indexI]);
848 neighbors.push_back(inBuffer[indexJ][indexI]);
849 if(difference>hThreshold)
860 sort(neighbors.begin(),neighbors.end());
861 assert(neighbors.size()>1);
862 inBuffer[(dimY-1)/2][x]=neighbors[1];
867 progress/=outputMask.nRows();
868 pfnProgress(progress,pszMessage,pProgressArg);
873 template<
class T>
unsigned long int Filter2d::dsm2dtm_nesw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
875 const char* pszMessage;
876 void* pProgressArg=NULL;
877 GDALProgressFunc pfnProgress=GDALTermProgress;
879 pfnProgress(progress,pszMessage,pProgressArg);
882 double noDataValue=0;
883 if(m_noDataValues.size())
884 noDataValue=m_noDataValues[0];
886 unsigned long int nchange=0;
893 if(outputMask.size()!=inputDSM.nRows())
894 outputMask.resize(inputDSM.nRows());
898 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
899 for(
int i=0;i<inputDSM.nCols();++i)
900 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
903 for(
int y=0;y<tmpDSM.nRows();++y){
906 inBuffer.erase(inBuffer.begin());
908 if(y+dimY/2<tmpDSM.nRows()){
910 inBuffer.push_back(inBuffer.back());
911 for(
int i=0;i<tmpDSM.nCols();++i)
912 inBuffer[inBuffer.size()-1][i]=tmpDSM[y+dimY/2][i];
915 int over=y+dimY/2-tmpDSM.nRows();
916 int index=(inBuffer.size()-1)-over;
918 assert(index<inBuffer.size());
919 inBuffer.push_back(inBuffer[index]);
922 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
923 double centerValue=inBuffer[(dimY-1)/2][x];
925 std::vector<T> neighbors;
926 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
927 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
932 else if(indexI>=tmpDSM.nCols())
933 indexI=tmpDSM.nCols()-i;
936 else if(y+j>=tmpDSM.nRows())
937 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
940 double difference=(centerValue-inBuffer[indexJ][indexI]);
942 neighbors.push_back(inBuffer[indexJ][indexI]);
943 if(difference>hThreshold)
954 sort(neighbors.begin(),neighbors.end());
955 assert(neighbors.size()>1);
956 inBuffer[(dimY-1)/2][x]=neighbors[1];
961 progress/=outputMask.nRows();
962 pfnProgress(progress,pszMessage,pProgressArg);
967 template<
class T>
unsigned long int Filter2d::dsm2dtm_senw(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
969 const char* pszMessage;
970 void* pProgressArg=NULL;
971 GDALProgressFunc pfnProgress=GDALTermProgress;
973 pfnProgress(progress,pszMessage,pProgressArg);
976 double noDataValue=0;
977 if(m_noDataValues.size())
978 noDataValue=m_noDataValues[0];
980 unsigned long int nchange=0;
987 if(outputMask.size()!=inputDSM.nRows())
988 outputMask.resize(inputDSM.nRows());
990 int indexJ=inputDSM.nRows()-1;
992 for(
int j=inputDSM.nRows()-dimY/2;j<inputDSM.nRows();--j){
993 for(
int i=0;i<inputDSM.nCols();++i)
994 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
997 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
998 if(y<tmpDSM.nRows()-1){
1000 inBuffer.erase(inBuffer.end()-1);
1004 inBuffer.insert(inBuffer.begin(),inBuffer.back());
1005 for(
int i=0;i<tmpDSM.nCols();++i)
1006 inBuffer[0][i]=tmpDSM[y-dimY/2][i];
1009 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1012 for(
int x=tmpDSM.nCols()-1;x>=0;--x){
1013 double centerValue=inBuffer[(dimY-1)/2][x];
1015 std::vector<T> neighbors;
1016 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
1017 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
1022 else if(indexI>=tmpDSM.nCols())
1023 indexI=tmpDSM.nCols()-i;
1026 else if(y+j>=tmpDSM.nRows())
1027 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1029 indexJ=(dimY-1)/2+j;
1030 double difference=(centerValue-inBuffer[indexJ][indexI]);
1032 neighbors.push_back(inBuffer[indexJ][indexI]);
1033 if(difference>hThreshold)
1037 if(nmasked<=nlimit){
1044 sort(neighbors.begin(),neighbors.end());
1045 assert(neighbors.size()>1);
1046 inBuffer[(dimY-1)/2][x]=neighbors[1];
1051 progress/=outputMask.nRows();
1052 pfnProgress(progress,pszMessage,pProgressArg);
1057 template<
class T>
unsigned long int Filter2d::dsm2dtm_swne(
const Vector2d<T>& inputDSM,
Vector2d<T>& outputMask,
double hThreshold,
int nlimit,
int dim)
1059 const char* pszMessage;
1060 void* pProgressArg=NULL;
1061 GDALProgressFunc pfnProgress=GDALTermProgress;
1063 pfnProgress(progress,pszMessage,pProgressArg);
1066 double noDataValue=0;
1067 if(m_noDataValues.size())
1068 noDataValue=m_noDataValues[0];
1070 unsigned long int nchange=0;
1077 if(outputMask.size()!=inputDSM.nRows())
1078 outputMask.resize(inputDSM.nRows());
1082 for(
int j=inputDSM.nRows()-dimY/2;j<inputDSM.nRows();--j){
1083 for(
int i=0;i<inputDSM.nCols();++i)
1084 inBuffer[indexJ][i]=tmpDSM[abs(j)][i];
1087 for(
int y=tmpDSM.nRows()-1;y>=0;--y){
1088 if(y<tmpDSM.nRows()-1){
1090 inBuffer.erase(inBuffer.end()-1);
1094 inBuffer.insert(inBuffer.begin(),inBuffer.back());
1095 for(
int i=0;i<tmpDSM.nCols();++i)
1096 inBuffer[0][i]=tmpDSM[y-dimY/2][i];
1099 inBuffer.insert(inBuffer.begin(),inBuffer[abs(y-dimY/2)]);
1102 for(
int x=0;x<tmpDSM.nCols();++x){
1103 double centerValue=inBuffer[(dimY-1)/2][x];
1105 std::vector<T> neighbors;
1106 for(
int j=-(dimY-1)/2;j<=dimY/2;++j){
1107 for(
int i=-(dimX-1)/2;i<=dimX/2;++i){
1112 else if(indexI>=tmpDSM.nCols())
1113 indexI=tmpDSM.nCols()-i;
1116 else if(y+j>=tmpDSM.nRows())
1117 indexJ=(dimY>2) ? (dimY-1)/2-j : 0;
1119 indexJ=(dimY-1)/2+j;
1120 double difference=(centerValue-inBuffer[indexJ][indexI]);
1122 neighbors.push_back(inBuffer[indexJ][indexI]);
1123 if(difference>hThreshold)
1127 if(nmasked<=nlimit){
1134 sort(neighbors.begin(),neighbors.end());
1135 assert(neighbors.size()>1);
1136 inBuffer[(dimY-1)/2][x]=neighbors[1];
1141 progress/=outputMask.nRows();
1142 pfnProgress(progress,pszMessage,pProgressArg);
1147 template<
class T>
void Filter2d::shadowDsm(
const Vector2d<T>& input,
Vector2d<T>& output,
double sza,
double saa,
double pixelSize,
short shadowFlag)
1149 unsigned int ncols=input.nCols();
1151 output.resize(input.nRows(),ncols);
1158 const char* pszMessage;
1159 void* pProgressArg=NULL;
1160 GDALProgressFunc pfnProgress=GDALTermProgress;
1162 pfnProgress(progress,pszMessage,pProgressArg);
1163 for(
int y=0;y<input.nRows();++y){
1164 for(
int x=0;x<input.nCols();++x){
1165 double currentValue=input[y][x];
1166 int theDist=
static_cast<int>(sqrt((currentValue*tan(DEG2RAD(sza))/pixelSize)*(currentValue*tan(DEG2RAD(sza))/pixelSize)));
1167 double theDir=DEG2RAD(saa)+PI/2.0;
1170 for(
int d=0;d<theDist;++d){
1171 indexI=x+d*cos(theDir);
1172 indexJ=y+d*sin(theDir);
1173 if(indexJ<0||indexJ>=input.size())
1175 if(indexI<0||indexI>=input[indexJ].size())
1177 if(input[indexJ][indexI]<currentValue-d*pixelSize/tan(DEG2RAD(sza))){
1178 output[indexJ][indexI]=shadowFlag;
1183 progress/=output.nRows();
1184 pfnProgress(progress,pszMessage,pProgressArg);
1188 template<
class T>
void Filter2d::dwtForward(
Vector2d<T>& theBuffer,
const std::string& wavelet_type,
int family){
1189 const char* pszMessage;
1190 void* pProgressArg=NULL;
1191 GDALProgressFunc pfnProgress=GDALTermProgress;
1193 pfnProgress(progress,pszMessage,pProgressArg);
1195 int nRow=theBuffer.size();
1197 int nCol=theBuffer[0].size();
1200 while(theBuffer.size()&(theBuffer.size()-1))
1201 theBuffer.push_back(theBuffer.back());
1202 for(
int irow=0;irow<theBuffer.size();++irow)
1203 while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1204 theBuffer[irow].push_back(theBuffer[irow].back());
1205 std::vector<double> vdata(theBuffer.size()*theBuffer[0].size());
1206 double* data=&(vdata[0]);
1207 for(
int irow=0;irow<theBuffer.size();++irow){
1208 for(
int icol=0;icol<theBuffer[0].size();++icol){
1209 int index=irow*theBuffer[0].size()+icol;
1210 data[index]=theBuffer[irow][icol];
1213 int nsize=theBuffer.size()*theBuffer[0].size();
1215 gsl_wavelet_workspace *work;
1217 w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1218 work=gsl_wavelet_workspace_alloc(nsize);
1219 gsl_wavelet2d_nstransform_forward (w, data, theBuffer.size(), theBuffer.size(),theBuffer[0].size(), work);
1220 theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1221 for(
int irow=0;irow<theBuffer.size();++irow){
1222 theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1223 for(
int icol=0;icol<theBuffer[irow].size();++icol){
1224 int index=irow*theBuffer[irow].size()+icol;
1225 theBuffer[irow][icol]=data[index];
1227 progress=(1.0+irow);
1228 progress/=theBuffer.nRows();
1229 pfnProgress(progress,pszMessage,pProgressArg);
1231 gsl_wavelet_free (w);
1232 gsl_wavelet_workspace_free (work);
1235 template<
class T>
void Filter2d::dwtInverse(
Vector2d<T>& theBuffer,
const std::string& wavelet_type,
int family){
1236 const char* pszMessage;
1237 void* pProgressArg=NULL;
1238 GDALProgressFunc pfnProgress=GDALTermProgress;
1240 pfnProgress(progress,pszMessage,pProgressArg);
1242 int nRow=theBuffer.size();
1244 int nCol=theBuffer[0].size();
1247 while(theBuffer.size()&(theBuffer.size()-1))
1248 theBuffer.push_back(theBuffer.back());
1249 for(
int irow=0;irow<theBuffer.size();++irow)
1250 while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1251 theBuffer[irow].push_back(theBuffer[irow].back());
1252 std::vector<double> vdata(theBuffer.size()*theBuffer[0].size());
1253 double* data=&(vdata[0]);
1255 for(
int irow=0;irow<theBuffer.size();++irow){
1256 for(
int icol=0;icol<theBuffer[0].size();++icol){
1257 int index=irow*theBuffer[0].size()+icol;
1258 data[index]=theBuffer[irow][icol];
1261 int nsize=theBuffer.size()*theBuffer[0].size();
1263 gsl_wavelet_workspace *work;
1265 w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1266 work=gsl_wavelet_workspace_alloc(nsize);
1267 gsl_wavelet2d_nstransform_inverse (w, data, theBuffer.size(), theBuffer.size(),theBuffer[0].size(), work);
1268 theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1269 for(
int irow=0;irow<theBuffer.size();++irow){
1270 theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1271 for(
int icol=0;icol<theBuffer[irow].size();++icol){
1272 int index=irow*theBuffer[irow].size()+icol;
1273 theBuffer[irow][icol]=data[index];
1275 progress=(1.0+irow);
1276 progress/=theBuffer.nRows();
1277 pfnProgress(progress,pszMessage,pProgressArg);
1279 gsl_wavelet_free (w);
1280 gsl_wavelet_workspace_free (work);
1283 template<
class T>
void Filter2d::dwtCut(
Vector2d<T>& theBuffer,
const std::string& wavelet_type,
int family,
double cut){
1284 const char* pszMessage;
1285 void* pProgressArg=NULL;
1286 GDALProgressFunc pfnProgress=GDALTermProgress;
1288 pfnProgress(progress,pszMessage,pProgressArg);
1290 int nRow=theBuffer.size();
1292 int nCol=theBuffer[0].size();
1295 while(theBuffer.size()&(theBuffer.size()-1))
1296 theBuffer.push_back(theBuffer.back());
1297 for(
int irow=0;irow<theBuffer.size();++irow)
1298 while(theBuffer[irow].size()&(theBuffer[irow].size()-1))
1299 theBuffer[irow].push_back(theBuffer[irow].back());
1300 double* data=
new double[theBuffer.size()*theBuffer[0].size()];
1301 double* abscoeff=
new double[theBuffer.size()*theBuffer[0].size()];
1302 size_t* p=
new size_t[theBuffer.size()*theBuffer[0].size()];
1303 for(
int irow=0;irow<theBuffer.size();++irow){
1304 for(
int icol=0;icol<theBuffer[0].size();++icol){
1305 int index=irow*theBuffer[0].size()+icol;
1306 assert(index<theBuffer.size()*theBuffer[0].size());
1307 data[index]=theBuffer[irow][icol];
1310 int nsize=theBuffer.size()*theBuffer[0].size();
1312 gsl_wavelet_workspace *work;
1314 w=gsl_wavelet_alloc(filter::Filter::getWaveletType(wavelet_type),family);
1315 work=gsl_wavelet_workspace_alloc(nsize);
1316 gsl_wavelet2d_nstransform_forward (w, data, theBuffer.size(), theBuffer[0].size(),theBuffer[0].size(), work);
1317 for(
int irow=0;irow<theBuffer.size();++irow){
1318 for(
int icol=0;icol<theBuffer[0].size();++icol){
1319 int index=irow*theBuffer[0].size()+icol;
1320 abscoeff[index]=fabs(data[index]);
1323 int nc=(100-cut)/100.0*nsize;
1324 gsl_sort_index(p,abscoeff,1,nsize);
1325 for(
int i=0;(i+nc)<nsize;i++)
1327 gsl_wavelet2d_nstransform_inverse (w, data, theBuffer.size(), theBuffer[0].size(),theBuffer[0].size(), work);
1328 for(
int irow=0;irow<theBuffer.size();++irow){
1329 for(
int icol=0;icol<theBuffer[irow].size();++icol){
1330 int index=irow*theBuffer[irow].size()+icol;
1331 theBuffer[irow][icol]=data[index];
1333 progress=(1.0+irow);
1334 progress/=theBuffer.nRows();
1335 pfnProgress(progress,pszMessage,pProgressArg);
1337 theBuffer.erase(theBuffer.begin()+nRow,theBuffer.end());
1338 for(
int irow=0;irow<theBuffer.size();++irow)
1339 theBuffer[irow].erase(theBuffer[irow].begin()+nCol,theBuffer[irow].end());
1343 gsl_wavelet_free (w);
1344 gsl_wavelet_workspace_free (work);