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;
220 for(
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
234 else if(m_noDataValues.size())
235 return m_noDataValues[0];
237 std::string errorString=
"Error: no valid data found";
242 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
245 typename std::vector<T>::iterator tmpIt;
246 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
260 else if(m_noDataValues.size())
261 return m_noDataValues[0];
263 std::string errorString=
"Error: no valid data found";
268 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
271 typename std::vector<T>::const_iterator tmpIt;
272 T minValue=minConstraint;
273 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
277 if((minConstraint<=*it)&&(*it<=minValue)){
290 else if(m_noDataValues.size())
291 return m_noDataValues[0];
293 std::string errorString=
"Error: no valid data found";
298 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
301 typename std::vector<T>::iterator tmpIt;
302 T minValue=minConstraint;
303 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
307 if((minConstraint<=*it)&&(*it<=minValue)){
320 else if(m_noDataValues.size())
321 return m_noDataValues[0];
323 std::string errorString=
"Error: no valid data found";
328 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
331 typename std::vector<T>::const_iterator tmpIt;
332 for (
typename std::vector<T>::iterator it = begin; it!=end; ++it){
346 else if(m_noDataValues.size())
347 return m_noDataValues[0];
349 std::string errorString=
"Error: no valid data found";
354 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
357 typename std::vector<T>::iterator tmpIt;
358 for (
typename std::vector<T>::iterator it = begin; it!=end; ++it){
376 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
379 typename std::vector<T>::const_iterator tmpIt;
380 T maxValue=maxConstraint;
381 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
385 if((maxConstraint>=*it)&&(*it>=maxValue)){
402 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
405 typename std::vector<T>::iterator tmpIt=v.end();
406 T maxValue=maxConstraint;
407 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
411 if((maxConstraint>=*it)&&(*it>=maxValue)){
428 template<
class T>
inline T StatFactory::mymin(
const std::vector<T>& v)
const
432 std::string errorString=
"Error: vector is empty";
436 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
450 else if(m_noDataValues.size())
451 return m_noDataValues[0];
453 std::string errorString=
"Error: no valid data found";
458 template<
class T>
inline T StatFactory::mymin(
const std::vector<T>& v, T minConstraint)
const
461 T minValue=minConstraint;
462 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
466 if((minConstraint<=*it)&&(*it<=minValue))
471 else if(m_noDataValues.size())
472 return m_noDataValues[0];
474 std::string errorString=
"Error: no valid data found";
479 template<
class T>
inline T StatFactory::mymax(
const std::vector<T>& v)
const
483 std::string errorString=
"Error: vector is empty";
487 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
501 else if(m_noDataValues.size())
502 return m_noDataValues[0];
504 std::string errorString=
"Error: no valid data found";
509 template<
class T>
inline T StatFactory::mymax(
const std::vector<T>& v, T maxConstraint)
const
512 T maxValue=maxConstraint;
513 for (
typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
516 if((maxValue<=*it)&&(*it<=maxConstraint))
521 else if(m_noDataValues.size())
522 return m_noDataValues[0];
524 std::string errorString=
"Error: no valid data found";
529 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
532 typename std::vector<T>::const_iterator tmpIt;
533 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
537 if(abs(*tmpIt)<abs(*it))
551 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
554 typename std::vector<T>::const_iterator tmpIt;
555 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
559 if(abs(*tmpIt)>abs(*it))
573 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
576 for (
typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
592 if(m_noDataValues.size()){
593 theMin=m_noDataValues[0];
594 theMax=m_noDataValues[0];
597 std::string errorString=
"Error: no valid data found";
603 template<
class T>
inline T StatFactory::sum(
const std::vector<T>& v)
const
606 typename std::vector<T>::const_iterator it;
608 for (it = v.begin(); it!= v.end(); ++it){
616 else if(m_noDataValues.size())
617 return m_noDataValues[0];
619 std::string errorString=
"Error: no valid data found";
624 template<
class T>
inline double StatFactory::mean(
const std::vector<T>& v)
const
626 typename std::vector<T>::const_iterator it;
628 unsigned int validSize=0;
629 for (it = v.begin(); it!= v.end(); ++it){
636 return static_cast<double>(tmpSum)/validSize;
637 else if(m_noDataValues.size())
638 return m_noDataValues[0];
640 std::string errorString=
"Error: no valid data found";
645 template<
class T>
inline void StatFactory::eraseNoData(std::vector<T>& v)
const
647 if(m_noDataValues.size()){
648 typename std::vector<T>::iterator it=v.begin();
658 template<
class T> T StatFactory::median(
const std::vector<T>& v)
const
660 std::vector<T> tmpV=v;
663 sort(tmpV.begin(),tmpV.end());
665 return tmpV[tmpV.size()/2];
667 return 0.5*(tmpV[tmpV.size()/2-1]+tmpV[tmpV.size()/2]);
669 else if(m_noDataValues.size())
670 return m_noDataValues[0];
672 std::string errorString=
"Error: no valid data found";
677 template<
class T>
double StatFactory::var(
const std::vector<T>& v)
const
679 typename std::vector<T>::const_iterator it;
680 unsigned int validSize=0;
683 for (it = v.begin(); it!= v.end(); ++it){
695 else if(m_noDataValues.size())
696 return m_noDataValues[0];
698 std::string errorString=
"Error: no valid data found";
703 template<
class T>
double StatFactory::moment(
const std::vector<T>& v,
int n)
const
705 unsigned int validSize=0;
706 typename std::vector<T>::const_iterator it;
709 for(it = v.begin(); it!= v.end(); ++it){
717 else if(m_noDataValues.size())
718 return m_noDataValues[0];
720 std::string errorString=
"Error: no valid data found";
726 template<
class T>
double StatFactory::cmoment(
const std::vector<T>& v,
int n)
const
728 unsigned int validSize=0;
729 typename std::vector<T>::const_iterator it;
732 for(it = v.begin(); it!= v.end(); ++it){
740 else if(m_noDataValues.size())
741 return m_noDataValues[0];
743 std::string errorString=
"Error: no valid data found";
748 template<
class T>
double StatFactory::skewness(
const std::vector<T>& v)
const
751 return cmoment(v,3)/pow(var(v),1.5);
754 template<
class T>
double StatFactory::kurtosis(
const std::vector<T>& v)
const
757 double m2=cmoment(v,2);
758 double m4=cmoment(v,4);
762 template<
class T>
void StatFactory::meanVar(
const std::vector<T>& v,
double& m1,
double& v1)
const
764 typename std::vector<T>::const_iterator it;
765 unsigned int validSize=0;
769 for (it = v.begin(); it!= v.end(); ++it){
781 else if(m_noDataValues.size()){
782 m1=m_noDataValues[0];
783 v1=m_noDataValues[0];
786 std::string errorString=
"Error: no valid data found";
791 template<
class T1,
class T2>
void StatFactory::scale2byte(
const std::vector<T1>& input, std::vector<T2>& output,
unsigned char lbound,
unsigned char ubound)
const
793 output.resize(input.size());
794 T1 minimum=mymin(input);
795 T1 maximum=mymax(input);
796 if(minimum>=maximum){
797 std::string errorString=
"Error: no valid data found";
800 double scale=(ubound-lbound)/(maximum-minimum);
802 for (
int i=0;i<input.size();++i){
803 output[i]=scale*(input[i]-(minimum))+lbound;
807 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
811 minmax(input,begin,end,minValue,maxValue);
812 if(minimum<maximum&&minimum>minValue)
814 if(minimum<maximum&&maximum<maxValue)
821 if(maximum<=minimum){
822 std::ostringstream s;
823 s<<
"Error: could not calculate distribution (min>=max)";
827 std::string errorString=
"Error: nbin not defined";
831 std::string errorString=
"Error: no valid data found";
834 if(output.size()!=nbin){
836 for(
int i=0;i<nbin;output[i++]=0);
839 typename std::vector<T>::const_iterator it;
840 for(it=begin;it!=end;++it){
854 for(
int ibin=0;ibin<nbin;++ibin){
855 double icenter=minimum+
static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin;
856 double thePdf=gsl_ran_gaussian_pdf(*it-icenter, sigma);
857 output[ibin]+=thePdf;
864 else if(*it>minimum && *it<maximum)
865 theBin=
static_cast<int>(
static_cast<double>((nbin-1)*(*it)-minimum)/(maximum-minimum));
874 std::string errorString=
"Error: no valid data found";
877 else if(!filename.empty()){
878 std::ofstream outputfile;
879 outputfile.open(filename.c_str());
881 std::ostringstream s;
882 s<<
"Error opening distribution file , " << filename;
885 for(
int ibin=0;ibin<nbin;++ibin)
886 outputfile << minimum+static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin <<
" " <<
static_cast<double>(output[ibin])/input.size() << std::endl;
891 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
894 std::ostringstream s;
895 s<<
"Error: inputX is empty";
898 if(inputX.size()!=inputY.size()){
899 std::ostringstream s;
900 s<<
"Error: inputX is empty";
903 unsigned long int npoint=inputX.size();
905 minmax(inputX,inputX.begin(),inputX.end(),minX,maxX);
907 std::ostringstream s;
908 s<<
"Error: could not calculate distribution (minX>=maxX)";
912 minmax(inputY,inputY.begin(),inputY.end(),minY,maxY);
914 std::ostringstream s;
915 s<<
"Error: could not calculate distribution (minY>=maxY)";
919 std::ostringstream s;
920 s<<
"Error: nbin must be larger than 1";
924 for(
int i=0;i<nbin;++i){
925 output[i].resize(nbin);
926 for(
int j=0;j<nbin;++j)
931 for(
unsigned long int ipoint=0;ipoint<npoint;++ipoint){
932 if(inputX[ipoint]==maxX)
935 binX=
static_cast<int>(
static_cast<double>(inputX[ipoint]-minX)/(maxX-minX)*nbin);
936 if(inputY[ipoint]==maxY)
939 binY=
static_cast<int>(
static_cast<double>(inputY[ipoint]-minY)/(maxY-minY)*nbin);
941 std::ostringstream s;
942 s<<
"Error: binX is smaller than 0";
945 if(output.size()<=binX){
946 std::ostringstream s;
947 s<<
"Error: output size must be larger than binX";
951 std::ostringstream s;
952 s<<
"Error: binY is smaller than 0";
955 if(output.size()<=binY){
956 std::ostringstream s;
957 s<<
"Error: output size must be larger than binY";
967 for(
int ibinX=0;ibinX<nbin;++ibinX){
968 double centerX=minX+
static_cast<double>(maxX-minX)*ibinX/nbin;
969 double pdfX=gsl_ran_gaussian_pdf(inputX[ipoint]-centerX, sigma);
970 for(
int ibinY=0;ibinY<nbin;++ibinY){
972 double centerY=minY+
static_cast<double>(maxY-minY)*ibinY/nbin;
973 double pdfY=gsl_ran_gaussian_pdf(inputY[ipoint]-centerY, sigma);
974 output[ibinX][binY]+=pdfX*pdfY;
979 ++output[binX][binY];
981 if(!filename.empty()){
982 std::ofstream outputfile;
983 outputfile.open(filename.c_str());
985 std::ostringstream s;
986 s<<
"Error opening distribution file , " << filename;
989 for(
int binX=0;binX<nbin;++binX){
990 outputfile << std::endl;
991 for(
int binY=0;binY<nbin;++binY){
993 if(nbin==maxX-minX+1)
996 binValueX=minX+
static_cast<double>(maxX-minX)*(binX+0.5)/nbin;
998 if(nbin==maxY-minY+1)
1001 binValueY=minY+
static_cast<double>(maxY-minY)*(binY+0.5)/nbin;
1003 value=
static_cast<double>(output[binX][binY])/npoint;
1004 outputfile << binValueX <<
" " << binValueY <<
" " << value << std::endl;
1014 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
1016 if(maximum<=minimum)
1017 minmax(input,begin,end,minimum,maximum);
1018 if(maximum<=minimum){
1019 std::ostringstream s;
1020 s<<
"Error: maximum must be at least minimum";
1024 std::ostringstream s;
1025 s<<
"Error: nbin must be larger than 1";
1029 std::ostringstream s;
1030 s<<
"Error: input is empty";
1033 output.resize(nbin);
1034 std::vector<T> inputSort;
1035 inputSort.assign(begin,end);
1036 typename std::vector<T>::iterator vit=inputSort.begin();
1037 while(vit!=inputSort.end()){
1038 if(*vit<minimum||*vit>maximum)
1039 inputSort.erase(vit);
1043 std::sort(inputSort.begin(),inputSort.end());
1044 vit=inputSort.begin();
1045 std::vector<T> inputBin;
1046 for(
int ibin=0;ibin<nbin;++ibin){
1048 while(inputBin.size()<inputSort.size()/nbin&&vit!=inputSort.end()){
1049 inputBin.push_back(*vit);
1052 if(inputBin.size()){
1053 output[ibin]=mymax(inputBin);
1056 if(!filename.empty()){
1057 std::ofstream outputfile;
1058 outputfile.open(filename.c_str());
1060 std::ostringstream s;
1061 s<<
"error opening distribution file , " << filename;
1064 for(
int ibin=0;ibin<nbin;++ibin)
1065 outputfile << ibin*100.0/nbin <<
" " << static_cast<double>(output[ibin])/input.size() << std::endl;
1071 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
1074 std::ostringstream s;
1075 s<<
"Error: input is empty";
1078 std::vector<T> inputSort;
1079 inputSort.assign(begin,end);
1080 typename std::vector<T>::iterator vit=inputSort.begin();
1081 while(vit!=inputSort.end()){
1082 if(maximum>minimum){
1083 if(*vit<minimum||*vit>maximum)
1084 inputSort.erase(vit);
1089 std::sort(inputSort.begin(),inputSort.end());
1090 return gsl_stats_quantile_from_sorted_data(&(inputSort[0]),1,inputSort.size(),percent/100.0);
1093 template<
class T>
void StatFactory::signature(
const std::vector<T>& input,
double&k,
double& alpha,
double& beta,
double e)
const
1095 double m1=moment(input,1);
1096 double m2=moment(input,2);
1097 signature(m1,m2,k,alpha,beta,e);
1101 template<
class T>
void StatFactory::normalize(
const std::vector<T>& input, std::vector<double>& output)
const{
1102 double total=sum(input);
1104 output.resize(input.size());
1105 for(
int index=0;index<input.size();++index)
1106 output[index]=input[index]/total;
1113 template<
class T>
void StatFactory::normalize_pct(std::vector<T>& input)
const{
1114 double total=sum(input);
1116 typename std::vector<T>::iterator it;
1117 for(it=input.begin();it!=input.end();++it)
1118 *it=100.0*(*it)/total;
1122 template<
class T>
double StatFactory::rmse(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1123 if(x.size()!=y.size()){
1124 std::ostringstream s;
1125 s<<
"Error: x and y not equal in size";
1129 std::ostringstream s;
1130 s<<
"Error: x is empty";
1134 for(
int isample=0;isample<x.size();++isample){
1135 if(isNoData(x[isample])||isNoData(y[isample]))
1137 double e=x[isample]-y[isample];
1147 template<
class T>
double StatFactory::gsl_covariance(
const std::vector<T>& x,
const std::vector<T>& y)
const{
1148 return(gsl_stats_covariance(&(x[0]),1,&(y[0]),1,x.size()));
1152 template<
class T>
double StatFactory::correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int delay)
const{
1158 meanVar(x,meanX,varX);
1159 meanVar(y,meanY,varY);
1160 double denom = sqrt(varX*varY);
1165 for (
int i=0;i<x.size();++i) {
1167 if (j < 0 || j >= y.size())
1169 else if(isNoData(x[i])||isNoData(y[j]))
1174 std::ostringstream s;
1175 s<<
"Error: i must be positive";
1179 std::ostringstream s;
1180 s<<
"Error: i must be smaller than x.size()";
1184 std::ostringstream s;
1185 s<<
"Error: j must be positive";
1189 std::ostringstream s;
1190 s<<
"Error: j must be smaller than y.size()";
1193 sXY += (x[i] - meanX) * (y[j] - meanY);
1197 double minSize=(x.size()<y.size())?x.size():y.size();
1198 return(sXY / denom / (minSize-1));
1200 else if(m_noDataValues.size())
1201 return m_noDataValues[0];
1203 std::string errorString=
"Error: no valid data found";
1212 template<
class T>
double StatFactory::cross_correlation(
const std::vector<T>& x,
const std::vector<T>& y,
int maxdelay, std::vector<T>& z)
const{
1214 double sumCorrelation=0;
1215 for (
int delay=-maxdelay;delay<maxdelay;delay++) {
1216 z.push_back(correlation(x,y,delay));
1217 sumCorrelation+=z.back();
1219 return sumCorrelation;
1223 template<
class T>
double StatFactory::linear_regression(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const{
1224 if(x.size()!=y.size()){
1225 std::ostringstream s;
1226 s<<
"Error: x and y not equal in size";
1230 std::ostringstream s;
1231 s<<
"Error: x is empty";
1238 gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1239 return (1-sumsq/var(y)/(y.size()-1));
1243 template<
class T>
double StatFactory::linear_regression_err(
const std::vector<T>& x,
const std::vector<T>& y,
double &c0,
double &c1)
const{
1244 if(x.size()!=y.size()){
1245 std::ostringstream s;
1246 s<<
"Error: x and y not equal in size";
1250 std::ostringstream s;
1251 s<<
"Error: x is empty";
1258 gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1259 return sqrt((sumsq)/(y.size()));
1265 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{
1266 if(wavelengthIn.empty()){
1267 std::ostringstream s;
1268 s<<
"Error: wavelengthIn is empty";
1271 std::vector<double> wavelengthOut=wavelengthIn;
1272 std::vector<T> validIn=input;
1273 if(input.size()!=wavelengthIn.size()){
1274 std::ostringstream s;
1275 s<<
"Error: x and y not equal in size";
1278 int nband=wavelengthIn.size();
1281 if(m_noDataValues.size()){
1282 typename std::vector<T>::iterator itValue=validIn.begin();
1283 typename std::vector<T>::iterator itWavelength=wavelengthOut.begin();
1284 while(itValue!=validIn.end()&&itWavelength!=wavelengthOut.end()){
1285 if(isNoData(*itValue)){
1286 validIn.erase(itValue);
1287 wavelengthOut.erase(itWavelength);
1294 if(validIn.size()>1){
1296 interpolateUp(wavelengthOut, validIn, wavelengthIn, type, output, verbose);
1309 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{
1310 if(wavelengthIn.empty()){
1311 std::ostringstream s;
1312 s<<
"Error: wavelengthIn is empty";
1315 if(input.size()!=wavelengthIn.size()){
1316 std::ostringstream s;
1317 s<<
"Error: input and wavelengthIn not equal in size";
1320 if(wavelengthOut.empty()){
1321 std::ostringstream s;
1322 s<<
"Error: wavelengthOut is empty";
1325 int nband=wavelengthIn.size();
1327 gsl_interp_accel *acc;
1330 getSpline(type,nband,spline);
1332 assert(&(wavelengthIn[0]));
1333 assert(&(input[0]));
1334 int status=initSpline(spline,&(wavelengthIn[0]),&(input[0]),nband);
1336 std::string errorString=
"Could not initialize spline";
1339 for(
int index=0;index<wavelengthOut.size();++index){
1340 if(wavelengthOut[index]<*wavelengthIn.begin()){
1341 output.push_back(*(input.begin()));
1344 else if(wavelengthOut[index]>wavelengthIn.back()){
1345 output.push_back(input.back());
1348 double dout=evalSpline(spline,wavelengthOut[index],acc);
1349 output.push_back(dout);
1351 gsl_spline_free(spline);
1352 gsl_interp_accel_free(acc);
1387 template<
class T>
void StatFactory::interpolateUp(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const
1390 std::ostringstream s;
1391 s<<
"Error: input is empty";
1395 std::ostringstream s;
1396 s<<
"Error: nbin must be larger than 0";
1400 int dim=input.size();
1401 for(
int i=0;i<dim;++i){
1403 double left=input[i];
1405 double right=(i<dim-1)? input[i+1]:input[i];
1406 deltaX=(right-left)/static_cast<double>(nbin);
1407 for(
int x=0;x<nbin;++x){
1408 output.push_back(left+x*deltaX);
1412 output.push_back(input.back());
1417 template<
class T>
void StatFactory::nearUp(
const std::vector<T>& input, std::vector<T>& output)
const
1420 std::ostringstream s;
1421 s<<
"Error: input is empty";
1424 if(output.size()<input.size()){
1425 std::ostringstream s;
1426 s<<
"Error: output size is smaller than input size";
1429 int dimInput=input.size();
1430 int dimOutput=output.size();
1432 for(
int iin=0;iin<dimInput;++iin){
1433 for(
int iout=0;iout<dimOutput/dimInput;++iout){
1434 int indexOutput=iin*dimOutput/dimInput+iout;
1435 if(indexOutput>=output.size()){
1436 std::ostringstream s;
1437 s<<
"Error: indexOutput must be smaller than output.size()";
1440 output[indexOutput]=input[iin];
1446 template<
class T>
void StatFactory::interpolateUp(
double* input,
int dim, std::vector<T>& output,
int nbin)
1449 std::ostringstream s;
1450 s<<
"Error: nbin must be larger than 0";
1454 for(
int i=0;i<dim;++i){
1456 double left=input[i];
1458 double right=(i<dim-1)? input[i+1]:input[i];
1459 deltaX=(right-left)/static_cast<double>(nbin);
1460 for(
int x=0;x<nbin;++x){
1461 output.push_back(left+x*deltaX);
1465 output.push_back(input[dim-1]);
1470 template<
class T>
void StatFactory::interpolateDown(
const std::vector<T>& input, std::vector<T>& output,
int nbin)
const
1473 std::ostringstream s;
1474 s<<
"Error: input is empty";
1478 std::ostringstream s;
1479 s<<
"Error: nbin must be larger than 0";
1483 int dim=input.size();
1485 output.push_back(input[0]);
1486 for(
int i=1;i<dim;++i){
1491 output.push_back(input[i]);
1497 template<
class T>
void StatFactory::interpolateDown(
double* input,
int dim, std::vector<T>& output,
int nbin)
1500 std::ostringstream s;
1501 s<<
"Error: nbin must be larger than 0";
1506 output.push_back(input[0]);
1507 for(
int i=1;i<dim;++i){
1512 output.push_back(input[i]);