20 #ifndef _STATFACTORY_H_
21 #define _STATFACTORY_H_
32 #include <gsl/gsl_fit.h>
33 #include <gsl/gsl_errno.h>
34 #include <gsl/gsl_spline.h>
35 #include <gsl/gsl_rng.h>
36 #include <gsl/gsl_randist.h>
37 #include <gsl/gsl_cdf.h>
38 #include <gsl/gsl_statistics.h>
46 enum INTERPOLATION_TYPE {undefined=0,linear=1,polynomial=2,cspline=3,cspline_periodic=4,akima=5,akima_periodic=6};
48 enum DISTRIBUTION_TYPE {uniform=1,gaussian=2};
52 INTERPOLATION_TYPE getInterpolationType(
const std::string interpolationType){
53 std::map<std::string, INTERPOLATION_TYPE> m_interpMap;
55 return m_interpMap[interpolationType];
57 DISTRIBUTION_TYPE getDistributionType(
const std::string distributionType){
58 std::map<std::string, DISTRIBUTION_TYPE> m_distMap;
60 return m_distMap[distributionType];
62 static void allocAcc(gsl_interp_accel *&acc){
63 acc = gsl_interp_accel_alloc ();
66 static void getSpline(
const std::string type,
int size, gsl_spline *& spline){
67 std::map<std::string, INTERPOLATION_TYPE> m_interpMap;
69 switch(m_interpMap[type]){
71 spline=gsl_spline_alloc(gsl_interp_polynomial,size);
74 spline=gsl_spline_alloc(gsl_interp_cspline,size);
76 case(cspline_periodic):
77 spline=gsl_spline_alloc(gsl_interp_cspline_periodic,size);
80 spline=gsl_spline_alloc(gsl_interp_akima,size);
83 spline=gsl_spline_alloc(gsl_interp_akima_periodic,size);
87 spline=gsl_spline_alloc(gsl_interp_linear,size);
92 static int initSpline(gsl_spline *spline,
const double *x,
const double *y,
int size){
93 return gsl_spline_init (spline, x, y, size);
95 static double evalSpline(gsl_spline *spline,
double x, gsl_interp_accel *acc){
96 return gsl_spline_eval (spline, x, acc);
99 static gsl_rng* getRandomGenerator(
unsigned long int theSeed){
100 const gsl_rng_type * T;
104 r = gsl_rng_alloc (gsl_rng_mt19937);
105 gsl_rng_set(r,theSeed);
108 void getNodataValues(std::vector<double>& nodatav)
const{nodatav=m_noDataValues;};
109 bool isNoData(
double value)
const{
110 if(m_noDataValues.empty())
113 return find(m_noDataValues.begin(),m_noDataValues.end(),value)!=m_noDataValues.end();
115 unsigned int pushNodDataValue(
double noDataValue){
116 if(find(m_noDataValues.begin(),m_noDataValues.end(),noDataValue)==m_noDataValues.end())
117 m_noDataValues.push_back(noDataValue);
118 return m_noDataValues.size();
120 unsigned int setNoDataValues(std::vector<double> vnodata){
121 m_noDataValues=vnodata;
122 return m_noDataValues.size();
124 double getRandomValue(
const gsl_rng* r,
const std::string type,
double a=0,
double b=1)
const{
125 std::map<std::string, DISTRIBUTION_TYPE> m_distMap;
128 switch(m_distMap[type]){
130 randValue = a+gsl_rng_uniform(r);
134 randValue = a+gsl_ran_gaussian(r, b);
141 template<
class T> T mymin(
const typename std::vector<T>& v)
const;
142 template<
class T> T mymax(
const typename std::vector<T>& v)
const;
143 template<
class T> T mymin(
const typename std::vector<T>& v, T minConstraint)
const;
144 template<
class T> T mymax(
const typename std::vector<T>& v, T maxConstraint)
const;
146 template<
class T>
typename std::vector<T>::const_iterator mymin(
const typename std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const;
147 template<
class T>
typename std::vector<T>::iterator mymin(
const typename std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end)
const;
148 template<
class T>
typename std::vector<T>::const_iterator mymin(
const typename std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, T minConstraint)
const;
149 template<
class T>
typename std::vector<T>::iterator mymin(
const typename std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end, T minConstraint)
const;
150 template<
class T>
typename std::vector<T>::const_iterator mymax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const;
151 template<
class T>
typename std::vector<T>::iterator mymax(
const std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end)
const;
152 template<
class T>
typename std::vector<T>::const_iterator mymax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, T maxConstraint)
const;
153 template<
class T>
typename std::vector<T>::iterator mymax(
const std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end, T maxConstraint)
const;
154 template<
class T>
typename std::vector<T>::const_iterator absmin(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const;
155 template<
class T>
typename std::vector<T>::const_iterator absmax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const;
157 template<
class T>
void minmax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, T& theMin, T& theMax)
const;
158 template<
class T> T sum(
const std::vector<T>& v)
const;
159 template<
class T>
double mean(
const std::vector<T>& v)
const;
160 template<
class T>
void eraseNoData(std::vector<T>& v)
const;
161 template<
class T> T median(
const std::vector<T>& v)
const;
162 template<
class T>
double var(
const std::vector<T>& v)
const;
163 template<
class T>
double moment(
const std::vector<T>& v,
int n)
const;
164 template<
class T>
double cmoment(
const std::vector<T>& v,
int n)
const;
165 template<
class T>
double skewness(
const std::vector<T>& v)
const;
166 template<
class T>
double kurtosis(
const std::vector<T>& v)
const;
167 template<
class T>
void meanVar(
const std::vector<T>& v,
double& m1,
double& v1)
const;
168 template<
class T1,
class T2>
void scale2byte(
const std::vector<T1>& input, std::vector<T2>& output,
unsigned char lbound=0,
unsigned char ubound=255)
const;
169 template<
class T>
void distribution(
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, std::vector<double>& output,
int nbin, T &minimum, T &maximum,
double sigma=0,
const std::string &filename=
"")
const;
170 template<
class T>
void distribution(
const std::vector<T>& input, std::vector<double>& output,
int nbin,
double sigma=0,
const std::string &filename=
"")
const{distribution(input,input.begin(),input.end(),output,nbin,0,0,sigma,filename);};
171 template<
class T>
void distribution2d(
const std::vector<T>& inputX,
const std::vector<T>& inputY, std::vector< std::vector<double> >& output,
int nbin, T& minX, T& maxX, T& minY, T& maxY,
double sigma=0,
const std::string& filename=
"")
const;
172 template<
class T>
void cumulative (
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, std::vector<int>& output,
int nbin, T &minimum, T &maximum)
const;
173 template<
class T>
void percentiles (
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, std::vector<T>& output,
int nbin, T &minimum, T &maximum,
const std::string &filename=
"")
const;
174 template<
class T> T percentile(
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end,
double percent, T minimum=0, T maximum=0)
const;
175 template<
class T>
void signature(
const std::vector<T>& input,
double& k,
double& alpha,
double& beta,
double e)
const;
176 void signature(
double m1,
double m2,
double& k,
double& alpha,
double& beta,
double e)
const;
177 template<
class T>
void normalize(
const std::vector<T>& input, std::vector<double>& output)
const;
178 template<
class T>
void normalize_pct(std::vector<T>& input)
const;
179 template<
class T>
double rmse(
const std::vector<T>& x,
const std::vector<T>& y)
const;
180 template<
class T>
double correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int delay=0)
const;
182 template<
class T>
double gsl_covariance(
const std::vector<T>& x,
const std::vector<T>& y)
const;
183 template<
class T>
double cross_correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int maxdelay, std::vector<T>& z)
const;
184 template<
class T>
double linear_regression(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const;
185 template<
class T>
double linear_regression_err(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const;
186 template<
class T>
void interpolateNoData(
const std::vector<double>& wavelengthIn,
const std::vector<T>& input,
const std::string& type, std::vector<T>& output,
bool verbose=
false)
const;
187 template<
class T>
void interpolateUp(
const std::vector<double>& wavelengthIn,
const std::vector<T>& input,
const std::vector<double>& wavelengthOut,
const std::string& type, std::vector<T>& output,
bool verbose=
false)
const;
188 template<
class T>
void interpolateUp(
const std::vector<double>& wavelengthIn,
const std::vector< std::vector<T> >& input,
const std::vector<double>& wavelengthOut,
const std::string& type, std::vector< std::vector<T> >& output,
bool verbose=
false)
const;
191 template<
class T>
void interpolateUp(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const;
192 template<
class T>
void nearUp(
const std::vector<T>& input, std::vector<T>& output)
const;
193 template<
class T>
void interpolateUp(
double* input,
int dim, std::vector<T>& output,
int nbin);
194 template<
class T>
void interpolateDown(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const;
195 template<
class T>
void interpolateDown(
double* input,
int dim, std::vector<T>& output,
int nbin);
198 static void initMap(std::map<std::string, INTERPOLATION_TYPE>& m_interpMap){
200 m_interpMap[
"linear"]=linear;
201 m_interpMap[
"polynomial"]=polynomial;
202 m_interpMap[
"cspline"]=cspline;
203 m_interpMap[
"cspline_periodic"]=cspline_periodic;
204 m_interpMap[
"akima"]=akima;
205 m_interpMap[
"akima_periodic"]=akima_periodic;
207 static void initDist(std::map<std::string, DISTRIBUTION_TYPE>& m_distMap){
209 m_distMap[
"gaussian"]=gaussian;
210 m_distMap[
"uniform"]=uniform;
212 std::vector<double> m_noDataValues;
216 template<
class T>
inline typename std::vector<T>::const_iterator StatFactory::mymin(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const
219 typename std::vector<T>::const_iterator tmpIt=begin;
220 for(
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
229 else if(m_noDataValues.size())
230 return m_noDataValues[0];
232 std::string errorString=
"Error: no valid data found";
237 template<
class T>
inline typename std::vector<T>::iterator StatFactory::mymin(
const std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end)
const
240 typename std::vector<T>::iterator tmpIt=begin;
241 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
250 else if(m_noDataValues.size())
251 return m_noDataValues[0];
253 std::string errorString=
"Error: no valid data found";
258 template<
class T>
inline typename std::vector<T>::const_iterator StatFactory::mymin(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, T minConstraint)
const
261 typename std::vector<T>::const_iterator tmpIt=v.end();
262 T minValue=minConstraint;
263 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
267 if((minConstraint<=*it)&&(*it<=minValue)){
274 else if(m_noDataValues.size())
275 return m_noDataValues[0];
277 std::string errorString=
"Error: no valid data found";
282 template<
class T>
inline typename std::vector<T>::iterator StatFactory::mymin(
const std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end, T minConstraint)
const
285 typename std::vector<T>::iterator tmpIt=v.end();
286 T minValue=minConstraint;
287 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
291 if((minConstraint<=*it)&&(*it<=minValue)){
298 else if(m_noDataValues.size())
299 return m_noDataValues[0];
301 std::string errorString=
"Error: no valid data found";
306 template<
class T>
inline typename std::vector<T>::const_iterator StatFactory::mymax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const
309 typename std::vector<T>::const_iterator tmpIt=begin;
310 for (
typename std::vector<T>::iterator it = begin; it!=end; ++it){
319 else if(m_noDataValues.size())
320 return m_noDataValues[0];
322 std::string errorString=
"Error: no valid data found";
327 template<
class T>
inline typename std::vector<T>::iterator StatFactory::mymax(
const std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end)
const
330 typename std::vector<T>::iterator tmpIt=begin;
331 for (
typename std::vector<T>::iterator it = begin; it!=end; ++it){
344 template<
class T>
inline typename std::vector<T>::const_iterator StatFactory::mymax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, T maxConstraint)
const
347 typename std::vector<T>::const_iterator tmpIt=v.end();
348 T maxValue=maxConstraint;
349 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
353 if((maxValue<=*it)&&(*it<=maxConstraint)){
364 template<
class T>
inline typename std::vector<T>::iterator StatFactory::mymax(
const std::vector<T>& v,
typename std::vector<T>::iterator begin,
typename std::vector<T>::iterator end, T maxConstraint)
const
367 typename std::vector<T>::iterator tmpIt=v.end();
368 T maxValue=maxConstraint;
369 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
373 if((maxValue<=*it)&&(*it<=maxConstraint)){
384 template<
class T>
inline T StatFactory::mymin(
const std::vector<T>& v)
const
388 std::string errorString=
"Error: vector is empty";
391 T minValue=*(v.begin());
392 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
401 else if(m_noDataValues.size())
402 return m_noDataValues[0];
404 std::string errorString=
"Error: no valid data found";
409 template<
class T>
inline T StatFactory::mymin(
const std::vector<T>& v, T minConstraint)
const
412 T minValue=minConstraint;
413 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
417 if((minConstraint<=*it)&&(*it<=minValue))
422 else if(m_noDataValues.size())
423 return m_noDataValues[0];
425 std::string errorString=
"Error: no valid data found";
430 template<
class T>
inline T StatFactory::mymax(
const std::vector<T>& v)
const
434 std::string errorString=
"Error: vector is empty";
437 T maxValue=*(v.begin());
438 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
447 else if(m_noDataValues.size())
448 return m_noDataValues[0];
450 std::string errorString=
"Error: no valid data found";
455 template<
class T>
inline T StatFactory::mymax(
const std::vector<T>& v, T maxConstraint)
const
458 T maxValue=maxConstraint;
459 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
462 if((maxValue<=*it)&&(*it<=maxConstraint))
467 else if(m_noDataValues.size())
468 return m_noDataValues[0];
470 std::string errorString=
"Error: no valid data found";
475 template<
class T>
inline typename std::vector<T>::const_iterator StatFactory::absmax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const
478 typename std::vector<T>::const_iterator tmpIt=begin;
479 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
483 if(abs(*tmpIt)<abs(*it))
492 template<
class T>
inline typename std::vector<T>::const_iterator StatFactory::absmin(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end)
const
495 typename std::vector<T>::const_iterator tmpIt=begin;
496 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
500 if(abs(*tmpIt)>abs(*it))
509 template<
class T>
inline void StatFactory::minmax(
const std::vector<T>& v,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, T& theMin, T& theMax)
const
514 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
524 if(m_noDataValues.size()){
525 theMin=m_noDataValues[0];
526 theMax=m_noDataValues[0];
529 std::string errorString=
"Error: no valid data found";
535 template<
class T>
inline T StatFactory::sum(
const std::vector<T>& v)
const
538 typename std::vector<T>::const_iterator it;
540 for (it = v.begin(); it!= v.end(); ++it){
548 else if(m_noDataValues.size())
549 return m_noDataValues[0];
551 std::string errorString=
"Error: no valid data found";
556 template<
class T>
inline double StatFactory::mean(
const std::vector<T>& v)
const
558 typename std::vector<T>::const_iterator it;
560 unsigned int validSize=0;
561 for (it = v.begin(); it!= v.end(); ++it){
568 return static_cast<double>(tmpSum)/validSize;
569 else if(m_noDataValues.size())
570 return m_noDataValues[0];
572 std::string errorString=
"Error: no valid data found";
577 template<
class T>
inline void StatFactory::eraseNoData(std::vector<T>& v)
const
579 if(m_noDataValues.size()){
580 typename std::vector<T>::iterator it=v.begin();
590 template<
class T> T StatFactory::median(
const std::vector<T>& v)
const
592 std::vector<T> tmpV=v;
595 sort(tmpV.begin(),tmpV.end());
597 return tmpV[tmpV.size()/2];
599 return 0.5*(tmpV[tmpV.size()/2-1]+tmpV[tmpV.size()/2]);
601 else if(m_noDataValues.size())
602 return m_noDataValues[0];
604 std::string errorString=
"Error: no valid data found";
609 template<
class T>
double StatFactory::var(
const std::vector<T>& v)
const
611 typename std::vector<T>::const_iterator it;
612 unsigned int validSize=0;
615 for (it = v.begin(); it!= v.end(); ++it){
627 else if(m_noDataValues.size())
628 return m_noDataValues[0];
630 std::string errorString=
"Error: no valid data found";
635 template<
class T>
double StatFactory::moment(
const std::vector<T>& v,
int n)
const
637 unsigned int validSize=0;
638 typename std::vector<T>::const_iterator it;
641 for(it = v.begin(); it!= v.end(); ++it){
649 else if(m_noDataValues.size())
650 return m_noDataValues[0];
652 std::string errorString=
"Error: no valid data found";
658 template<
class T>
double StatFactory::cmoment(
const std::vector<T>& v,
int n)
const
660 unsigned int validSize=0;
661 typename std::vector<T>::const_iterator it;
664 for(it = v.begin(); it!= v.end(); ++it){
672 else if(m_noDataValues.size())
673 return m_noDataValues[0];
675 std::string errorString=
"Error: no valid data found";
680 template<
class T>
double StatFactory::skewness(
const std::vector<T>& v)
const
683 return cmoment(v,3)/pow(var(v),1.5);
686 template<
class T>
double StatFactory::kurtosis(
const std::vector<T>& v)
const
689 double m2=cmoment(v,2);
690 double m4=cmoment(v,4);
694 template<
class T>
void StatFactory::meanVar(
const std::vector<T>& v,
double& m1,
double& v1)
const
696 typename std::vector<T>::const_iterator it;
697 unsigned int validSize=0;
701 for (it = v.begin(); it!= v.end(); ++it){
713 else if(m_noDataValues.size()){
714 m1=m_noDataValues[0];
715 v1=m_noDataValues[0];
718 std::string errorString=
"Error: no valid data found";
723 template<
class T1,
class T2>
void StatFactory::scale2byte(
const std::vector<T1>& input, std::vector<T2>& output,
unsigned char lbound,
unsigned char ubound)
const
725 output.resize(input.size());
726 T1 minimum=mymin(input);
727 T1 maximum=mymax(input);
728 if(minimum>=maximum){
729 std::string errorString=
"Error: no valid data found";
732 double scale=(ubound-lbound)/(maximum-minimum);
734 for (
int i=0;i<input.size();++i){
735 output[i]=scale*(input[i]-(minimum))+lbound;
739 template<
class T>
void StatFactory::distribution(
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, std::vector<double>& output,
int nbin, T &minimum, T &maximum,
double sigma,
const std::string &filename)
const
743 minmax(input,begin,end,minValue,maxValue);
744 if(minimum<maximum&&minimum>minValue)
746 if(minimum<maximum&&maximum<maxValue)
753 if(maximum<=minimum){
754 std::ostringstream s;
755 s<<
"Error: could not calculate distribution (min>=max)";
759 std::string errorString=
"Error: nbin not defined";
763 std::string errorString=
"Error: no valid data found";
766 if(output.size()!=nbin){
768 for(
int i=0;i<nbin;output[i++]=0);
771 typename std::vector<T>::const_iterator it;
772 for(it=begin;it!=end;++it){
786 for(
int ibin=0;ibin<nbin;++ibin){
787 double icenter=minimum+
static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin;
788 double thePdf=gsl_ran_gaussian_pdf(*it-icenter, sigma);
789 output[ibin]+=thePdf;
796 else if(*it>minimum && *it<maximum)
797 theBin=
static_cast<int>(
static_cast<double>((nbin-1)*(*it)-minimum)/(maximum-minimum));
806 std::string errorString=
"Error: no valid data found";
809 else if(!filename.empty()){
810 std::ofstream outputfile;
811 outputfile.open(filename.c_str());
813 std::ostringstream s;
814 s<<
"Error opening distribution file , " << filename;
817 for(
int ibin=0;ibin<nbin;++ibin)
818 outputfile << minimum+static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin <<
" " <<
static_cast<double>(output[ibin])/input.size() << std::endl;
823 template<
class T>
void StatFactory::distribution2d(
const std::vector<T>& inputX,
const std::vector<T>& inputY, std::vector< std::vector<double> >& output,
int nbin, T& minX, T& maxX, T& minY, T& maxY,
double sigma,
const std::string& filename)
const
826 std::ostringstream s;
827 s<<
"Error: inputX is empty";
830 if(inputX.size()!=inputY.size()){
831 std::ostringstream s;
832 s<<
"Error: inputX is empty";
835 unsigned long int npoint=inputX.size();
837 minmax(inputX,inputX.begin(),inputX.end(),minX,maxX);
839 std::ostringstream s;
840 s<<
"Error: could not calculate distribution (minX>=maxX)";
844 minmax(inputY,inputY.begin(),inputY.end(),minY,maxY);
846 std::ostringstream s;
847 s<<
"Error: could not calculate distribution (minY>=maxY)";
851 std::ostringstream s;
852 s<<
"Error: nbin must be larger than 1";
856 for(
int i=0;i<nbin;++i){
857 output[i].resize(nbin);
858 for(
int j=0;j<nbin;++j)
863 for(
unsigned long int ipoint=0;ipoint<npoint;++ipoint){
864 if(inputX[ipoint]==maxX)
867 binX=
static_cast<int>(
static_cast<double>(inputX[ipoint]-minX)/(maxX-minX)*nbin);
868 if(inputY[ipoint]==maxY)
871 binY=
static_cast<int>(
static_cast<double>(inputY[ipoint]-minY)/(maxY-minY)*nbin);
873 std::ostringstream s;
874 s<<
"Error: binX is smaller than 0";
877 if(output.size()<=binX){
878 std::ostringstream s;
879 s<<
"Error: output size must be larger than binX";
883 std::ostringstream s;
884 s<<
"Error: binY is smaller than 0";
887 if(output.size()<=binY){
888 std::ostringstream s;
889 s<<
"Error: output size must be larger than binY";
899 for(
int ibinX=0;ibinX<nbin;++ibinX){
900 double centerX=minX+
static_cast<double>(maxX-minX)*ibinX/nbin;
901 double pdfX=gsl_ran_gaussian_pdf(inputX[ipoint]-centerX, sigma);
902 for(
int ibinY=0;ibinY<nbin;++ibinY){
904 double centerY=minY+
static_cast<double>(maxY-minY)*ibinY/nbin;
905 double pdfY=gsl_ran_gaussian_pdf(inputY[ipoint]-centerY, sigma);
906 output[ibinX][binY]+=pdfX*pdfY;
911 ++output[binX][binY];
913 if(!filename.empty()){
914 std::ofstream outputfile;
915 outputfile.open(filename.c_str());
917 std::ostringstream s;
918 s<<
"Error opening distribution file , " << filename;
921 for(
int binX=0;binX<nbin;++binX){
922 outputfile << std::endl;
923 for(
int binY=0;binY<nbin;++binY){
925 if(nbin==maxX-minX+1)
928 binValueX=minX+
static_cast<double>(maxX-minX)*(binX+0.5)/nbin;
930 if(nbin==maxY-minY+1)
933 binValueY=minY+
static_cast<double>(maxY-minY)*(binY+0.5)/nbin;
935 value=
static_cast<double>(output[binX][binY])/npoint;
936 outputfile << binValueX <<
" " << binValueY <<
" " << value << std::endl;
946 template<
class T>
void StatFactory::percentiles (
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end, std::vector<T>& output,
int nbin, T &minimum, T &maximum,
const std::string &filename)
const
949 minmax(input,begin,end,minimum,maximum);
950 if(maximum<=minimum){
951 std::ostringstream s;
952 s<<
"Error: maximum must be at least minimum";
956 std::ostringstream s;
957 s<<
"Error: nbin must be larger than 1";
961 std::ostringstream s;
962 s<<
"Error: input is empty";
966 std::vector<T> inputSort;
967 inputSort.assign(begin,end);
968 typename std::vector<T>::iterator vit=inputSort.begin();
969 while(vit!=inputSort.end()){
970 if(*vit<minimum||*vit>maximum)
971 inputSort.erase(vit);
975 std::sort(inputSort.begin(),inputSort.end());
976 vit=inputSort.begin();
977 std::vector<T> inputBin;
978 for(
int ibin=0;ibin<nbin;++ibin){
980 while(inputBin.size()<inputSort.size()/nbin&&vit!=inputSort.end()){
981 inputBin.push_back(*vit);
985 output[ibin]=mymax(inputBin);
988 if(!filename.empty()){
989 std::ofstream outputfile;
990 outputfile.open(filename.c_str());
992 std::ostringstream s;
993 s<<
"error opening distribution file , " << filename;
996 for(
int ibin=0;ibin<nbin;++ibin)
997 outputfile << ibin*100.0/nbin <<
" " << static_cast<double>(output[ibin])/input.size() << std::endl;
1003 template<
class T> T StatFactory::percentile(
const std::vector<T>& input,
typename std::vector<T>::const_iterator begin,
typename std::vector<T>::const_iterator end,
double percent, T minimum, T maximum)
const
1006 std::ostringstream s;
1007 s<<
"Error: input is empty";
1010 std::vector<T> inputSort;
1011 inputSort.assign(begin,end);
1012 typename std::vector<T>::iterator vit=inputSort.begin();
1013 while(vit!=inputSort.end()){
1014 if(maximum>minimum){
1015 if(*vit<minimum||*vit>maximum)
1016 inputSort.erase(vit);
1021 std::sort(inputSort.begin(),inputSort.end());
1022 return gsl_stats_quantile_from_sorted_data(&(inputSort[0]),1,inputSort.size(),percent/100.0);
1025 template<
class T>
void StatFactory::signature(
const std::vector<T>& input,
double&k,
double& alpha,
double& beta,
double e)
const
1027 double m1=moment(input,1);
1028 double m2=moment(input,2);
1029 signature(m1,m2,k,alpha,beta,e);
1033 template<
class T>
void StatFactory::normalize(
const std::vector<T>& input, std::vector<double>& output)
const{
1034 double total=sum(input);
1036 output.resize(input.size());
1037 for(
int index=0;index<input.size();++index)
1038 output[index]=input[index]/total;
1045 template<
class T>
void StatFactory::normalize_pct(std::vector<T>& input)
const{
1046 double total=sum(input);
1048 typename std::vector<T>::iterator it;
1049 for(it=input.begin();it!=input.end();++it)
1050 *it=100.0*(*it)/total;
1054 template<
class T>
double StatFactory::rmse(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1055 if(x.size()!=y.size()){
1056 std::ostringstream s;
1057 s<<
"Error: x and y not equal in size";
1061 std::ostringstream s;
1062 s<<
"Error: x is empty";
1066 for(
int isample=0;isample<x.size();++isample){
1067 if(isNoData(x[isample])||isNoData(y[isample]))
1069 double e=x[isample]-y[isample];
1079 template<
class T>
double StatFactory::gsl_covariance(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1080 return(gsl_stats_covariance(&(x[0]),1,&(y[0]),1,x.size()));
1084 template<
class T>
double StatFactory::correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int delay)
const{
1090 meanVar(x,meanX,varX);
1091 meanVar(y,meanY,varY);
1092 double denom = sqrt(varX*varY);
1097 for (
int i=0;i<x.size();++i) {
1099 if (j < 0 || j >= y.size())
1101 else if(isNoData(x[i])||isNoData(y[j]))
1106 std::ostringstream s;
1107 s<<
"Error: i must be positive";
1111 std::ostringstream s;
1112 s<<
"Error: i must be smaller than x.size()";
1116 std::ostringstream s;
1117 s<<
"Error: j must be positive";
1121 std::ostringstream s;
1122 s<<
"Error: j must be smaller than y.size()";
1125 sXY += (x[i] - meanX) * (y[j] - meanY);
1129 double minSize=(x.size()<y.size())?x.size():y.size();
1130 return(sXY / denom / (minSize-1));
1132 else if(m_noDataValues.size())
1133 return m_noDataValues[0];
1135 std::string errorString=
"Error: no valid data found";
1144 template<
class T>
double StatFactory::cross_correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int maxdelay, std::vector<T>& z)
const{
1146 double sumCorrelation=0;
1147 for (
int delay=-maxdelay;delay<maxdelay;delay++) {
1148 z.push_back(correlation(x,y,delay));
1149 sumCorrelation+=z.back();
1151 return sumCorrelation;
1155 template<
class T>
double StatFactory::linear_regression(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const{
1156 if(x.size()!=y.size()){
1157 std::ostringstream s;
1158 s<<
"Error: x and y not equal in size";
1162 std::ostringstream s;
1163 s<<
"Error: x is empty";
1170 gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1171 return (1-sumsq/var(y)/(y.size()-1));
1175 template<
class T>
double StatFactory::linear_regression_err(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const{
1176 if(x.size()!=y.size()){
1177 std::ostringstream s;
1178 s<<
"Error: x and y not equal in size";
1182 std::ostringstream s;
1183 s<<
"Error: x is empty";
1190 gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1191 return sqrt((sumsq)/(y.size()));
1197 template<
class T>
void StatFactory::interpolateNoData(
const std::vector<double>& wavelengthIn,
const std::vector<T>& input,
const std::string& type, std::vector<T>& output,
bool verbose)
const{
1198 if(wavelengthIn.empty()){
1199 std::ostringstream s;
1200 s<<
"Error: wavelengthIn is empty";
1203 std::vector<double> wavelengthOut=wavelengthIn;
1204 std::vector<T> validIn=input;
1205 if(input.size()!=wavelengthIn.size()){
1206 std::ostringstream s;
1207 s<<
"Error: x and y not equal in size";
1210 int nband=wavelengthIn.size();
1213 if(m_noDataValues.size()){
1214 typename std::vector<T>::iterator itValue=validIn.begin();
1215 typename std::vector<T>::iterator itWavelength=wavelengthOut.begin();
1216 while(itValue!=validIn.end()&&itWavelength!=wavelengthOut.end()){
1217 if(isNoData(*itValue)){
1218 validIn.erase(itValue);
1219 wavelengthOut.erase(itWavelength);
1226 if(validIn.size()>1){
1228 interpolateUp(wavelengthOut, validIn, wavelengthIn, type, output, verbose);
1241 template<
class T>
void StatFactory::interpolateUp(
const std::vector<double>& wavelengthIn,
const std::vector<T>& input,
const std::vector<double>& wavelengthOut,
const std::string& type, std::vector<T>& output,
bool verbose)
const{
1242 if(wavelengthIn.empty()){
1243 std::ostringstream s;
1244 s<<
"Error: wavelengthIn is empty";
1247 if(input.size()!=wavelengthIn.size()){
1248 std::ostringstream s;
1249 s<<
"Error: input and wavelengthIn not equal in size";
1252 if(wavelengthOut.empty()){
1253 std::ostringstream s;
1254 s<<
"Error: wavelengthOut is empty";
1257 int nband=wavelengthIn.size();
1259 gsl_interp_accel *acc;
1262 getSpline(type,nband,spline);
1264 assert(&(wavelengthIn[0]));
1265 assert(&(input[0]));
1266 int status=initSpline(spline,&(wavelengthIn[0]),&(input[0]),nband);
1268 std::string errorString=
"Could not initialize spline";
1271 for(
int index=0;index<wavelengthOut.size();++index){
1272 if(wavelengthOut[index]<*wavelengthIn.begin()){
1273 output.push_back(*(input.begin()));
1276 else if(wavelengthOut[index]>wavelengthIn.back()){
1277 output.push_back(input.back());
1280 double dout=evalSpline(spline,wavelengthOut[index],acc);
1281 output.push_back(dout);
1283 gsl_spline_free(spline);
1284 gsl_interp_accel_free(acc);
1319 template<
class T>
void StatFactory::interpolateUp(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const
1322 std::ostringstream s;
1323 s<<
"Error: input is empty";
1327 std::ostringstream s;
1328 s<<
"Error: nbin must be larger than 0";
1332 int dim=input.size();
1333 for(
int i=0;i<dim;++i){
1335 double left=input[i];
1337 double right=(i<dim-1)? input[i+1]:input[i];
1338 deltaX=(right-left)/static_cast<double>(nbin);
1339 for(
int x=0;x<nbin;++x){
1340 output.push_back(left+x*deltaX);
1344 output.push_back(input.back());
1349 template<
class T>
void StatFactory::nearUp(
const std::vector<T>& input, std::vector<T>& output)
const
1352 std::ostringstream s;
1353 s<<
"Error: input is empty";
1356 if(output.size()<input.size()){
1357 std::ostringstream s;
1358 s<<
"Error: output size is smaller than input size";
1361 int dimInput=input.size();
1362 int dimOutput=output.size();
1364 for(
int iin=0;iin<dimInput;++iin){
1365 for(
int iout=0;iout<dimOutput/dimInput;++iout){
1366 int indexOutput=iin*dimOutput/dimInput+iout;
1367 if(indexOutput>=output.size()){
1368 std::ostringstream s;
1369 s<<
"Error: indexOutput must be smaller than output.size()";
1372 output[indexOutput]=input[iin];
1378 template<
class T>
void StatFactory::interpolateUp(
double* input,
int dim, std::vector<T>& output,
int nbin)
1381 std::ostringstream s;
1382 s<<
"Error: nbin must be larger than 0";
1386 for(
int i=0;i<dim;++i){
1388 double left=input[i];
1390 double right=(i<dim-1)? input[i+1]:input[i];
1391 deltaX=(right-left)/static_cast<double>(nbin);
1392 for(
int x=0;x<nbin;++x){
1393 output.push_back(left+x*deltaX);
1397 output.push_back(input[dim-1]);
1402 template<
class T>
void StatFactory::interpolateDown(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const
1405 std::ostringstream s;
1406 s<<
"Error: input is empty";
1410 std::ostringstream s;
1411 s<<
"Error: nbin must be larger than 0";
1415 int dim=input.size();
1417 output.push_back(input[0]);
1418 for(
int i=1;i<dim;++i){
1423 output.push_back(input[i]);
1429 template<
class T>
void StatFactory::interpolateDown(
double* input,
int dim, std::vector<T>& output,
int nbin)
1432 std::ostringstream s;
1433 s<<
"Error: nbin must be larger than 0";
1438 output.push_back(input[0]);
1439 for(
int i=1;i<dim;++i){
1444 output.push_back(input[i]);