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>
unsigned int nvalid(
const std::vector<T>& v)
const;
162 template<
class T> T median(
const std::vector<T>& v)
const;
163 template<
class T>
double var(
const std::vector<T>& v)
const;
164 template<
class T>
double moment(
const std::vector<T>& v,
int n)
const;
165 template<
class T>
double cmoment(
const std::vector<T>& v,
int n)
const;
166 template<
class T>
double skewness(
const std::vector<T>& v)
const;
167 template<
class T>
double kurtosis(
const std::vector<T>& v)
const;
168 template<
class T>
void meanVar(
const std::vector<T>& v,
double& m1,
double& v1)
const;
169 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;
170 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;
171 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);};
172 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;
173 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;
174 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;
175 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;
176 template<
class T>
void signature(
const std::vector<T>& input,
double& k,
double& alpha,
double& beta,
double e)
const;
177 void signature(
double m1,
double m2,
double& k,
double& alpha,
double& beta,
double e)
const;
178 template<
class T>
void normalize(
const std::vector<T>& input, std::vector<double>& output)
const;
179 template<
class T>
void normalize_pct(std::vector<T>& input)
const;
180 template<
class T>
double rmse(
const std::vector<T>& x,
const std::vector<T>& y)
const;
181 template<
class T>
double nrmse(
const std::vector<T>& x,
const std::vector<T>& y)
const;
182 template<
class T>
double cvrmse(
const std::vector<T>& x,
const std::vector<T>& y)
const;
183 template<
class T>
double correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int delay=0)
const;
185 template<
class T>
double gsl_covariance(
const std::vector<T>& x,
const std::vector<T>& y)
const;
186 template<
class T>
double cross_correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int maxdelay, std::vector<T>& z)
const;
187 template<
class T>
double linear_regression(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const;
188 template<
class T>
double linear_regression_err(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const;
189 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;
190 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;
191 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;
194 template<
class T>
void interpolateUp(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const;
195 template<
class T>
void nearUp(
const std::vector<T>& input, std::vector<T>& output)
const;
196 template<
class T>
void interpolateUp(
double* input,
int dim, std::vector<T>& output,
int nbin);
197 template<
class T>
void interpolateDown(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const;
198 template<
class T>
void interpolateDown(
double* input,
int dim, std::vector<T>& output,
int nbin);
201 static void initMap(std::map<std::string, INTERPOLATION_TYPE>& m_interpMap){
203 m_interpMap[
"linear"]=linear;
204 m_interpMap[
"polynomial"]=polynomial;
205 m_interpMap[
"cspline"]=cspline;
206 m_interpMap[
"cspline_periodic"]=cspline_periodic;
207 m_interpMap[
"akima"]=akima;
208 m_interpMap[
"akima_periodic"]=akima_periodic;
210 static void initDist(std::map<std::string, DISTRIBUTION_TYPE>& m_distMap){
212 m_distMap[
"gaussian"]=gaussian;
213 m_distMap[
"uniform"]=uniform;
215 std::vector<double> m_noDataValues;
219 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
222 typename std::vector<T>::const_iterator tmpIt;
223 for(
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
237 else if(m_noDataValues.size())
238 return m_noDataValues[0];
240 std::string errorString=
"Error: no valid data found";
245 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
248 typename std::vector<T>::iterator tmpIt;
249 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
263 else if(m_noDataValues.size())
264 return m_noDataValues[0];
266 std::string errorString=
"Error: no valid data found";
271 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
274 typename std::vector<T>::const_iterator tmpIt;
275 T minValue=minConstraint;
276 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
280 if((minConstraint<=*it)&&(*it<minValue)){
295 else if(m_noDataValues.size())
296 return m_noDataValues[0];
298 std::string errorString=
"Error: no valid data found";
303 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
306 typename std::vector<T>::iterator tmpIt;
307 T minValue=minConstraint;
308 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
312 if((minConstraint<=*it)&&(*it<minValue)){
318 if(*it<minConstraint)
327 else if(m_noDataValues.size())
328 return m_noDataValues[0];
330 std::string errorString=
"Error: no valid data found";
335 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
338 typename std::vector<T>::const_iterator tmpIt;
339 for (
typename std::vector<T>::iterator it = begin; it!=end; ++it){
353 else if(m_noDataValues.size())
354 return m_noDataValues[0];
356 std::string errorString=
"Error: no valid data found";
361 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
364 typename std::vector<T>::iterator tmpIt;
365 for (
typename std::vector<T>::iterator it = begin; it!=end; ++it){
383 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
386 typename std::vector<T>::const_iterator tmpIt;
387 T maxValue=maxConstraint;
388 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
392 if((maxConstraint>=*it)&&(*it>maxValue)){
398 if(*it>maxConstraint)
411 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
414 typename std::vector<T>::iterator tmpIt=v.end();
415 T maxValue=maxConstraint;
416 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
420 if((maxConstraint>=*it)&&(*it>maxValue)){
439 template<
class T>
inline T StatFactory::mymin(
const std::vector<T>& v)
const
443 std::string errorString=
"Error: vector is empty";
447 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
461 else if(m_noDataValues.size())
462 return m_noDataValues[0];
464 std::string errorString=
"Error: no valid data found";
469 template<
class T>
inline T StatFactory::mymin(
const std::vector<T>& v, T minConstraint)
const
472 T minValue=minConstraint;
473 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
477 if((minConstraint<=*it)&&(*it<minValue))
489 else if(m_noDataValues.size())
490 return m_noDataValues[0];
492 std::string errorString=
"Error: no valid data found";
497 template<
class T>
inline T StatFactory::mymax(
const std::vector<T>& v)
const
501 std::string errorString=
"Error: vector is empty";
505 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
519 else if(m_noDataValues.size())
520 return m_noDataValues[0];
522 std::string errorString=
"Error: no valid data found";
527 template<
class T>
inline T StatFactory::mymax(
const std::vector<T>& v, T maxConstraint)
const
530 T maxValue=maxConstraint;
531 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
535 if((*it<=maxConstraint)&&(*it>maxValue))
547 else if(m_noDataValues.size())
548 return m_noDataValues[0];
550 std::string errorString=
"Error: no valid data found";
555 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
558 typename std::vector<T>::const_iterator tmpIt;
559 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
563 if(abs(*tmpIt)<abs(*it))
577 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
580 typename std::vector<T>::const_iterator tmpIt;
581 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
585 if(abs(*tmpIt)>abs(*it))
599 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
601 bool isConstraint=(theMax>theMin);
602 double minConstraint=theMin;
603 double maxConstraint=theMax;
605 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
610 if(*it<minConstraint)
612 if(*it>maxConstraint)
622 if(*it<minConstraint)
624 if(*it>maxConstraint)
633 if(m_noDataValues.size()){
634 theMin=m_noDataValues[0];
635 theMax=m_noDataValues[0];
638 std::string errorString=
"Error: no valid data found";
644 template<
class T>
inline T StatFactory::sum(
const std::vector<T>& v)
const
647 typename std::vector<T>::const_iterator it;
649 for (it = v.begin(); it!= v.end(); ++it){
657 else if(m_noDataValues.size())
658 return m_noDataValues[0];
660 std::string errorString=
"Error: no valid data found";
665 template<
class T>
inline double StatFactory::mean(
const std::vector<T>& v)
const
667 typename std::vector<T>::const_iterator it;
669 unsigned int validSize=0;
670 for (it = v.begin(); it!= v.end(); ++it){
677 return static_cast<double>(tmpSum)/validSize;
678 else if(m_noDataValues.size())
679 return m_noDataValues[0];
681 std::string errorString=
"Error: no valid data found";
686 template<
class T>
inline void StatFactory::eraseNoData(std::vector<T>& v)
const
688 if(m_noDataValues.size()){
689 typename std::vector<T>::iterator it=v.begin();
699 template<
class T>
unsigned int StatFactory::nvalid(
const std::vector<T>& v)
const{
700 std::vector<T> tmpV=v;
705 template<
class T> T StatFactory::median(
const std::vector<T>& v)
const
707 std::vector<T> tmpV=v;
710 sort(tmpV.begin(),tmpV.end());
712 return tmpV[tmpV.size()/2];
714 return 0.5*(tmpV[tmpV.size()/2-1]+tmpV[tmpV.size()/2]);
716 else if(m_noDataValues.size())
717 return m_noDataValues[0];
719 std::string errorString=
"Error: no valid data found";
724 template<
class T>
double StatFactory::var(
const std::vector<T>& v)
const
726 typename std::vector<T>::const_iterator it;
727 unsigned int validSize=0;
730 for (it = v.begin(); it!= v.end(); ++it){
742 else if(m_noDataValues.size())
743 return m_noDataValues[0];
745 std::string errorString=
"Error: no valid data found";
750 template<
class T>
double StatFactory::moment(
const std::vector<T>& v,
int n)
const
752 unsigned int validSize=0;
753 typename std::vector<T>::const_iterator it;
756 for(it = v.begin(); it!= v.end(); ++it){
764 else if(m_noDataValues.size())
765 return m_noDataValues[0];
767 std::string errorString=
"Error: no valid data found";
773 template<
class T>
double StatFactory::cmoment(
const std::vector<T>& v,
int n)
const
775 unsigned int validSize=0;
776 typename std::vector<T>::const_iterator it;
779 for(it = v.begin(); it!= v.end(); ++it){
787 else if(m_noDataValues.size())
788 return m_noDataValues[0];
790 std::string errorString=
"Error: no valid data found";
795 template<
class T>
double StatFactory::skewness(
const std::vector<T>& v)
const
798 return cmoment(v,3)/pow(var(v),1.5);
801 template<
class T>
double StatFactory::kurtosis(
const std::vector<T>& v)
const
804 double m2=cmoment(v,2);
805 double m4=cmoment(v,4);
809 template<
class T>
void StatFactory::meanVar(
const std::vector<T>& v,
double& m1,
double& v1)
const
811 typename std::vector<T>::const_iterator it;
812 unsigned int validSize=0;
816 for (it = v.begin(); it!= v.end(); ++it){
828 else if(m_noDataValues.size()){
829 m1=m_noDataValues[0];
830 v1=m_noDataValues[0];
833 std::string errorString=
"Error: no valid data found";
838 template<
class T1,
class T2>
void StatFactory::scale2byte(
const std::vector<T1>& input, std::vector<T2>& output,
unsigned char lbound,
unsigned char ubound)
const
840 output.resize(input.size());
841 T1 minimum=mymin(input);
842 T1 maximum=mymax(input);
843 if(minimum>=maximum){
844 std::string errorString=
"Error: no valid data found";
847 double scale=(ubound-lbound)/(maximum-minimum);
849 for (
int i=0;i<input.size();++i){
850 output[i]=scale*(input[i]-(minimum))+lbound;
854 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
858 minmax(input,begin,end,minimum,maximum);
859 if(minimum<maximum&&minimum>minValue)
861 if(minimum<maximum&&maximum<maxValue)
868 if(maximum<=minimum){
869 std::ostringstream s;
870 s<<
"Error: could not calculate distribution (min>=max)";
874 std::string errorString=
"Error: nbin not defined";
878 std::string errorString=
"Error: no valid data found";
881 if(output.size()!=nbin){
883 for(
int i=0;i<nbin;output[i++]=0);
886 typename std::vector<T>::const_iterator it;
887 for(it=begin;it!=end;++it){
901 for(
int ibin=0;ibin<nbin;++ibin){
902 double icenter=minimum+
static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin;
903 double thePdf=gsl_ran_gaussian_pdf(*it-icenter, sigma);
904 output[ibin]+=thePdf;
911 else if(*it>minimum && *it<maximum)
912 theBin=
static_cast<int>(
static_cast<double>((nbin-1)*(*it)-minimum)/(maximum-minimum));
921 std::string errorString=
"Error: no valid data found";
924 else if(!filename.empty()){
925 std::ofstream outputfile;
926 outputfile.open(filename.c_str());
928 std::ostringstream s;
929 s<<
"Error opening distribution file , " << filename;
932 for(
int ibin=0;ibin<nbin;++ibin)
933 outputfile << minimum+static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin <<
" " <<
static_cast<double>(output[ibin])/input.size() << std::endl;
938 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
941 std::ostringstream s;
942 s<<
"Error: inputX is empty";
945 if(inputX.size()!=inputY.size()){
946 std::ostringstream s;
947 s<<
"Error: inputX is empty";
950 unsigned long int npoint=inputX.size();
952 minmax(inputX,inputX.begin(),inputX.end(),minX,maxX);
954 std::ostringstream s;
955 s<<
"Error: could not calculate distribution (minX>=maxX)";
959 minmax(inputY,inputY.begin(),inputY.end(),minY,maxY);
961 std::ostringstream s;
962 s<<
"Error: could not calculate distribution (minY>=maxY)";
966 std::ostringstream s;
967 s<<
"Error: nbin must be larger than 1";
971 for(
int i=0;i<nbin;++i){
972 output[i].resize(nbin);
973 for(
int j=0;j<nbin;++j)
978 for(
unsigned long int ipoint=0;ipoint<npoint;++ipoint){
979 if(inputX[ipoint]==maxX)
982 binX=
static_cast<int>(
static_cast<double>(inputX[ipoint]-minX)/(maxX-minX)*nbin);
983 if(inputY[ipoint]==maxY)
986 binY=
static_cast<int>(
static_cast<double>(inputY[ipoint]-minY)/(maxY-minY)*nbin);
988 std::ostringstream s;
989 s<<
"Error: binX is smaller than 0";
992 if(output.size()<=binX){
993 std::ostringstream s;
994 s<<
"Error: output size must be larger than binX";
998 std::ostringstream s;
999 s<<
"Error: binY is smaller than 0";
1002 if(output.size()<=binY){
1003 std::ostringstream s;
1004 s<<
"Error: output size must be larger than binY";
1014 for(
int ibinX=0;ibinX<nbin;++ibinX){
1015 double centerX=minX+
static_cast<double>(maxX-minX)*ibinX/nbin;
1016 double pdfX=gsl_ran_gaussian_pdf(inputX[ipoint]-centerX, sigma);
1017 for(
int ibinY=0;ibinY<nbin;++ibinY){
1019 double centerY=minY+
static_cast<double>(maxY-minY)*ibinY/nbin;
1020 double pdfY=gsl_ran_gaussian_pdf(inputY[ipoint]-centerY, sigma);
1021 output[ibinX][binY]+=pdfX*pdfY;
1026 ++output[binX][binY];
1028 if(!filename.empty()){
1029 std::ofstream outputfile;
1030 outputfile.open(filename.c_str());
1032 std::ostringstream s;
1033 s<<
"Error opening distribution file , " << filename;
1036 for(
int binX=0;binX<nbin;++binX){
1037 outputfile << std::endl;
1038 for(
int binY=0;binY<nbin;++binY){
1040 if(nbin==maxX-minX+1)
1041 binValueX=minX+binX;
1043 binValueX=minX+
static_cast<double>(maxX-minX)*(binX+0.5)/nbin;
1045 if(nbin==maxY-minY+1)
1046 binValueY=minY+binY;
1048 binValueY=minY+
static_cast<double>(maxY-minY)*(binY+0.5)/nbin;
1050 value=
static_cast<double>(output[binX][binY])/npoint;
1051 outputfile << binValueX <<
" " << binValueY <<
" " << value << std::endl;
1061 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
1063 if(maximum<=minimum)
1064 minmax(input,begin,end,minimum,maximum);
1065 if(maximum<=minimum){
1066 std::ostringstream s;
1067 s<<
"Error: maximum must be at least minimum";
1071 std::ostringstream s;
1072 s<<
"Error: nbin must be larger than 1";
1076 std::ostringstream s;
1077 s<<
"Error: input is empty";
1080 output.resize(nbin);
1081 std::vector<T> inputSort;
1082 inputSort.assign(begin,end);
1083 typename std::vector<T>::iterator vit=inputSort.begin();
1084 while(vit!=inputSort.end()){
1085 if(*vit<minimum||*vit>maximum)
1086 inputSort.erase(vit);
1090 eraseNoData(inputSort);
1091 std::sort(inputSort.begin(),inputSort.end());
1092 vit=inputSort.begin();
1093 std::vector<T> inputBin;
1094 for(
int ibin=0;ibin<nbin;++ibin){
1096 while(inputBin.size()<inputSort.size()/nbin&&vit!=inputSort.end()){
1097 inputBin.push_back(*vit);
1100 if(inputBin.size()){
1101 output[ibin]=mymax(inputBin);
1104 if(!filename.empty()){
1105 std::ofstream outputfile;
1106 outputfile.open(filename.c_str());
1108 std::ostringstream s;
1109 s<<
"error opening distribution file , " << filename;
1112 for(
int ibin=0;ibin<nbin;++ibin)
1113 outputfile << ibin*100.0/nbin <<
" " << static_cast<double>(output[ibin])/input.size() << std::endl;
1119 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
1122 std::ostringstream s;
1123 s<<
"Error: input is empty";
1126 std::vector<T> inputSort;
1127 inputSort.assign(begin,end);
1128 typename std::vector<T>::iterator vit=inputSort.begin();
1129 while(vit!=inputSort.end()){
1130 if(maximum>minimum){
1131 if(*vit<minimum||*vit>maximum)
1132 inputSort.erase(vit);
1137 eraseNoData(inputSort);
1138 std::sort(inputSort.begin(),inputSort.end());
1139 return gsl_stats_quantile_from_sorted_data(&(inputSort[0]),1,inputSort.size(),percent/100.0);
1142 template<
class T>
void StatFactory::signature(
const std::vector<T>& input,
double&k,
double& alpha,
double& beta,
double e)
const
1144 double m1=moment(input,1);
1145 double m2=moment(input,2);
1146 signature(m1,m2,k,alpha,beta,e);
1150 template<
class T>
void StatFactory::normalize(
const std::vector<T>& input, std::vector<double>& output)
const{
1151 double total=sum(input);
1153 output.resize(input.size());
1154 for(
int index=0;index<input.size();++index)
1155 output[index]=input[index]/total;
1162 template<
class T>
void StatFactory::normalize_pct(std::vector<T>& input)
const{
1163 double total=sum(input);
1165 typename std::vector<T>::iterator it;
1166 for(it=input.begin();it!=input.end();++it)
1167 *it=100.0*(*it)/total;
1171 template<
class T>
double StatFactory::rmse(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1172 if(x.size()!=y.size()){
1173 std::ostringstream s;
1174 s<<
"Error: x and y not equal in size";
1178 std::ostringstream s;
1179 s<<
"Error: x is empty";
1183 for(
int isample=0;isample<x.size();++isample){
1184 if(isNoData(x[isample])||isNoData(y[isample]))
1186 double e=x[isample]-y[isample];
1193 template<
class T>
double StatFactory::nrmse(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1194 if(x.size()!=y.size()){
1195 std::ostringstream s;
1196 s<<
"Error: x and y not equal in size";
1200 std::ostringstream s;
1201 s<<
"Error: x is empty";
1204 std::vector<T> tmpX=x;
1206 std::vector<T> tmpY=y;
1208 double maxY=mymax(y);
1209 double minY=mymin(y);
1210 double rangeY=maxY-minY;
1212 for(
int isample=0;isample<x.size();++isample){
1213 double e=x[isample]-y[isample];
1216 return sqrt(mse)/rangeY;
1220 template<
class T>
double StatFactory::cvrmse(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1221 if(x.size()!=y.size()){
1222 std::ostringstream s;
1223 s<<
"Error: x and y not equal in size";
1227 std::ostringstream s;
1228 s<<
"Error: x is empty";
1231 std::vector<T> tmpX=x;
1233 std::vector<T> tmpY=y;
1235 double maxY=mymax(tmpY);
1236 double minY=mymin(tmpY);
1237 double rangeY=maxY-minY;
1239 for(
int isample=0;isample<x.size();++isample){
1240 double e=x[isample]-y[isample];
1243 return sqrt(mse)/mean(tmpY);
1250 template<
class T>
double StatFactory::gsl_covariance(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1251 return(gsl_stats_covariance(&(x[0]),1,&(y[0]),1,x.size()));
1255 template<
class T>
double StatFactory::correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int delay)
const{
1261 meanVar(x,meanX,varX);
1262 meanVar(y,meanY,varY);
1263 double denom = sqrt(varX*varY);
1268 for (
int i=0;i<x.size();++i) {
1270 if (j < 0 || j >= y.size())
1272 else if(isNoData(x[i])||isNoData(y[j]))
1277 std::ostringstream s;
1278 s<<
"Error: i must be positive";
1282 std::ostringstream s;
1283 s<<
"Error: i must be smaller than x.size()";
1287 std::ostringstream s;
1288 s<<
"Error: j must be positive";
1292 std::ostringstream s;
1293 s<<
"Error: j must be smaller than y.size()";
1296 sXY += (x[i] - meanX) * (y[j] - meanY);
1300 double minSize=(x.size()<y.size())?x.size():y.size();
1301 return(sXY / denom / (minSize-1));
1303 else if(m_noDataValues.size())
1304 return m_noDataValues[0];
1306 std::string errorString=
"Error: no valid data found";
1315 template<
class T>
double StatFactory::cross_correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int maxdelay, std::vector<T>& z)
const{
1317 double sumCorrelation=0;
1318 for (
int delay=-maxdelay;delay<maxdelay;delay++) {
1319 z.push_back(correlation(x,y,delay));
1320 sumCorrelation+=z.back();
1322 return sumCorrelation;
1326 template<
class T>
double StatFactory::linear_regression(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const{
1327 if(x.size()!=y.size()){
1328 std::ostringstream s;
1329 s<<
"Error: x and y not equal in size";
1333 std::ostringstream s;
1334 s<<
"Error: x is empty";
1341 gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1342 return (1-sumsq/var(y)/(y.size()-1));
1346 template<
class T>
double StatFactory::linear_regression_err(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const{
1347 if(x.size()!=y.size()){
1348 std::ostringstream s;
1349 s<<
"Error: x and y not equal in size";
1353 std::ostringstream s;
1354 s<<
"Error: x is empty";
1361 gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1362 return sqrt((sumsq)/(y.size()));
1368 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{
1369 if(wavelengthIn.empty()){
1370 std::ostringstream s;
1371 s<<
"Error: wavelengthIn is empty";
1374 std::vector<double> wavelengthOut=wavelengthIn;
1375 std::vector<T> validIn=input;
1376 if(input.size()!=wavelengthIn.size()){
1377 std::ostringstream s;
1378 s<<
"Error: x and y not equal in size";
1381 int nband=wavelengthIn.size();
1384 if(m_noDataValues.size()){
1385 typename std::vector<T>::iterator itValue=validIn.begin();
1386 typename std::vector<T>::iterator itWavelength=wavelengthOut.begin();
1387 while(itValue!=validIn.end()&&itWavelength!=wavelengthOut.end()){
1388 if(isNoData(*itValue)){
1389 validIn.erase(itValue);
1390 wavelengthOut.erase(itWavelength);
1397 if(validIn.size()>1){
1399 interpolateUp(wavelengthOut, validIn, wavelengthIn, type, output, verbose);
1412 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{
1413 if(wavelengthIn.empty()){
1414 std::ostringstream s;
1415 s<<
"Error: wavelengthIn is empty";
1418 if(input.size()!=wavelengthIn.size()){
1419 std::ostringstream s;
1420 s<<
"Error: input and wavelengthIn not equal in size";
1423 if(wavelengthOut.empty()){
1424 std::ostringstream s;
1425 s<<
"Error: wavelengthOut is empty";
1428 int nband=wavelengthIn.size();
1430 gsl_interp_accel *acc;
1433 getSpline(type,nband,spline);
1435 assert(&(wavelengthIn[0]));
1436 assert(&(input[0]));
1437 int status=initSpline(spline,&(wavelengthIn[0]),&(input[0]),nband);
1439 std::string errorString=
"Could not initialize spline";
1442 for(
int index=0;index<wavelengthOut.size();++index){
1443 if(wavelengthOut[index]<*wavelengthIn.begin()){
1444 output.push_back(*(input.begin()));
1447 else if(wavelengthOut[index]>wavelengthIn.back()){
1448 output.push_back(input.back());
1451 double dout=evalSpline(spline,wavelengthOut[index],acc);
1452 output.push_back(dout);
1454 gsl_spline_free(spline);
1455 gsl_interp_accel_free(acc);
1490 template<
class T>
void StatFactory::interpolateUp(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const
1493 std::ostringstream s;
1494 s<<
"Error: input is empty";
1498 std::ostringstream s;
1499 s<<
"Error: nbin must be larger than 0";
1503 int dim=input.size();
1504 for(
int i=0;i<dim;++i){
1506 double left=input[i];
1508 double right=(i<dim-1)? input[i+1]:input[i];
1509 deltaX=(right-left)/static_cast<double>(nbin);
1510 for(
int x=0;x<nbin;++x){
1511 output.push_back(left+x*deltaX);
1515 output.push_back(input.back());
1520 template<
class T>
void StatFactory::nearUp(
const std::vector<T>& input, std::vector<T>& output)
const
1523 std::ostringstream s;
1524 s<<
"Error: input is empty";
1527 if(output.size()<input.size()){
1528 std::ostringstream s;
1529 s<<
"Error: output size is smaller than input size";
1532 int dimInput=input.size();
1533 int dimOutput=output.size();
1535 for(
int iin=0;iin<dimInput;++iin){
1536 for(
int iout=0;iout<dimOutput/dimInput;++iout){
1537 int indexOutput=iin*dimOutput/dimInput+iout;
1538 if(indexOutput>=output.size()){
1539 std::ostringstream s;
1540 s<<
"Error: indexOutput must be smaller than output.size()";
1543 output[indexOutput]=input[iin];
1549 template<
class T>
void StatFactory::interpolateUp(
double* input,
int dim, std::vector<T>& output,
int nbin)
1552 std::ostringstream s;
1553 s<<
"Error: nbin must be larger than 0";
1557 for(
int i=0;i<dim;++i){
1559 double left=input[i];
1561 double right=(i<dim-1)? input[i+1]:input[i];
1562 deltaX=(right-left)/static_cast<double>(nbin);
1563 for(
int x=0;x<nbin;++x){
1564 output.push_back(left+x*deltaX);
1568 output.push_back(input[dim-1]);
1573 template<
class T>
void StatFactory::interpolateDown(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const
1576 std::ostringstream s;
1577 s<<
"Error: input is empty";
1581 std::ostringstream s;
1582 s<<
"Error: nbin must be larger than 0";
1586 int dim=input.size();
1588 output.push_back(input[0]);
1589 for(
int i=1;i<dim;++i){
1594 output.push_back(input[i]);
1600 template<
class T>
void StatFactory::interpolateDown(
double* input,
int dim, std::vector<T>& output,
int nbin)
1603 std::ostringstream s;
1604 s<<
"Error: nbin must be larger than 0";
1609 output.push_back(input[0]);
1610 for(
int i=1;i<dim;++i){
1615 output.push_back(input[i]);