pktools  2.6.5
Processing Kernel for geospatial data
StatFactory.h
1 /**********************************************************************
2 StatFactory.h: class for statistical operations on vectors
3 Copyright (C) 2008-2013 Pieter Kempeneers
4 
5 This file is part of pktools
6 
7 pktools is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11 
12 pktools is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16 
17 You should have received a copy of the GNU General Public License
18 along with pktools. If not, see <http://www.gnu.org/licenses/>.
19 ***********************************************************************/
20 #ifndef _STATFACTORY_H_
21 #define _STATFACTORY_H_
22 
23 #include <iostream>
24 #include <vector>
25 #include <map>
26 #include <math.h>
27 #include <assert.h>
28 #include <string>
29 #include <sstream>
30 #include <fstream>
31 #include <algorithm>
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>
39 
40 namespace statfactory
41 {
42 
44 
45 public:
46  enum INTERPOLATION_TYPE {undefined=0,linear=1,polynomial=2,cspline=3,cspline_periodic=4,akima=5,akima_periodic=6};
47  //todo: expand with other distributions as in http://www.gnu.org/software/gsl/manual/gsl-ref_toc.html#TOC320
48  enum DISTRIBUTION_TYPE {uniform=1,gaussian=2};
49 
50  StatFactory(void){};
51  virtual ~StatFactory(void){};
52  INTERPOLATION_TYPE getInterpolationType(const std::string interpolationType){
53  std::map<std::string, INTERPOLATION_TYPE> m_interpMap;
54  initMap(m_interpMap);
55  return m_interpMap[interpolationType];
56  };
57  DISTRIBUTION_TYPE getDistributionType(const std::string distributionType){
58  std::map<std::string, DISTRIBUTION_TYPE> m_distMap;
59  initDist(m_distMap);
60  return m_distMap[distributionType];
61  };
62  static void allocAcc(gsl_interp_accel *&acc){
63  acc = gsl_interp_accel_alloc ();
64  };
65 
66  static void getSpline(const std::string type, int size, gsl_spline *& spline){
67  std::map<std::string, INTERPOLATION_TYPE> m_interpMap;
68  initMap(m_interpMap);
69  switch(m_interpMap[type]){
70  case(polynomial):
71  spline=gsl_spline_alloc(gsl_interp_polynomial,size);
72  break;
73  case(cspline):
74  spline=gsl_spline_alloc(gsl_interp_cspline,size);
75  break;
76  case(cspline_periodic):
77  spline=gsl_spline_alloc(gsl_interp_cspline_periodic,size);
78  break;
79  case(akima):
80  spline=gsl_spline_alloc(gsl_interp_akima,size);
81  break;
82  case(akima_periodic):
83  spline=gsl_spline_alloc(gsl_interp_akima_periodic,size);
84  break;
85  case(linear):
86  default:
87  spline=gsl_spline_alloc(gsl_interp_linear,size);
88  break;
89  }
90  assert(spline);
91  };
92  static int initSpline(gsl_spline *spline, const double *x, const double *y, int size){
93  return gsl_spline_init (spline, x, y, size);
94  };
95  static double evalSpline(gsl_spline *spline, double x, gsl_interp_accel *acc){
96  return gsl_spline_eval (spline, x, acc);
97  };
98 
99  static gsl_rng* getRandomGenerator(unsigned long int theSeed){
100  const gsl_rng_type * T;
101  gsl_rng * r;
102 
103  // select random number generator
104  r = gsl_rng_alloc (gsl_rng_mt19937);
105  gsl_rng_set(r,theSeed);
106  return r;
107  }
108  void getNodataValues(std::vector<double>& nodatav) const{nodatav=m_noDataValues;};
109  bool isNoData(double value) const{
110  if(m_noDataValues.empty())
111  return false;
112  else
113  return find(m_noDataValues.begin(),m_noDataValues.end(),value)!=m_noDataValues.end();
114  };
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();
119  };
120  unsigned int setNoDataValues(std::vector<double> vnodata){
121  m_noDataValues=vnodata;
122  return m_noDataValues.size();
123  };
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;
126  initDist(m_distMap);
127  double randValue=0;
128  switch(m_distMap[type]){
129  case(uniform):
130  randValue = a+gsl_rng_uniform(r);
131  break;
132  case(gaussian):
133  default:
134  randValue = a+gsl_ran_gaussian(r, b);
135  break;
136  }
137  return randValue;
138  };
139 
140 
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;
145 // 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;
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;
156 
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;
184  // template<class T> double gsl_correlation(const std::vector<T>& x, const std::vector<T>& y) 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;
192  // template<class T> void interpolateUp(const std::vector< std::vector<T> >& input, std::vector< std::vector<T> >& output, double start, double end, double step, const gsl_interp_type* type);
193  // template<class T> void interpolateUp(const std::vector< std::vector<T> >& input, const std::vector<double>& wavelengthIn, std::vector< std::vector<T> >& output, std::vector<double>& wavelengthOut, double start, double end, double step, const gsl_interp_type* type);
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);
199 
200 private:
201  static void initMap(std::map<std::string, INTERPOLATION_TYPE>& m_interpMap){
202  //initialize selMap
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;
209  }
210  static void initDist(std::map<std::string, DISTRIBUTION_TYPE>& m_distMap){
211  //initialize distMap
212  m_distMap["gaussian"]=gaussian;
213  m_distMap["uniform"]=uniform;
214  }
215  std::vector<double> m_noDataValues;
216 };
217 
218 
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
220 {
221  bool isValid=false;
222  typename std::vector<T>::const_iterator tmpIt;
223  for(typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
224  if(!isNoData(*it)){
225  if(isValid){
226  if(*tmpIt>*it)
227  tmpIt=it;
228  }
229  else{
230  tmpIt=it;
231  isValid=true;
232  }
233  }
234  }
235  if(isValid)
236  return tmpIt;
237  else if(m_noDataValues.size())
238  return m_noDataValues[0];
239  else{
240  std::string errorString="Error: no valid data found";
241  throw(errorString);
242  }
243 }
244 
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
246 {
247  bool isValid=false;
248  typename std::vector<T>::iterator tmpIt;
249  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
250  if(!isNoData(*it)){
251  if(isValid){
252  if(*tmpIt>*it)
253  tmpIt=it;
254  }
255  else{
256  tmpIt=it;
257  isValid=true;
258  }
259  }
260  }
261  if(isValid)
262  return tmpIt;
263  else if(m_noDataValues.size())
264  return m_noDataValues[0];
265  else{
266  std::string errorString="Error: no valid data found";
267  throw(errorString);
268  }
269 }
270 
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
272 {
273  bool isValid=false;
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){
277  if(isNoData(*it))
278  continue;
279  if(isValid){
280  if((minConstraint<=*it)&&(*it<minValue)){
281  tmpIt=it;
282  minValue=*it;
283  }
284  }
285  else{
286  if(*it<minValue)
287  continue;
288  tmpIt=it;
289  minValue=*it;
290  isValid=true;
291  }
292  }
293  if(isValid)
294  return tmpIt;
295  else if(m_noDataValues.size())
296  return m_noDataValues[0];
297  else{
298  std::string errorString="Error: no valid data found";
299  throw(errorString);
300  }
301 }
302 
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
304 {
305  bool isValid=false;
306  typename std::vector<T>::iterator tmpIt;
307  T minValue=minConstraint;
308  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
309  if(isNoData(*it))
310  continue;
311  if(isValid){
312  if((minConstraint<=*it)&&(*it<minValue)){
313  tmpIt=it;
314  minValue=*it;
315  }
316  }
317  else{
318  if(*it<minConstraint)
319  continue;
320  tmpIt=it;
321  minValue=*it;
322  isValid=true;
323  }
324  }
325  if(isValid)
326  return tmpIt;
327  else if(m_noDataValues.size())
328  return m_noDataValues[0];
329  else{
330  std::string errorString="Error: no valid data found";
331  throw(errorString);
332  }
333 }
334 
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
336 {
337  bool isValid=false;
338  typename std::vector<T>::const_iterator tmpIt;
339  for (typename std::vector<T>::iterator it = begin; it!=end; ++it){
340  if(isNoData(*it))
341  continue;
342  if(isValid){
343  if(*tmpIt<*it)
344  tmpIt=it;
345  }
346  else{
347  tmpIt=it;
348  isValid=true;
349  }
350  }
351  if(isValid)
352  return tmpIt;
353  else if(m_noDataValues.size())
354  return m_noDataValues[0];
355  else{
356  std::string errorString="Error: no valid data found";
357  throw(errorString);
358  }
359 }
360 
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
362 {
363  bool isValid=false;
364  typename std::vector<T>::iterator tmpIt;
365  for (typename std::vector<T>::iterator it = begin; it!=end; ++it){
366  if(isNoData(*it))
367  continue;
368  if(isValid){
369  if(*tmpIt<*it)
370  tmpIt=it;
371  }
372  else{
373  tmpIt=it;
374  isValid=true;
375  }
376  }
377  if(isValid)
378  return tmpIt;
379  else
380  return end;
381 }
382 
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
384 {
385  bool isValid=false;
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){
389  if(isNoData(*it))
390  continue;
391  if(isValid){
392  if((maxConstraint>=*it)&&(*it>maxValue)){
393  tmpIt=it;
394  maxValue=*it;
395  }
396  }
397  else{
398  if(*it>maxConstraint)
399  continue;
400  tmpIt=it;
401  maxValue=*it;
402  isValid=true;
403  }
404  }
405  if(isValid)
406  return tmpIt;
407  else
408  return end;
409 }
410 
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
412 {
413  bool isValid=false;
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){
417  if(isNoData(*it))
418  continue;
419  if(isValid){
420  if((maxConstraint>=*it)&&(*it>maxValue)){
421  tmpIt=it;
422  maxValue=*it;
423  }
424  }
425  else{
426  if(*it>maxValue)
427  continue;
428  tmpIt=it;
429  maxValue=*it;
430  isValid=true;
431  }
432  }
433  if(isValid)
434  return tmpIt;
435  else
436  return end;
437 }
438 
439 template<class T> inline T StatFactory::mymin(const std::vector<T>& v) const
440 {
441  bool isValid=false;
442  if(v.empty()){
443  std::string errorString="Error: vector is empty";
444  throw(errorString);
445  }
446  T minValue;
447  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
448  if(isNoData(*it))
449  continue;
450  if(isValid){
451  if(minValue>*it)
452  minValue=*it;
453  }
454  else{
455  minValue=*it;
456  isValid=true;
457  }
458  }
459  if(isValid)
460  return minValue;
461  else if(m_noDataValues.size())
462  return m_noDataValues[0];
463  else{
464  std::string errorString="Error: no valid data found";
465  throw(errorString);
466  }
467 }
468 
469  template<class T> inline T StatFactory::mymin(const std::vector<T>& v, T minConstraint) const
470 {
471  bool isValid=false;
472  T minValue=minConstraint;
473  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
474  if(isNoData(*it))
475  continue;
476  if(isValid){
477  if((minConstraint<=*it)&&(*it<minValue))
478  minValue=*it;
479  }
480  else{
481  if(*it<minValue)
482  continue;
483  minValue=*it;
484  isValid=true;
485  }
486  }
487  if(isValid)
488  return minValue;
489  else if(m_noDataValues.size())
490  return m_noDataValues[0];
491  else{
492  std::string errorString="Error: no valid data found";
493  throw(errorString);
494  }
495 }
496 
497 template<class T> inline T StatFactory::mymax(const std::vector<T>& v) const
498 {
499  bool isValid=false;
500  if(v.empty()){
501  std::string errorString="Error: vector is empty";
502  throw(errorString);
503  }
504  T maxValue;
505  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
506  if(isNoData(*it))
507  continue;
508  if(isValid){
509  if(maxValue<*it)
510  maxValue=*it;
511  }
512  else{
513  maxValue=*it;
514  isValid=true;
515  }
516  }
517  if(isValid)
518  return maxValue;
519  else if(m_noDataValues.size())
520  return m_noDataValues[0];
521  else{
522  std::string errorString="Error: no valid data found";
523  throw(errorString);
524  }
525 }
526 
527 template<class T> inline T StatFactory::mymax(const std::vector<T>& v, T maxConstraint) const
528 {
529  bool isValid=false;
530  T maxValue=maxConstraint;
531  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
532  if(isNoData(*it))
533  continue;
534  if(isValid){
535  if((*it<=maxConstraint)&&(*it>maxValue))
536  maxValue=*it;
537  }
538  else{
539  if(*it>maxValue)
540  continue;
541  maxValue=*it;
542  isValid=true;
543  }
544  }
545  if(isValid)
546  return maxValue;
547  else if(m_noDataValues.size())
548  return m_noDataValues[0];
549  else{
550  std::string errorString="Error: no valid data found";
551  throw(errorString);
552  }
553 }
554 
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
556 {
557  bool isValid=false;
558  typename std::vector<T>::const_iterator tmpIt;
559  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
560  if(isNoData(*it))
561  continue;
562  if(isValid){
563  if(abs(*tmpIt)<abs(*it))
564  tmpIt=it;
565  }
566  else{
567  tmpIt=it;
568  isValid=true;
569  }
570  }
571  if(isValid)
572  return tmpIt;
573  else
574  return end;
575 }
576 
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
578 {
579  bool isValid=false;
580  typename std::vector<T>::const_iterator tmpIt;
581  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
582  if(isNoData(*it))
583  continue;
584  if(isValid){
585  if(abs(*tmpIt)>abs(*it))
586  tmpIt=it;
587  }
588  else{
589  tmpIt=it;
590  isValid=true;
591  }
592  }
593  if(isValid)
594  return tmpIt;
595  else
596  return end;
597 }
598 
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
600 {
601  bool isConstraint=(theMax>theMin);
602  double minConstraint=theMin;
603  double maxConstraint=theMax;
604  bool isValid=false;
605  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
606  if(isNoData(*it))
607  continue;
608  if(isValid){
609  if(isConstraint){
610  if(*it<minConstraint)
611  continue;
612  if(*it>maxConstraint)
613  continue;
614  }
615  if(*it<theMin)
616  theMin=*it;
617  if(*it>theMax)
618  theMax=*it;
619  }
620  else{
621  if(isConstraint){
622  if(*it<minConstraint)
623  continue;
624  if(*it>maxConstraint)
625  continue;
626  }
627  theMin=*it;
628  theMax=*it;
629  isValid=true;
630  }
631  }
632  if(!isValid){
633  if(m_noDataValues.size()){
634  theMin=m_noDataValues[0];
635  theMax=m_noDataValues[0];
636  }
637  else{
638  std::string errorString="Error: no valid data found";
639  throw(errorString);
640  }
641  }
642 }
643 
644 template<class T> inline T StatFactory::sum(const std::vector<T>& v) const
645 {
646  bool isValid=false;
647  typename std::vector<T>::const_iterator it;
648  T tmpSum=0;
649  for (it = v.begin(); it!= v.end(); ++it){
650  if(isNoData(*it))
651  continue;
652  isValid=true;
653  tmpSum+=*it;
654  }
655  if(isValid)
656  return tmpSum;
657  else if(m_noDataValues.size())
658  return m_noDataValues[0];
659  else{
660  std::string errorString="Error: no valid data found";
661  throw(errorString);
662  }
663 }
664 
665 template<class T> inline double StatFactory::mean(const std::vector<T>& v) const
666 {
667  typename std::vector<T>::const_iterator it;
668  T tmpSum=0;
669  unsigned int validSize=0;
670  for (it = v.begin(); it!= v.end(); ++it){
671  if(isNoData(*it))
672  continue;
673  ++validSize;
674  tmpSum+=*it;
675  }
676  if(validSize)
677  return static_cast<double>(tmpSum)/validSize;
678  else if(m_noDataValues.size())
679  return m_noDataValues[0];
680  else{
681  std::string errorString="Error: no valid data found";
682  throw(errorString);
683  }
684 }
685 
686 template<class T> inline void StatFactory::eraseNoData(std::vector<T>& v) const
687 {
688  if(m_noDataValues.size()){
689  typename std::vector<T>::iterator it=v.begin();
690  while(it!=v.end()){
691  if(isNoData(*it))
692  v.erase(it);
693  else
694  ++it;
695  }
696  }
697 }
698 
699  template<class T> unsigned int StatFactory::nvalid(const std::vector<T>& v) const{
700  std::vector<T> tmpV=v;
701  eraseNoData(tmpV);
702  return(tmpV.size());
703  }
704 
705 template<class T> T StatFactory::median(const std::vector<T>& v) const
706 {
707  std::vector<T> tmpV=v;
708  eraseNoData(tmpV);
709  if(tmpV.size()){
710  sort(tmpV.begin(),tmpV.end());
711  if(tmpV.size()%2)
712  return tmpV[tmpV.size()/2];
713  else
714  return 0.5*(tmpV[tmpV.size()/2-1]+tmpV[tmpV.size()/2]);
715  }
716  else if(m_noDataValues.size())
717  return m_noDataValues[0];
718  else{
719  std::string errorString="Error: no valid data found";
720  throw(errorString);
721  }
722 }
723 
724 template<class T> double StatFactory::var(const std::vector<T>& v) const
725 {
726  typename std::vector<T>::const_iterator it;
727  unsigned int validSize=0;
728  double m1=0;
729  double m2=0;
730  for (it = v.begin(); it!= v.end(); ++it){
731  if(isNoData(*it))
732  continue;
733  m1+=*it;
734  m2+=(*it)*(*it);
735  ++validSize;
736  }
737  if(validSize){
738  m2/=validSize;
739  m1/=validSize;
740  return m2-m1*m1;
741  }
742  else if(m_noDataValues.size())
743  return m_noDataValues[0];
744  else{
745  std::string errorString="Error: no valid data found";
746  throw(errorString);
747  }
748 }
749 
750 template<class T> double StatFactory::moment(const std::vector<T>& v, int n) const
751 {
752  unsigned int validSize=0;
753  typename std::vector<T>::const_iterator it;
754  double m=0;
755 // double m1=mean(v);
756  for(it = v.begin(); it!= v.end(); ++it){
757  if(isNoData(*it))
758  continue;
759  m+=pow((*it),n);
760  ++validSize;
761  }
762  if(validSize)
763  return m/validSize;
764  else if(m_noDataValues.size())
765  return m_noDataValues[0];
766  else{
767  std::string errorString="Error: no valid data found";
768  throw(errorString);
769  }
770 }
771 
772  //central moment
773 template<class T> double StatFactory::cmoment(const std::vector<T>& v, int n) const
774 {
775  unsigned int validSize=0;
776  typename std::vector<T>::const_iterator it;
777  double m=0;
778  double m1=mean(v);
779  for(it = v.begin(); it!= v.end(); ++it){
780  if(isNoData(*it))
781  continue;
782  m+=pow((*it-m1),n);
783  ++validSize;
784  }
785  if(validSize)
786  return m/validSize;
787  else if(m_noDataValues.size())
788  return m_noDataValues[0];
789  else{
790  std::string errorString="Error: no valid data found";
791  throw(errorString);
792  }
793 }
794 
795 template<class T> double StatFactory::skewness(const std::vector<T>& v) const
796 {
797  //todo: what if nodata value?
798  return cmoment(v,3)/pow(var(v),1.5);
799 }
800 
801 template<class T> double StatFactory::kurtosis(const std::vector<T>& v) const
802 {
803  //todo: what if nodata value?
804  double m2=cmoment(v,2);
805  double m4=cmoment(v,4);
806  return m4/m2/m2-3.0;
807 }
808 
809 template<class T> void StatFactory::meanVar(const std::vector<T>& v, double& m1, double& v1) const
810 {
811  typename std::vector<T>::const_iterator it;
812  unsigned int validSize=0;
813  m1=0;
814  v1=0;
815  double m2=0;
816  for (it = v.begin(); it!= v.end(); ++it){
817  if(isNoData(*it))
818  continue;
819  m1+=*it;
820  m2+=(*it)*(*it);
821  ++validSize;
822  }
823  if(validSize){
824  m2/=validSize;
825  m1/=validSize;
826  v1=m2-m1*m1;
827  }
828  else if(m_noDataValues.size()){
829  m1=m_noDataValues[0];
830  v1=m_noDataValues[0];
831  }
832  else{
833  std::string errorString="Error: no valid data found";
834  throw(errorString);
835  }
836 }
837 
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
839 {
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";
845  throw(errorString);
846  }
847  double scale=(ubound-lbound)/(maximum-minimum);
848  //todo: what if nodata value?
849  for (int i=0;i<input.size();++i){
850  output[i]=scale*(input[i]-(minimum))+lbound;
851  }
852 }
853 
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
855 {
856  double minValue=0;
857  double maxValue=0;
858  minmax(input,begin,end,minimum,maximum);
859  if(minimum<maximum&&minimum>minValue)
860  minValue=minimum;
861  if(minimum<maximum&&maximum<maxValue)
862  maxValue=maximum;
863 
864  //todo: check...
865  minimum=minValue;
866  maximum=maxValue;
867 
868  if(maximum<=minimum){
869  std::ostringstream s;
870  s<<"Error: could not calculate distribution (min>=max)";
871  throw(s.str());
872  }
873  if(!nbin){
874  std::string errorString="Error: nbin not defined";
875  throw(errorString);
876  }
877  if(!input.size()){
878  std::string errorString="Error: no valid data found";
879  throw(errorString);
880  }
881  if(output.size()!=nbin){
882  output.resize(nbin);
883  for(int i=0;i<nbin;output[i++]=0);
884  }
885  bool isValid=false;
886  typename std::vector<T>::const_iterator it;
887  for(it=begin;it!=end;++it){
888  if(*it<minimum)
889  continue;
890  if(*it>maximum)
891  continue;
892  if(isNoData(*it))
893  continue;
894  isValid=true;
895  if(sigma>0){
896  // minimum-=2*sigma;
897  // maximum+=2*sigma;
898  //create kde for Gaussian basis function
899  //todo: speed up by calculating first and last bin with non-zero contriubtion...
900  //todo: calculate real surface below pdf by using gsl_cdf_gaussian_P(x-mean+binsize,sigma)-gsl_cdf_gaussian_P(x-mean,sigma)
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;
905  }
906  }
907  else{
908  int theBin=0;
909  if(*it==maximum)
910  theBin=nbin-1;
911  else if(*it>minimum && *it<maximum)
912  theBin=static_cast<int>(static_cast<double>((nbin-1)*(*it)-minimum)/(maximum-minimum));
913  ++output[theBin];
914  // if(*it==maximum)
915  // ++output[nbin-1];
916  // else if(*it>=minimum && *it<maximum)
917  // ++output[static_cast<int>(static_cast<double>((*it)-minimum)/(maximum-minimum)*nbin)];
918  }
919  }
920  if(!isValid){
921  std::string errorString="Error: no valid data found";
922  throw(errorString);
923  }
924  else if(!filename.empty()){
925  std::ofstream outputfile;
926  outputfile.open(filename.c_str());
927  if(!outputfile){
928  std::ostringstream s;
929  s<<"Error opening distribution file , " << filename;
930  throw(s.str());
931  }
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;
934  outputfile.close();
935  }
936 }
937 
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
939 {
940  if(inputX.empty()){
941  std::ostringstream s;
942  s<<"Error: inputX is empty";
943  throw(s.str());
944  }
945  if(inputX.size()!=inputY.size()){
946  std::ostringstream s;
947  s<<"Error: inputX is empty";
948  throw(s.str());
949  }
950  unsigned long int npoint=inputX.size();
951  if(maxX<=minX)
952  minmax(inputX,inputX.begin(),inputX.end(),minX,maxX);
953  if(maxX<=minX){
954  std::ostringstream s;
955  s<<"Error: could not calculate distribution (minX>=maxX)";
956  throw(s.str());
957  }
958  if(maxY<=minY)
959  minmax(inputY,inputY.begin(),inputY.end(),minY,maxY);
960  if(maxY<=minY){
961  std::ostringstream s;
962  s<<"Error: could not calculate distribution (minY>=maxY)";
963  throw(s.str());
964  }
965  if(nbin<=1){
966  std::ostringstream s;
967  s<<"Error: nbin must be larger than 1";
968  throw(s.str());
969  }
970  output.resize(nbin);
971  for(int i=0;i<nbin;++i){
972  output[i].resize(nbin);
973  for(int j=0;j<nbin;++j)
974  output[i][j]=0;
975  }
976  int binX=0;
977  int binY=0;
978  for(unsigned long int ipoint=0;ipoint<npoint;++ipoint){
979  if(inputX[ipoint]==maxX)
980  binX=nbin-1;
981  else
982  binX=static_cast<int>(static_cast<double>(inputX[ipoint]-minX)/(maxX-minX)*nbin);
983  if(inputY[ipoint]==maxY)
984  binY=nbin-1;
985  else
986  binY=static_cast<int>(static_cast<double>(inputY[ipoint]-minY)/(maxY-minY)*nbin);
987  if(binX<0){
988  std::ostringstream s;
989  s<<"Error: binX is smaller than 0";
990  throw(s.str());
991  }
992  if(output.size()<=binX){
993  std::ostringstream s;
994  s<<"Error: output size must be larger than binX";
995  throw(s.str());
996  }
997  if(binY<0){
998  std::ostringstream s;
999  s<<"Error: binY is smaller than 0";
1000  throw(s.str());
1001  }
1002  if(output.size()<=binY){
1003  std::ostringstream s;
1004  s<<"Error: output size must be larger than binY";
1005  throw(s.str());
1006  }
1007  if(sigma>0){
1008  // minX-=2*sigma;
1009  // maxX+=2*sigma;
1010  // minY-=2*sigma;
1011  // maxY+=2*sigma;
1012  //create kde for Gaussian basis function
1013  //todo: speed up by calculating first and last bin with non-zero contriubtion...
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){
1018  //calculate \integral_ibinX^(ibinX+1)
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;
1022  }
1023  }
1024  }
1025  else
1026  ++output[binX][binY];
1027  }
1028  if(!filename.empty()){
1029  std::ofstream outputfile;
1030  outputfile.open(filename.c_str());
1031  if(!outputfile){
1032  std::ostringstream s;
1033  s<<"Error opening distribution file , " << filename;
1034  throw(s.str());
1035  }
1036  for(int binX=0;binX<nbin;++binX){
1037  outputfile << std::endl;
1038  for(int binY=0;binY<nbin;++binY){
1039  double binValueX=0;
1040  if(nbin==maxX-minX+1)
1041  binValueX=minX+binX;
1042  else
1043  binValueX=minX+static_cast<double>(maxX-minX)*(binX+0.5)/nbin;
1044  double binValueY=0;
1045  if(nbin==maxY-minY+1)
1046  binValueY=minY+binY;
1047  else
1048  binValueY=minY+static_cast<double>(maxY-minY)*(binY+0.5)/nbin;
1049  double value=0;
1050  value=static_cast<double>(output[binX][binY])/npoint;
1051  outputfile << binValueX << " " << binValueY << " " << value << std::endl;
1052  /* double value=static_cast<double>(output[binX][binY])/npoint; */
1053  /* outputfile << (maxX-minX)*bin/(nbin-1)+minX << " " << (maxY-minY)*bin/(nbin-1)+minY << " " << value << std::endl; */
1054  }
1055  }
1056  outputfile.close();
1057  }
1058 }
1059 
1060 //todo: what with nodata values?
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
1062 {
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";
1068  throw(s.str());
1069  }
1070  if(nbin<=1){
1071  std::ostringstream s;
1072  s<<"Error: nbin must be larger than 1";
1073  throw(s.str());
1074  }
1075  if(input.empty()){
1076  std::ostringstream s;
1077  s<<"Error: input is empty";
1078  throw(s.str());
1079  }
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);
1087  else
1088  ++vit;
1089  }
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){
1095  inputBin.clear();
1096  while(inputBin.size()<inputSort.size()/nbin&&vit!=inputSort.end()){
1097  inputBin.push_back(*vit);
1098  ++vit;
1099  }
1100  if(inputBin.size()){
1101  output[ibin]=mymax(inputBin);
1102  }
1103  }
1104  if(!filename.empty()){
1105  std::ofstream outputfile;
1106  outputfile.open(filename.c_str());
1107  if(!outputfile){
1108  std::ostringstream s;
1109  s<<"error opening distribution file , " << filename;
1110  throw(s.str());
1111  }
1112  for(int ibin=0;ibin<nbin;++ibin)
1113  outputfile << ibin*100.0/nbin << " " << static_cast<double>(output[ibin])/input.size() << std::endl;
1114  outputfile.close();
1115  }
1116 }
1117 
1118 //todo: what with nodata values?
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
1120 {
1121  if(input.empty()){
1122  std::ostringstream s;
1123  s<<"Error: input is empty";
1124  throw(s.str());
1125  }
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);
1133  }
1134  else
1135  ++vit;
1136  }
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);
1140 }
1141 
1142 template<class T> void StatFactory::signature(const std::vector<T>& input, double&k, double& alpha, double& beta, double e) const
1143 {
1144  double m1=moment(input,1);
1145  double m2=moment(input,2);
1146  signature(m1,m2,k,alpha,beta,e);
1147 }
1148 
1149 //todo: what with nodata values?
1150 template<class T> void StatFactory::normalize(const std::vector<T>& input, std::vector<double>& output) const{
1151  double total=sum(input);
1152  if(total){
1153  output.resize(input.size());
1154  for(int index=0;index<input.size();++index)
1155  output[index]=input[index]/total;
1156  }
1157  else
1158  output=input;
1159 }
1160 
1161 //todo: what with nodata values?
1162 template<class T> void StatFactory::normalize_pct(std::vector<T>& input) const{
1163  double total=sum(input);
1164  if(total){
1165  typename std::vector<T>::iterator it;
1166  for(it=input.begin();it!=input.end();++it)
1167  *it=100.0*(*it)/total;
1168  }
1169 }
1170 
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";
1175  throw(s.str());
1176  }
1177  if(x.empty()){
1178  std::ostringstream s;
1179  s<<"Error: x is empty";
1180  throw(s.str());
1181  }
1182  double mse=0;
1183  for(int isample=0;isample<x.size();++isample){
1184  if(isNoData(x[isample])||isNoData(y[isample]))
1185  continue;
1186  double e=x[isample]-y[isample];
1187  mse+=e*e/x.size();
1188  }
1189  return sqrt(mse);
1190 }
1191 
1192 //normalized root mean square error
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";
1197  throw(s.str());
1198  }
1199  if(x.empty()){
1200  std::ostringstream s;
1201  s<<"Error: x is empty";
1202  throw(s.str());
1203  }
1204  std::vector<T> tmpX=x;
1205  eraseNoData(tmpX);
1206  std::vector<T> tmpY=y;
1207  eraseNoData(tmpY);
1208  double maxY=mymax(y);
1209  double minY=mymin(y);
1210  double rangeY=maxY-minY;
1211  double mse=0;
1212  for(int isample=0;isample<x.size();++isample){
1213  double e=x[isample]-y[isample];
1214  mse+=e*e/x.size();
1215  }
1216  return sqrt(mse)/rangeY;
1217 }
1218 
1219 // coefficient of variation root mean square error
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";
1224  throw(s.str());
1225  }
1226  if(x.empty()){
1227  std::ostringstream s;
1228  s<<"Error: x is empty";
1229  throw(s.str());
1230  }
1231  std::vector<T> tmpX=x;
1232  eraseNoData(tmpX);
1233  std::vector<T> tmpY=y;
1234  eraseNoData(tmpY);
1235  double maxY=mymax(tmpY);
1236  double minY=mymin(tmpY);
1237  double rangeY=maxY-minY;
1238  double mse=0;
1239  for(int isample=0;isample<x.size();++isample){
1240  double e=x[isample]-y[isample];
1241  mse+=e*e/x.size();
1242  }
1243  return sqrt(mse)/mean(tmpY);
1244 }
1245 
1246 // template<class T> double StatFactory::gsl_correlation(const std::vector<T>& x, const std::vector<T>& y) const{
1247 // return(gsl_stats_correlation(&(x[0]),1,&(y[0]),1,x.size()));
1248 // }
1249 
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()));
1252  }
1253 
1254 
1255 template<class T> double StatFactory::correlation(const std::vector<T>& x, const std::vector<T>& y, int delay) const{
1256  double meanX=0;
1257  double meanY=0;
1258  double varX=0;
1259  double varY=0;
1260  double sXY=0;
1261  meanVar(x,meanX,varX);
1262  meanVar(y,meanY,varY);
1263  double denom = sqrt(varX*varY);
1264  bool isValid=false;
1265  if(denom){
1266  //Calculate the correlation series
1267  sXY = 0;
1268  for (int i=0;i<x.size();++i) {
1269  int j = i + delay;
1270  if (j < 0 || j >= y.size())
1271  continue;
1272  else if(isNoData(x[i])||isNoData(y[j]))
1273  continue;
1274  else{
1275  isValid=true;
1276  if(i<0){
1277  std::ostringstream s;
1278  s<<"Error: i must be positive";
1279  throw(s.str());
1280  }
1281  if(i>=x.size()){
1282  std::ostringstream s;
1283  s<<"Error: i must be smaller than x.size()";
1284  throw(s.str());
1285  }
1286  if(j<0){
1287  std::ostringstream s;
1288  s<<"Error: j must be positive";
1289  throw(s.str());
1290  }
1291  if(j>=y.size()){
1292  std::ostringstream s;
1293  s<<"Error: j must be smaller than y.size()";
1294  throw(s.str());
1295  }
1296  sXY += (x[i] - meanX) * (y[j] - meanY);
1297  }
1298  }
1299  if(isValid){
1300  double minSize=(x.size()<y.size())?x.size():y.size();
1301  return(sXY / denom / (minSize-1));
1302  }
1303  else if(m_noDataValues.size())
1304  return m_noDataValues[0];
1305  else{
1306  std::string errorString="Error: no valid data found";
1307  throw(errorString);
1308  }
1309  }
1310  else
1311  return 0;
1312 }
1313 
1314 //todo: what if no valid data?
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{
1316  z.clear();
1317  double sumCorrelation=0;
1318  for (int delay=-maxdelay;delay<maxdelay;delay++) {
1319  z.push_back(correlation(x,y,delay));
1320  sumCorrelation+=z.back();
1321  }
1322  return sumCorrelation;
1323 }
1324 
1325 //todo: nodata?
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";
1330  throw(s.str());
1331  }
1332  if(x.empty()){
1333  std::ostringstream s;
1334  s<<"Error: x is empty";
1335  throw(s.str());
1336  }
1337  double cov00;
1338  double cov01;
1339  double cov11;
1340  double sumsq;
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));
1343 }
1344 
1345 //todo: nodata?
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";
1350  throw(s.str());
1351  }
1352  if(x.empty()){
1353  std::ostringstream s;
1354  s<<"Error: x is empty";
1355  throw(s.str());
1356  }
1357  double cov00;
1358  double cov01;
1359  double cov11;
1360  double sumsq;
1361  gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1362  return sqrt((sumsq)/(y.size()));
1363 }
1364 
1365 //alternatively: use GNU scientific library:
1366 // gsl_stats_correlation (const double data1[], const size_t stride1, const double data2[], const size_t stride2, const size_t n)
1367 
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";
1372  throw(s.str());
1373  }
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";
1379  throw(s.str());
1380  }
1381  int nband=wavelengthIn.size();
1382  output.clear();
1383  //remove nodata from input and corresponding wavelengthIn
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);
1391  }
1392  else{
1393  ++itValue;
1394  ++itWavelength;
1395  }
1396  }
1397  if(validIn.size()>1){
1398  try{
1399  interpolateUp(wavelengthOut, validIn, wavelengthIn, type, output, verbose);
1400  }
1401  catch(...){
1402  output=input;
1403  }
1404  }
1405  else//we can not interpolate if no valid data
1406  output=input;
1407  }
1408  else//no nodata values to interpolate
1409  output=input;
1410 }
1411 
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";
1416  throw(s.str());
1417  }
1418  if(input.size()!=wavelengthIn.size()){
1419  std::ostringstream s;
1420  s<<"Error: input and wavelengthIn not equal in size";
1421  throw(s.str());
1422  }
1423  if(wavelengthOut.empty()){
1424  std::ostringstream s;
1425  s<<"Error: wavelengthOut is empty";
1426  throw(s.str());
1427  }
1428  int nband=wavelengthIn.size();
1429  output.clear();
1430  gsl_interp_accel *acc;
1431  allocAcc(acc);
1432  gsl_spline *spline;
1433  getSpline(type,nband,spline);
1434  assert(spline);
1435  assert(&(wavelengthIn[0]));
1436  assert(&(input[0]));
1437  int status=initSpline(spline,&(wavelengthIn[0]),&(input[0]),nband);
1438  if(status){
1439  std::string errorString="Could not initialize spline";
1440  throw(errorString);
1441  }
1442  for(int index=0;index<wavelengthOut.size();++index){
1443  if(wavelengthOut[index]<*wavelengthIn.begin()){
1444  output.push_back(*(input.begin()));
1445  continue;
1446  }
1447  else if(wavelengthOut[index]>wavelengthIn.back()){
1448  output.push_back(input.back());
1449  continue;
1450  }
1451  double dout=evalSpline(spline,wavelengthOut[index],acc);
1452  output.push_back(dout);
1453  }
1454  gsl_spline_free(spline);
1455  gsl_interp_accel_free(acc);
1456 }
1457 
1458 // template<class T> void StatFactory::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){
1459 // assert(wavelengthIn.size());
1460 // assert(wavelengthOut.size());
1461 // int nsample=input.size();
1462 // int nband=wavelengthIn.size();
1463 // output.clear();
1464 // output.resize(nsample);
1465 // gsl_interp_accel *acc;
1466 // allocAcc(acc);
1467 // gsl_spline *spline;
1468 // getSpline(type,nband,spline);
1469 // for(int isample=0;isample<nsample;++isample){
1470 // assert(input[isample].size()==wavelengthIn.size());
1471 // initSpline(spline,&(wavelengthIn[0]),&(input[isample][0]),nband);
1472 // for(int index=0;index<wavelengthOut.size();++index){
1473 // if(type=="linear"){
1474 // if(wavelengthOut[index]<wavelengthIn.back())
1475 // output[isample].push_back(*(input.begin()));
1476 // else if(wavelengthOut[index]>wavelengthIn.back())
1477 // output[isample].push_back(input.back());
1478 // }
1479 // else{
1480 // double dout=evalSpline(spline,wavelengthOut[index],acc);
1481 // output[isample].push_back(dout);
1482 // }
1483 // }
1484 // }
1485 // gsl_spline_free(spline);
1486 // gsl_interp_accel_free(acc);
1487 // }
1488 
1489 //todo: nodata?
1490 template<class T> void StatFactory::interpolateUp(const std::vector<T>& input, std::vector<T>& output, int nbin) const
1491 {
1492  if(input.empty()){
1493  std::ostringstream s;
1494  s<<"Error: input is empty";
1495  throw(s.str());
1496  }
1497  if(!nbin){
1498  std::ostringstream s;
1499  s<<"Error: nbin must be larger than 0";
1500  throw(s.str());
1501  }
1502  output.clear();
1503  int dim=input.size();
1504  for(int i=0;i<dim;++i){
1505  double deltaX=0;
1506  double left=input[i];
1507  if(i<dim-1){
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);
1512  }
1513  }
1514  else
1515  output.push_back(input.back());
1516  }
1517 }
1518 
1519 //todo: nodata?
1520 template<class T> void StatFactory::nearUp(const std::vector<T>& input, std::vector<T>& output) const
1521 {
1522  if(input.empty()){
1523  std::ostringstream s;
1524  s<<"Error: input is empty";
1525  throw(s.str());
1526  }
1527  if(output.size()<input.size()){
1528  std::ostringstream s;
1529  s<<"Error: output size is smaller than input size";
1530  throw(s.str());
1531  }
1532  int dimInput=input.size();
1533  int dimOutput=output.size();
1534 
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()";
1541  throw(s.str());
1542  }
1543  output[indexOutput]=input[iin];
1544  }
1545  }
1546 }
1547 
1548 //todo: nodata?
1549 template<class T> void StatFactory::interpolateUp(double* input, int dim, std::vector<T>& output, int nbin)
1550 {
1551  if(!nbin){
1552  std::ostringstream s;
1553  s<<"Error: nbin must be larger than 0";
1554  throw(s.str());
1555  }
1556  output.clear();
1557  for(int i=0;i<dim;++i){
1558  double deltaX=0;
1559  double left=input[i];
1560  if(i<dim-1){
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);
1565  }
1566  }
1567  else
1568  output.push_back(input[dim-1]);
1569  }
1570 }
1571 
1572 //todo: nodata?
1573 template<class T> void StatFactory::interpolateDown(const std::vector<T>& input, std::vector<T>& output, int nbin) const
1574 {
1575  if(input.empty()){
1576  std::ostringstream s;
1577  s<<"Error: input is empty";
1578  throw(s.str());
1579  }
1580  if(!nbin){
1581  std::ostringstream s;
1582  s<<"Error: nbin must be larger than 0";
1583  throw(s.str());
1584  }
1585  output.clear();
1586  int dim=input.size();
1587  int x=0;
1588  output.push_back(input[0]);
1589  for(int i=1;i<dim;++i){
1590  if(i%nbin)
1591  continue;
1592  else{
1593  x=(i-1)/nbin+1;
1594  output.push_back(input[i]);
1595  }
1596  }
1597 }
1598 
1599 //todo: nodata?
1600 template<class T> void StatFactory::interpolateDown(double* input, int dim, std::vector<T>& output, int nbin)
1601 {
1602  if(!nbin){
1603  std::ostringstream s;
1604  s<<"Error: nbin must be larger than 0";
1605  throw(s.str());
1606  }
1607  output.clear();
1608  int x=0;
1609  output.push_back(input[0]);
1610  for(int i=1;i<dim;++i){
1611  if(i%nbin)
1612  continue;
1613  else{
1614  x=(i-1)/nbin+1;
1615  output.push_back(input[i]);
1616  }
1617  }
1618 }
1619 }
1620 
1621 #endif /* _STATFACTORY_H_ */
1622 
1623 // void Histogram::signature(double m1, double m2, double& k, double& alpha, double& beta, double e)
1624 // {
1625 // double y=m1*m1/m2;
1626 // beta=F_1(y,0.1,10.0,e);
1627 // double fb=F(beta);
1628 // double g=exp(lgamma(1.0/beta));
1629 // alpha=m1*g/exp(lgamma(2.0/beta));
1630 // k=beta/(2*alpha*g);
1631 // // std::cout << "y, alpha, beta: " << y << ", " << alpha << ", " << beta << std::endl;
1632 // }
1633 
1634 // double Histogram::F(double x)
1635 // {
1636 // double g2=exp(lgamma(2.0/x));
1637 // return(g2*g2/exp(lgamma(3.0/x))/exp(lgamma(1.0/x)));
1638 // }
1639 
1640 // //x1 is under estimate, x2 is over estimate, e is error
1641 // double Histogram::F_1(double y, double x1, double x2, double e)
1642 // {
1643 // double f1=F(x1);
1644 // double f2=F(x2);
1645 // assert(f1!=f2);
1646 // double x=x1+(x2-x1)*(y-f1)/(f2-f1);
1647 // double f=F(x);
1648 // while(f-y>=e||y-f>=e){
1649 // if(f<y)
1650 // x1=x;
1651 // else
1652 // x2=x;
1653 // if(x1==x2)
1654 // return x1;
1655 // assert(f1!=f2);
1656 // x=x1+(x2-x1)*(y-f1)/(f2-f1);
1657 // f=F(x);
1658 // }
1659 // return x;
1660 // }