pktools  2.6.3
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> 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;
181  // template<class T> double gsl_correlation(const std::vector<T>& x, const std::vector<T>& y) 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;
189  // 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);
190  // 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);
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);
196 
197 private:
198  static void initMap(std::map<std::string, INTERPOLATION_TYPE>& m_interpMap){
199  //initialize selMap
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;
206  }
207  static void initDist(std::map<std::string, DISTRIBUTION_TYPE>& m_distMap){
208  //initialize distMap
209  m_distMap["gaussian"]=gaussian;
210  m_distMap["uniform"]=uniform;
211  }
212  std::vector<double> m_noDataValues;
213 };
214 
215 
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
217 {
218  bool isValid=false;
219  typename std::vector<T>::const_iterator tmpIt=begin;
220  for(typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
221  if(!isNoData(*it)){
222  isValid=true;
223  if(*tmpIt>*it)
224  tmpIt=it;
225  }
226  }
227  if(isValid)
228  return tmpIt;
229  else if(m_noDataValues.size())
230  return m_noDataValues[0];
231  else{
232  std::string errorString="Error: no valid data found";
233  throw(errorString);
234  }
235 }
236 
237 template<class T> inline typename std::vector<T>::iterator StatFactory::mymin(const std::vector<T>& v, typename std::vector<T>::iterator begin, typename std::vector<T>::iterator end) const
238 {
239  bool isValid=false;
240  typename std::vector<T>::iterator tmpIt=begin;
241  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
242  if(!isNoData(*it)){
243  isValid=true;
244  if(*tmpIt>*it)
245  tmpIt=it;
246  }
247  }
248  if(isValid)
249  return tmpIt;
250  else if(m_noDataValues.size())
251  return m_noDataValues[0];
252  else{
253  std::string errorString="Error: no valid data found";
254  throw(errorString);
255  }
256 }
257 
258 template<class T> inline typename std::vector<T>::const_iterator StatFactory::mymin(const std::vector<T>& v, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end, T minConstraint) const
259 {
260  bool isValid=false;
261  typename std::vector<T>::const_iterator tmpIt=v.end();
262  T minValue=minConstraint;
263  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
264  if(isNoData(*it))
265  continue;
266  isValid=true;
267  if((minConstraint<=*it)&&(*it<=minValue)){
268  tmpIt=it;
269  minValue=*it;
270  }
271  }
272  if(isValid)
273  return tmpIt;
274  else if(m_noDataValues.size())
275  return m_noDataValues[0];
276  else{
277  std::string errorString="Error: no valid data found";
278  throw(errorString);
279  }
280 }
281 
282 template<class T> inline typename std::vector<T>::iterator StatFactory::mymin(const std::vector<T>& v, typename std::vector<T>::iterator begin, typename std::vector<T>::iterator end, T minConstraint) const
283 {
284  bool isValid=false;
285  typename std::vector<T>::iterator tmpIt=v.end();
286  T minValue=minConstraint;
287  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
288  if(isNoData(*it))
289  continue;
290  isValid=true;
291  if((minConstraint<=*it)&&(*it<=minValue)){
292  tmpIt=it;
293  minValue=*it;
294  }
295  }
296  if(isValid)
297  return tmpIt;
298  else if(m_noDataValues.size())
299  return m_noDataValues[0];
300  else{
301  std::string errorString="Error: no valid data found";
302  throw(errorString);
303  }
304 }
305 
306 template<class T> inline typename std::vector<T>::const_iterator StatFactory::mymax(const std::vector<T>& v, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end) const
307 {
308  bool isValid=false;
309  typename std::vector<T>::const_iterator tmpIt=begin;
310  for (typename std::vector<T>::iterator it = begin; it!=end; ++it){
311  if(isNoData(*it))
312  continue;
313  isValid=true;
314  if(*tmpIt<*it)
315  tmpIt=it;
316  }
317  if(isValid)
318  return tmpIt;
319  else if(m_noDataValues.size())
320  return m_noDataValues[0];
321  else{
322  std::string errorString="Error: no valid data found";
323  throw(errorString);
324  }
325 }
326 
327 template<class T> inline typename std::vector<T>::iterator StatFactory::mymax(const std::vector<T>& v, typename std::vector<T>::iterator begin, typename std::vector<T>::iterator end) const
328 {
329  bool isValid=false;
330  typename std::vector<T>::iterator tmpIt=begin;
331  for (typename std::vector<T>::iterator it = begin; it!=end; ++it){
332  if(isNoData(*it))
333  continue;
334  isValid=true;
335  if(*tmpIt<*it)
336  tmpIt=it;
337  }
338  if(isValid)
339  return tmpIt;
340  else
341  return end;
342 }
343 
344 template<class T> inline typename std::vector<T>::const_iterator StatFactory::mymax(const std::vector<T>& v, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end, T maxConstraint) const
345 {
346  bool isValid=false;
347  typename std::vector<T>::const_iterator tmpIt=v.end();
348  T maxValue=maxConstraint;
349  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
350  if(isNoData(*it))
351  continue;
352  isValid=true;
353  if((maxValue<=*it)&&(*it<=maxConstraint)){
354  tmpIt=it;
355  maxValue=*it;
356  }
357  }
358  if(isValid)
359  return tmpIt;
360  else
361  return end;
362 }
363 
364 template<class T> inline typename std::vector<T>::iterator StatFactory::mymax(const std::vector<T>& v, typename std::vector<T>::iterator begin, typename std::vector<T>::iterator end, T maxConstraint) const
365 {
366  bool isValid=false;
367  typename std::vector<T>::iterator tmpIt=v.end();
368  T maxValue=maxConstraint;
369  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
370  if(isNoData(*it))
371  continue;
372  isValid=true;
373  if((maxValue<=*it)&&(*it<=maxConstraint)){
374  tmpIt=it;
375  maxValue=*it;
376  }
377  }
378  if(isValid)
379  return tmpIt;
380  else
381  return end;
382 }
383 
384 template<class T> inline T StatFactory::mymin(const std::vector<T>& v) const
385 {
386  bool isValid=false;
387  if(v.empty()){
388  std::string errorString="Error: vector is empty";
389  throw(errorString);
390  }
391  T minValue=*(v.begin());
392  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
393  if(isNoData(*it))
394  continue;
395  isValid=true;
396  if(minValue>*it)
397  minValue=*it;
398  }
399  if(isValid)
400  return minValue;
401  else if(m_noDataValues.size())
402  return m_noDataValues[0];
403  else{
404  std::string errorString="Error: no valid data found";
405  throw(errorString);
406  }
407 }
408 
409  template<class T> inline T StatFactory::mymin(const std::vector<T>& v, T minConstraint) const
410 {
411  bool isValid=false;
412  T minValue=minConstraint;
413  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
414  if(isNoData(*it))
415  continue;
416  isValid=true;
417  if((minConstraint<=*it)&&(*it<=minValue))
418  minValue=*it;
419  }
420  if(isValid)
421  return minValue;
422  else if(m_noDataValues.size())
423  return m_noDataValues[0];
424  else{
425  std::string errorString="Error: no valid data found";
426  throw(errorString);
427  }
428 }
429 
430 template<class T> inline T StatFactory::mymax(const std::vector<T>& v) const
431 {
432  bool isValid=false;
433  if(v.empty()){
434  std::string errorString="Error: vector is empty";
435  throw(errorString);
436  }
437  T maxValue=*(v.begin());
438  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
439  if(isNoData(*it))
440  continue;
441  isValid=true;
442  if(maxValue<*it)
443  maxValue=*it;
444  }
445  if(isValid)
446  return maxValue;
447  else if(m_noDataValues.size())
448  return m_noDataValues[0];
449  else{
450  std::string errorString="Error: no valid data found";
451  throw(errorString);
452  }
453 }
454 
455 template<class T> inline T StatFactory::mymax(const std::vector<T>& v, T maxConstraint) const
456 {
457  bool isValid=false;
458  T maxValue=maxConstraint;
459  for (typename std::vector<T>::const_iterator it = v.begin(); it!=v.end(); ++it){
460  if(isNoData(*it))
461  continue;
462  if((maxValue<=*it)&&(*it<=maxConstraint))
463  maxValue=*it;
464  }
465  if(isValid)
466  return maxValue;
467  else if(m_noDataValues.size())
468  return m_noDataValues[0];
469  else{
470  std::string errorString="Error: no valid data found";
471  throw(errorString);
472  }
473 }
474 
475 template<class T> inline typename std::vector<T>::const_iterator StatFactory::absmax(const std::vector<T>& v, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end) const
476 {
477  bool isValid=false;
478  typename std::vector<T>::const_iterator tmpIt=begin;
479  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
480  if(isNoData(*it))
481  continue;
482  isValid=true;
483  if(abs(*tmpIt)<abs(*it))
484  tmpIt=it;
485  }
486  if(isValid)
487  return tmpIt;
488  else
489  return end;
490 }
491 
492 template<class T> inline typename std::vector<T>::const_iterator StatFactory::absmin(const std::vector<T>& v, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end) const
493 {
494  bool isValid=false;
495  typename std::vector<T>::const_iterator tmpIt=begin;
496  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
497  if(isNoData(*it))
498  continue;
499  isValid=true;
500  if(abs(*tmpIt)>abs(*it))
501  tmpIt=it;
502  }
503  if(isValid)
504  return tmpIt;
505  else
506  return end;
507 }
508 
509 template<class T> inline void StatFactory::minmax(const std::vector<T>& v, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end, T& theMin, T& theMax) const
510 {
511  bool isValid=false;
512  theMin=*begin;
513  theMax=*begin;
514  for (typename std::vector<T>::const_iterator it = begin; it!=end; ++it){
515  if(isNoData(*it))
516  continue;
517  isValid=true;
518  if(theMin>*it)
519  theMin=*it;
520  if(theMax<*it)
521  theMax=*it;
522  }
523  if(!isValid){
524  if(m_noDataValues.size()){
525  theMin=m_noDataValues[0];
526  theMax=m_noDataValues[0];
527  }
528  else{
529  std::string errorString="Error: no valid data found";
530  throw(errorString);
531  }
532  }
533 }
534 
535 template<class T> inline T StatFactory::sum(const std::vector<T>& v) const
536 {
537  bool isValid=false;
538  typename std::vector<T>::const_iterator it;
539  T tmpSum=0;
540  for (it = v.begin(); it!= v.end(); ++it){
541  if(isNoData(*it))
542  continue;
543  isValid=true;
544  tmpSum+=*it;
545  }
546  if(isValid)
547  return tmpSum;
548  else if(m_noDataValues.size())
549  return m_noDataValues[0];
550  else{
551  std::string errorString="Error: no valid data found";
552  throw(errorString);
553  }
554 }
555 
556 template<class T> inline double StatFactory::mean(const std::vector<T>& v) const
557 {
558  typename std::vector<T>::const_iterator it;
559  T tmpSum=0;
560  unsigned int validSize=0;
561  for (it = v.begin(); it!= v.end(); ++it){
562  if(isNoData(*it))
563  continue;
564  ++validSize;
565  tmpSum+=*it;
566  }
567  if(validSize)
568  return static_cast<double>(tmpSum)/validSize;
569  else if(m_noDataValues.size())
570  return m_noDataValues[0];
571  else{
572  std::string errorString="Error: no valid data found";
573  throw(errorString);
574  }
575 }
576 
577 template<class T> inline void StatFactory::eraseNoData(std::vector<T>& v) const
578 {
579  if(m_noDataValues.size()){
580  typename std::vector<T>::iterator it=v.begin();
581  while(it!=v.end()){
582  if(isNoData(*it))
583  v.erase(it);
584  else
585  ++it;
586  }
587  }
588 }
589 
590 template<class T> T StatFactory::median(const std::vector<T>& v) const
591 {
592  std::vector<T> tmpV=v;
593  eraseNoData(tmpV);
594  if(tmpV.size()){
595  sort(tmpV.begin(),tmpV.end());
596  if(tmpV.size()%2)
597  return tmpV[tmpV.size()/2];
598  else
599  return 0.5*(tmpV[tmpV.size()/2-1]+tmpV[tmpV.size()/2]);
600  }
601  else if(m_noDataValues.size())
602  return m_noDataValues[0];
603  else{
604  std::string errorString="Error: no valid data found";
605  throw(errorString);
606  }
607 }
608 
609 template<class T> double StatFactory::var(const std::vector<T>& v) const
610 {
611  typename std::vector<T>::const_iterator it;
612  unsigned int validSize=0;
613  double m1=0;
614  double m2=0;
615  for (it = v.begin(); it!= v.end(); ++it){
616  if(isNoData(*it))
617  continue;
618  m1+=*it;
619  m2+=(*it)*(*it);
620  ++validSize;
621  }
622  if(validSize){
623  m2/=validSize;
624  m1/=validSize;
625  return m2-m1*m1;
626  }
627  else if(m_noDataValues.size())
628  return m_noDataValues[0];
629  else{
630  std::string errorString="Error: no valid data found";
631  throw(errorString);
632  }
633 }
634 
635 template<class T> double StatFactory::moment(const std::vector<T>& v, int n) const
636 {
637  unsigned int validSize=0;
638  typename std::vector<T>::const_iterator it;
639  double m=0;
640 // double m1=mean(v);
641  for(it = v.begin(); it!= v.end(); ++it){
642  if(isNoData(*it))
643  continue;
644  m+=pow((*it),n);
645  ++validSize;
646  }
647  if(validSize)
648  return m/validSize;
649  else if(m_noDataValues.size())
650  return m_noDataValues[0];
651  else{
652  std::string errorString="Error: no valid data found";
653  throw(errorString);
654  }
655 }
656 
657  //central moment
658 template<class T> double StatFactory::cmoment(const std::vector<T>& v, int n) const
659 {
660  unsigned int validSize=0;
661  typename std::vector<T>::const_iterator it;
662  double m=0;
663  double m1=mean(v);
664  for(it = v.begin(); it!= v.end(); ++it){
665  if(isNoData(*it))
666  continue;
667  m+=pow((*it-m1),n);
668  ++validSize;
669  }
670  if(validSize)
671  return m/validSize;
672  else if(m_noDataValues.size())
673  return m_noDataValues[0];
674  else{
675  std::string errorString="Error: no valid data found";
676  throw(errorString);
677  }
678 }
679 
680 template<class T> double StatFactory::skewness(const std::vector<T>& v) const
681 {
682  //todo: what if nodata value?
683  return cmoment(v,3)/pow(var(v),1.5);
684 }
685 
686 template<class T> double StatFactory::kurtosis(const std::vector<T>& v) const
687 {
688  //todo: what if nodata value?
689  double m2=cmoment(v,2);
690  double m4=cmoment(v,4);
691  return m4/m2/m2-3.0;
692 }
693 
694 template<class T> void StatFactory::meanVar(const std::vector<T>& v, double& m1, double& v1) const
695 {
696  typename std::vector<T>::const_iterator it;
697  unsigned int validSize=0;
698  m1=0;
699  v1=0;
700  double m2=0;
701  for (it = v.begin(); it!= v.end(); ++it){
702  if(isNoData(*it))
703  continue;
704  m1+=*it;
705  m2+=(*it)*(*it);
706  ++validSize;
707  }
708  if(validSize){
709  m2/=validSize;
710  m1/=validSize;
711  v1=m2-m1*m1;
712  }
713  else if(m_noDataValues.size()){
714  m1=m_noDataValues[0];
715  v1=m_noDataValues[0];
716  }
717  else{
718  std::string errorString="Error: no valid data found";
719  throw(errorString);
720  }
721 }
722 
723 template<class T1, class T2> void StatFactory::scale2byte(const std::vector<T1>& input, std::vector<T2>& output, unsigned char lbound, unsigned char ubound) const
724 {
725  output.resize(input.size());
726  T1 minimum=mymin(input);
727  T1 maximum=mymax(input);
728  if(minimum>=maximum){
729  std::string errorString="Error: no valid data found";
730  throw(errorString);
731  }
732  double scale=(ubound-lbound)/(maximum-minimum);
733  //todo: what if nodata value?
734  for (int i=0;i<input.size();++i){
735  output[i]=scale*(input[i]-(minimum))+lbound;
736  }
737 }
738 
739 template<class T> void StatFactory::distribution(const std::vector<T>& input, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end, std::vector<double>& output, int nbin, T &minimum, T &maximum, double sigma, const std::string &filename) const
740 {
741  double minValue=0;
742  double maxValue=0;
743  minmax(input,begin,end,minValue,maxValue);
744  if(minimum<maximum&&minimum>minValue)
745  minValue=minimum;
746  if(minimum<maximum&&maximum<maxValue)
747  maxValue=maximum;
748 
749  //todo: check...
750  minimum=minValue;
751  maximum=maxValue;
752 
753  if(maximum<=minimum){
754  std::ostringstream s;
755  s<<"Error: could not calculate distribution (min>=max)";
756  throw(s.str());
757  }
758  if(!nbin){
759  std::string errorString="Error: nbin not defined";
760  throw(errorString);
761  }
762  if(!input.size()){
763  std::string errorString="Error: no valid data found";
764  throw(errorString);
765  }
766  if(output.size()!=nbin){
767  output.resize(nbin);
768  for(int i=0;i<nbin;output[i++]=0);
769  }
770  bool isValid=false;
771  typename std::vector<T>::const_iterator it;
772  for(it=begin;it!=end;++it){
773  if(*it<minimum)
774  continue;
775  if(*it>maximum)
776  continue;
777  if(isNoData(*it))
778  continue;
779  isValid=true;
780  if(sigma>0){
781  // minimum-=2*sigma;
782  // maximum+=2*sigma;
783  //create kde for Gaussian basis function
784  //todo: speed up by calculating first and last bin with non-zero contriubtion...
785  //todo: calculate real surface below pdf by using gsl_cdf_gaussian_P(x-mean+binsize,sigma)-gsl_cdf_gaussian_P(x-mean,sigma)
786  for(int ibin=0;ibin<nbin;++ibin){
787  double icenter=minimum+static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin;
788  double thePdf=gsl_ran_gaussian_pdf(*it-icenter, sigma);
789  output[ibin]+=thePdf;
790  }
791  }
792  else{
793  int theBin=0;
794  if(*it==maximum)
795  theBin=nbin-1;
796  else if(*it>minimum && *it<maximum)
797  theBin=static_cast<int>(static_cast<double>((nbin-1)*(*it)-minimum)/(maximum-minimum));
798  ++output[theBin];
799  // if(*it==maximum)
800  // ++output[nbin-1];
801  // else if(*it>=minimum && *it<maximum)
802  // ++output[static_cast<int>(static_cast<double>((*it)-minimum)/(maximum-minimum)*nbin)];
803  }
804  }
805  if(!isValid){
806  std::string errorString="Error: no valid data found";
807  throw(errorString);
808  }
809  else if(!filename.empty()){
810  std::ofstream outputfile;
811  outputfile.open(filename.c_str());
812  if(!outputfile){
813  std::ostringstream s;
814  s<<"Error opening distribution file , " << filename;
815  throw(s.str());
816  }
817  for(int ibin=0;ibin<nbin;++ibin)
818  outputfile << minimum+static_cast<double>(maximum-minimum)*(ibin+0.5)/nbin << " " << static_cast<double>(output[ibin])/input.size() << std::endl;
819  outputfile.close();
820  }
821 }
822 
823 template<class T> void StatFactory::distribution2d(const std::vector<T>& inputX, const std::vector<T>& inputY, std::vector< std::vector<double> >& output, int nbin, T& minX, T& maxX, T& minY, T& maxY, double sigma, const std::string& filename) const
824 {
825  if(inputX.empty()){
826  std::ostringstream s;
827  s<<"Error: inputX is empty";
828  throw(s.str());
829  }
830  if(inputX.size()!=inputY.size()){
831  std::ostringstream s;
832  s<<"Error: inputX is empty";
833  throw(s.str());
834  }
835  unsigned long int npoint=inputX.size();
836  if(maxX<=minX)
837  minmax(inputX,inputX.begin(),inputX.end(),minX,maxX);
838  if(maxX<=minX){
839  std::ostringstream s;
840  s<<"Error: could not calculate distribution (minX>=maxX)";
841  throw(s.str());
842  }
843  if(maxY<=minY)
844  minmax(inputY,inputY.begin(),inputY.end(),minY,maxY);
845  if(maxY<=minY){
846  std::ostringstream s;
847  s<<"Error: could not calculate distribution (minY>=maxY)";
848  throw(s.str());
849  }
850  if(nbin<=1){
851  std::ostringstream s;
852  s<<"Error: nbin must be larger than 1";
853  throw(s.str());
854  }
855  output.resize(nbin);
856  for(int i=0;i<nbin;++i){
857  output[i].resize(nbin);
858  for(int j=0;j<nbin;++j)
859  output[i][j]=0;
860  }
861  int binX=0;
862  int binY=0;
863  for(unsigned long int ipoint=0;ipoint<npoint;++ipoint){
864  if(inputX[ipoint]==maxX)
865  binX=nbin-1;
866  else
867  binX=static_cast<int>(static_cast<double>(inputX[ipoint]-minX)/(maxX-minX)*nbin);
868  if(inputY[ipoint]==maxY)
869  binY=nbin-1;
870  else
871  binY=static_cast<int>(static_cast<double>(inputY[ipoint]-minY)/(maxY-minY)*nbin);
872  if(binX<0){
873  std::ostringstream s;
874  s<<"Error: binX is smaller than 0";
875  throw(s.str());
876  }
877  if(output.size()<=binX){
878  std::ostringstream s;
879  s<<"Error: output size must be larger than binX";
880  throw(s.str());
881  }
882  if(binY<0){
883  std::ostringstream s;
884  s<<"Error: binY is smaller than 0";
885  throw(s.str());
886  }
887  if(output.size()<=binY){
888  std::ostringstream s;
889  s<<"Error: output size must be larger than binY";
890  throw(s.str());
891  }
892  if(sigma>0){
893  // minX-=2*sigma;
894  // maxX+=2*sigma;
895  // minY-=2*sigma;
896  // maxY+=2*sigma;
897  //create kde for Gaussian basis function
898  //todo: speed up by calculating first and last bin with non-zero contriubtion...
899  for(int ibinX=0;ibinX<nbin;++ibinX){
900  double centerX=minX+static_cast<double>(maxX-minX)*ibinX/nbin;
901  double pdfX=gsl_ran_gaussian_pdf(inputX[ipoint]-centerX, sigma);
902  for(int ibinY=0;ibinY<nbin;++ibinY){
903  //calculate \integral_ibinX^(ibinX+1)
904  double centerY=minY+static_cast<double>(maxY-minY)*ibinY/nbin;
905  double pdfY=gsl_ran_gaussian_pdf(inputY[ipoint]-centerY, sigma);
906  output[ibinX][binY]+=pdfX*pdfY;
907  }
908  }
909  }
910  else
911  ++output[binX][binY];
912  }
913  if(!filename.empty()){
914  std::ofstream outputfile;
915  outputfile.open(filename.c_str());
916  if(!outputfile){
917  std::ostringstream s;
918  s<<"Error opening distribution file , " << filename;
919  throw(s.str());
920  }
921  for(int binX=0;binX<nbin;++binX){
922  outputfile << std::endl;
923  for(int binY=0;binY<nbin;++binY){
924  double binValueX=0;
925  if(nbin==maxX-minX+1)
926  binValueX=minX+binX;
927  else
928  binValueX=minX+static_cast<double>(maxX-minX)*(binX+0.5)/nbin;
929  double binValueY=0;
930  if(nbin==maxY-minY+1)
931  binValueY=minY+binY;
932  else
933  binValueY=minY+static_cast<double>(maxY-minY)*(binY+0.5)/nbin;
934  double value=0;
935  value=static_cast<double>(output[binX][binY])/npoint;
936  outputfile << binValueX << " " << binValueY << " " << value << std::endl;
937  /* double value=static_cast<double>(output[binX][binY])/npoint; */
938  /* outputfile << (maxX-minX)*bin/(nbin-1)+minX << " " << (maxY-minY)*bin/(nbin-1)+minY << " " << value << std::endl; */
939  }
940  }
941  outputfile.close();
942  }
943 }
944 
945 //todo: what with nodata values?
946 template<class T> void StatFactory::percentiles (const std::vector<T>& input, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end, std::vector<T>& output, int nbin, T &minimum, T &maximum, const std::string &filename) const
947 {
948  if(maximum<=minimum)
949  minmax(input,begin,end,minimum,maximum);
950  if(maximum<=minimum){
951  std::ostringstream s;
952  s<<"Error: maximum must be at least minimum";
953  throw(s.str());
954  }
955  if(nbin<=1){
956  std::ostringstream s;
957  s<<"Error: nbin must be larger than 1";
958  throw(s.str());
959  }
960  if(input.empty()){
961  std::ostringstream s;
962  s<<"Error: input is empty";
963  throw(s.str());
964  }
965  output.resize(nbin);
966  std::vector<T> inputSort;
967  inputSort.assign(begin,end);
968  typename std::vector<T>::iterator vit=inputSort.begin();
969  while(vit!=inputSort.end()){
970  if(*vit<minimum||*vit>maximum)
971  inputSort.erase(vit);
972  else
973  ++vit;
974  }
975  std::sort(inputSort.begin(),inputSort.end());
976  vit=inputSort.begin();
977  std::vector<T> inputBin;
978  for(int ibin=0;ibin<nbin;++ibin){
979  inputBin.clear();
980  while(inputBin.size()<inputSort.size()/nbin&&vit!=inputSort.end()){
981  inputBin.push_back(*vit);
982  ++vit;
983  }
984  if(inputBin.size()){
985  output[ibin]=mymax(inputBin);
986  }
987  }
988  if(!filename.empty()){
989  std::ofstream outputfile;
990  outputfile.open(filename.c_str());
991  if(!outputfile){
992  std::ostringstream s;
993  s<<"error opening distribution file , " << filename;
994  throw(s.str());
995  }
996  for(int ibin=0;ibin<nbin;++ibin)
997  outputfile << ibin*100.0/nbin << " " << static_cast<double>(output[ibin])/input.size() << std::endl;
998  outputfile.close();
999  }
1000 }
1001 
1002 //todo: what with nodata values?
1003 template<class T> T StatFactory::percentile(const std::vector<T>& input, typename std::vector<T>::const_iterator begin, typename std::vector<T>::const_iterator end, double percent, T minimum, T maximum) const
1004 {
1005  if(input.empty()){
1006  std::ostringstream s;
1007  s<<"Error: input is empty";
1008  throw(s.str());
1009  }
1010  std::vector<T> inputSort;
1011  inputSort.assign(begin,end);
1012  typename std::vector<T>::iterator vit=inputSort.begin();
1013  while(vit!=inputSort.end()){
1014  if(maximum>minimum){
1015  if(*vit<minimum||*vit>maximum)
1016  inputSort.erase(vit);
1017  }
1018  else
1019  ++vit;
1020  }
1021  std::sort(inputSort.begin(),inputSort.end());
1022  return gsl_stats_quantile_from_sorted_data(&(inputSort[0]),1,inputSort.size(),percent/100.0);
1023 }
1024 
1025 template<class T> void StatFactory::signature(const std::vector<T>& input, double&k, double& alpha, double& beta, double e) const
1026 {
1027  double m1=moment(input,1);
1028  double m2=moment(input,2);
1029  signature(m1,m2,k,alpha,beta,e);
1030 }
1031 
1032 //todo: what with nodata values?
1033 template<class T> void StatFactory::normalize(const std::vector<T>& input, std::vector<double>& output) const{
1034  double total=sum(input);
1035  if(total){
1036  output.resize(input.size());
1037  for(int index=0;index<input.size();++index)
1038  output[index]=input[index]/total;
1039  }
1040  else
1041  output=input;
1042 }
1043 
1044 //todo: what with nodata values?
1045 template<class T> void StatFactory::normalize_pct(std::vector<T>& input) const{
1046  double total=sum(input);
1047  if(total){
1048  typename std::vector<T>::iterator it;
1049  for(it=input.begin();it!=input.end();++it)
1050  *it=100.0*(*it)/total;
1051  }
1052 }
1053 
1054 template<class T> double StatFactory::rmse(const std::vector<T>& x, const std::vector<T>& y) const{
1055  if(x.size()!=y.size()){
1056  std::ostringstream s;
1057  s<<"Error: x and y not equal in size";
1058  throw(s.str());
1059  }
1060  if(x.empty()){
1061  std::ostringstream s;
1062  s<<"Error: x is empty";
1063  throw(s.str());
1064  }
1065  double mse=0;
1066  for(int isample=0;isample<x.size();++isample){
1067  if(isNoData(x[isample])||isNoData(y[isample]))
1068  continue;
1069  double e=x[isample]-y[isample];
1070  mse+=e*e/x.size();
1071  }
1072  return sqrt(mse);
1073 }
1074 
1075 // template<class T> double StatFactory::gsl_correlation(const std::vector<T>& x, const std::vector<T>& y) const{
1076 // return(gsl_stats_correlation(&(x[0]),1,&(y[0]),1,x.size()));
1077 // }
1078 
1079  template<class T> double StatFactory::gsl_covariance(const std::vector<T>& x, const std::vector<T>& y) const{
1080  return(gsl_stats_covariance(&(x[0]),1,&(y[0]),1,x.size()));
1081  }
1082 
1083 
1084 template<class T> double StatFactory::correlation(const std::vector<T>& x, const std::vector<T>& y, int delay) const{
1085  double meanX=0;
1086  double meanY=0;
1087  double varX=0;
1088  double varY=0;
1089  double sXY=0;
1090  meanVar(x,meanX,varX);
1091  meanVar(y,meanY,varY);
1092  double denom = sqrt(varX*varY);
1093  bool isValid=false;
1094  if(denom){
1095  //Calculate the correlation series
1096  sXY = 0;
1097  for (int i=0;i<x.size();++i) {
1098  int j = i + delay;
1099  if (j < 0 || j >= y.size())
1100  continue;
1101  else if(isNoData(x[i])||isNoData(y[j]))
1102  continue;
1103  else{
1104  isValid=true;
1105  if(i<0){
1106  std::ostringstream s;
1107  s<<"Error: i must be positive";
1108  throw(s.str());
1109  }
1110  if(i>=x.size()){
1111  std::ostringstream s;
1112  s<<"Error: i must be smaller than x.size()";
1113  throw(s.str());
1114  }
1115  if(j<0){
1116  std::ostringstream s;
1117  s<<"Error: j must be positive";
1118  throw(s.str());
1119  }
1120  if(j>=y.size()){
1121  std::ostringstream s;
1122  s<<"Error: j must be smaller than y.size()";
1123  throw(s.str());
1124  }
1125  sXY += (x[i] - meanX) * (y[j] - meanY);
1126  }
1127  }
1128  if(isValid){
1129  double minSize=(x.size()<y.size())?x.size():y.size();
1130  return(sXY / denom / (minSize-1));
1131  }
1132  else if(m_noDataValues.size())
1133  return m_noDataValues[0];
1134  else{
1135  std::string errorString="Error: no valid data found";
1136  throw(errorString);
1137  }
1138  }
1139  else
1140  return 0;
1141 }
1142 
1143 //todo: what if no valid data?
1144 template<class T> double StatFactory::cross_correlation(const std::vector<T>& x, const std::vector<T>& y, int maxdelay, std::vector<T>& z) const{
1145  z.clear();
1146  double sumCorrelation=0;
1147  for (int delay=-maxdelay;delay<maxdelay;delay++) {
1148  z.push_back(correlation(x,y,delay));
1149  sumCorrelation+=z.back();
1150  }
1151  return sumCorrelation;
1152 }
1153 
1154 //todo: nodata?
1155 template<class T> double StatFactory::linear_regression(const std::vector<T>& x, const std::vector<T>& y, double &c0, double &c1) const{
1156  if(x.size()!=y.size()){
1157  std::ostringstream s;
1158  s<<"Error: x and y not equal in size";
1159  throw(s.str());
1160  }
1161  if(x.empty()){
1162  std::ostringstream s;
1163  s<<"Error: x is empty";
1164  throw(s.str());
1165  }
1166  double cov00;
1167  double cov01;
1168  double cov11;
1169  double sumsq;
1170  gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1171  return (1-sumsq/var(y)/(y.size()-1));
1172 }
1173 
1174 //todo: nodata?
1175 template<class T> double StatFactory::linear_regression_err(const std::vector<T>& x, const std::vector<T>& y, double &c0, double &c1) const{
1176  if(x.size()!=y.size()){
1177  std::ostringstream s;
1178  s<<"Error: x and y not equal in size";
1179  throw(s.str());
1180  }
1181  if(x.empty()){
1182  std::ostringstream s;
1183  s<<"Error: x is empty";
1184  throw(s.str());
1185  }
1186  double cov00;
1187  double cov01;
1188  double cov11;
1189  double sumsq;
1190  gsl_fit_linear(&(x[0]),1,&(y[0]),1,x.size(),&c0,&c1,&cov00,&cov01,&cov11,&sumsq);
1191  return sqrt((sumsq)/(y.size()));
1192 }
1193 
1194 //alternatively: use GNU scientific library:
1195 // gsl_stats_correlation (const double data1[], const size_t stride1, const double data2[], const size_t stride2, const size_t n)
1196 
1197 template<class T> void StatFactory::interpolateNoData(const std::vector<double>& wavelengthIn, const std::vector<T>& input, const std::string& type, std::vector<T>& output, bool verbose) const{
1198  if(wavelengthIn.empty()){
1199  std::ostringstream s;
1200  s<<"Error: wavelengthIn is empty";
1201  throw(s.str());
1202  }
1203  std::vector<double> wavelengthOut=wavelengthIn;
1204  std::vector<T> validIn=input;
1205  if(input.size()!=wavelengthIn.size()){
1206  std::ostringstream s;
1207  s<<"Error: x and y not equal in size";
1208  throw(s.str());
1209  }
1210  int nband=wavelengthIn.size();
1211  output.clear();
1212  //remove nodata from input and corresponding wavelengthIn
1213  if(m_noDataValues.size()){
1214  typename std::vector<T>::iterator itValue=validIn.begin();
1215  typename std::vector<T>::iterator itWavelength=wavelengthOut.begin();
1216  while(itValue!=validIn.end()&&itWavelength!=wavelengthOut.end()){
1217  if(isNoData(*itValue)){
1218  validIn.erase(itValue);
1219  wavelengthOut.erase(itWavelength);
1220  }
1221  else{
1222  ++itValue;
1223  ++itWavelength;
1224  }
1225  }
1226  if(validIn.size()>1){
1227  try{
1228  interpolateUp(wavelengthOut, validIn, wavelengthIn, type, output, verbose);
1229  }
1230  catch(...){
1231  output=input;
1232  }
1233  }
1234  else//we can not interpolate if no valid data
1235  output=input;
1236  }
1237  else//no nodata values to interpolate
1238  output=input;
1239 }
1240 
1241 template<class T> void StatFactory::interpolateUp(const std::vector<double>& wavelengthIn, const std::vector<T>& input, const std::vector<double>& wavelengthOut, const std::string& type, std::vector<T>& output, bool verbose) const{
1242  if(wavelengthIn.empty()){
1243  std::ostringstream s;
1244  s<<"Error: wavelengthIn is empty";
1245  throw(s.str());
1246  }
1247  if(input.size()!=wavelengthIn.size()){
1248  std::ostringstream s;
1249  s<<"Error: input and wavelengthIn not equal in size";
1250  throw(s.str());
1251  }
1252  if(wavelengthOut.empty()){
1253  std::ostringstream s;
1254  s<<"Error: wavelengthOut is empty";
1255  throw(s.str());
1256  }
1257  int nband=wavelengthIn.size();
1258  output.clear();
1259  gsl_interp_accel *acc;
1260  allocAcc(acc);
1261  gsl_spline *spline;
1262  getSpline(type,nband,spline);
1263  assert(spline);
1264  assert(&(wavelengthIn[0]));
1265  assert(&(input[0]));
1266  int status=initSpline(spline,&(wavelengthIn[0]),&(input[0]),nband);
1267  if(status){
1268  std::string errorString="Could not initialize spline";
1269  throw(errorString);
1270  }
1271  for(int index=0;index<wavelengthOut.size();++index){
1272  if(wavelengthOut[index]<*wavelengthIn.begin()){
1273  output.push_back(*(input.begin()));
1274  continue;
1275  }
1276  else if(wavelengthOut[index]>wavelengthIn.back()){
1277  output.push_back(input.back());
1278  continue;
1279  }
1280  double dout=evalSpline(spline,wavelengthOut[index],acc);
1281  output.push_back(dout);
1282  }
1283  gsl_spline_free(spline);
1284  gsl_interp_accel_free(acc);
1285 }
1286 
1287 // 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){
1288 // assert(wavelengthIn.size());
1289 // assert(wavelengthOut.size());
1290 // int nsample=input.size();
1291 // int nband=wavelengthIn.size();
1292 // output.clear();
1293 // output.resize(nsample);
1294 // gsl_interp_accel *acc;
1295 // allocAcc(acc);
1296 // gsl_spline *spline;
1297 // getSpline(type,nband,spline);
1298 // for(int isample=0;isample<nsample;++isample){
1299 // assert(input[isample].size()==wavelengthIn.size());
1300 // initSpline(spline,&(wavelengthIn[0]),&(input[isample][0]),nband);
1301 // for(int index=0;index<wavelengthOut.size();++index){
1302 // if(type=="linear"){
1303 // if(wavelengthOut[index]<wavelengthIn.back())
1304 // output[isample].push_back(*(input.begin()));
1305 // else if(wavelengthOut[index]>wavelengthIn.back())
1306 // output[isample].push_back(input.back());
1307 // }
1308 // else{
1309 // double dout=evalSpline(spline,wavelengthOut[index],acc);
1310 // output[isample].push_back(dout);
1311 // }
1312 // }
1313 // }
1314 // gsl_spline_free(spline);
1315 // gsl_interp_accel_free(acc);
1316 // }
1317 
1318 //todo: nodata?
1319 template<class T> void StatFactory::interpolateUp(const std::vector<T>& input, std::vector<T>& output, int nbin) const
1320 {
1321  if(input.empty()){
1322  std::ostringstream s;
1323  s<<"Error: input is empty";
1324  throw(s.str());
1325  }
1326  if(!nbin){
1327  std::ostringstream s;
1328  s<<"Error: nbin must be larger than 0";
1329  throw(s.str());
1330  }
1331  output.clear();
1332  int dim=input.size();
1333  for(int i=0;i<dim;++i){
1334  double deltaX=0;
1335  double left=input[i];
1336  if(i<dim-1){
1337  double right=(i<dim-1)? input[i+1]:input[i];
1338  deltaX=(right-left)/static_cast<double>(nbin);
1339  for(int x=0;x<nbin;++x){
1340  output.push_back(left+x*deltaX);
1341  }
1342  }
1343  else
1344  output.push_back(input.back());
1345  }
1346 }
1347 
1348 //todo: nodata?
1349 template<class T> void StatFactory::nearUp(const std::vector<T>& input, std::vector<T>& output) const
1350 {
1351  if(input.empty()){
1352  std::ostringstream s;
1353  s<<"Error: input is empty";
1354  throw(s.str());
1355  }
1356  if(output.size()<input.size()){
1357  std::ostringstream s;
1358  s<<"Error: output size is smaller than input size";
1359  throw(s.str());
1360  }
1361  int dimInput=input.size();
1362  int dimOutput=output.size();
1363 
1364  for(int iin=0;iin<dimInput;++iin){
1365  for(int iout=0;iout<dimOutput/dimInput;++iout){
1366  int indexOutput=iin*dimOutput/dimInput+iout;
1367  if(indexOutput>=output.size()){
1368  std::ostringstream s;
1369  s<<"Error: indexOutput must be smaller than output.size()";
1370  throw(s.str());
1371  }
1372  output[indexOutput]=input[iin];
1373  }
1374  }
1375 }
1376 
1377 //todo: nodata?
1378 template<class T> void StatFactory::interpolateUp(double* input, int dim, std::vector<T>& output, int nbin)
1379 {
1380  if(!nbin){
1381  std::ostringstream s;
1382  s<<"Error: nbin must be larger than 0";
1383  throw(s.str());
1384  }
1385  output.clear();
1386  for(int i=0;i<dim;++i){
1387  double deltaX=0;
1388  double left=input[i];
1389  if(i<dim-1){
1390  double right=(i<dim-1)? input[i+1]:input[i];
1391  deltaX=(right-left)/static_cast<double>(nbin);
1392  for(int x=0;x<nbin;++x){
1393  output.push_back(left+x*deltaX);
1394  }
1395  }
1396  else
1397  output.push_back(input[dim-1]);
1398  }
1399 }
1400 
1401 //todo: nodata?
1402 template<class T> void StatFactory::interpolateDown(const std::vector<T>& input, std::vector<T>& output, int nbin) const
1403 {
1404  if(input.empty()){
1405  std::ostringstream s;
1406  s<<"Error: input is empty";
1407  throw(s.str());
1408  }
1409  if(!nbin){
1410  std::ostringstream s;
1411  s<<"Error: nbin must be larger than 0";
1412  throw(s.str());
1413  }
1414  output.clear();
1415  int dim=input.size();
1416  int x=0;
1417  output.push_back(input[0]);
1418  for(int i=1;i<dim;++i){
1419  if(i%nbin)
1420  continue;
1421  else{
1422  x=(i-1)/nbin+1;
1423  output.push_back(input[i]);
1424  }
1425  }
1426 }
1427 
1428 //todo: nodata?
1429 template<class T> void StatFactory::interpolateDown(double* input, int dim, std::vector<T>& output, int nbin)
1430 {
1431  if(!nbin){
1432  std::ostringstream s;
1433  s<<"Error: nbin must be larger than 0";
1434  throw(s.str());
1435  }
1436  output.clear();
1437  int x=0;
1438  output.push_back(input[0]);
1439  for(int i=1;i<dim;++i){
1440  if(i%nbin)
1441  continue;
1442  else{
1443  x=(i-1)/nbin+1;
1444  output.push_back(input[i]);
1445  }
1446  }
1447 }
1448 }
1449 
1450 #endif /* _STATFACTORY_H_ */
1451 
1452 // void Histogram::signature(double m1, double m2, double& k, double& alpha, double& beta, double e)
1453 // {
1454 // double y=m1*m1/m2;
1455 // beta=F_1(y,0.1,10.0,e);
1456 // double fb=F(beta);
1457 // double g=exp(lgamma(1.0/beta));
1458 // alpha=m1*g/exp(lgamma(2.0/beta));
1459 // k=beta/(2*alpha*g);
1460 // // std::cout << "y, alpha, beta: " << y << ", " << alpha << ", " << beta << std::endl;
1461 // }
1462 
1463 // double Histogram::F(double x)
1464 // {
1465 // double g2=exp(lgamma(2.0/x));
1466 // return(g2*g2/exp(lgamma(3.0/x))/exp(lgamma(1.0/x)));
1467 // }
1468 
1469 // //x1 is under estimate, x2 is over estimate, e is error
1470 // double Histogram::F_1(double y, double x1, double x2, double e)
1471 // {
1472 // double f1=F(x1);
1473 // double f2=F(x2);
1474 // assert(f1!=f2);
1475 // double x=x1+(x2-x1)*(y-f1)/(f2-f1);
1476 // double f=F(x);
1477 // while(f-y>=e||y-f>=e){
1478 // if(f<y)
1479 // x1=x;
1480 // else
1481 // x2=x;
1482 // if(x1==x2)
1483 // return x1;
1484 // assert(f1!=f2);
1485 // x=x1+(x2-x1)*(y-f1)/(f2-f1);
1486 // f=F(x);
1487 // }
1488 // return x;
1489 // }