00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_ARRAY_H 00024 #define O2SCL_ARRAY_H 00025 00026 /** \file array.h 00027 \brief Various array classes 00028 00029 This file contains the three alloction classes 00030 - array_alloc 00031 - array_2d_alloc 00032 - pointer_alloc 00033 the classes for the manipulation of array in smart_interp 00034 - array_reverse 00035 - array_subvector 00036 - array_subvector_reverse 00037 - array_const_reverse 00038 - array_const_subvector 00039 - array_const_subvector_reverse 00040 the array equivalent of omatrix_row or umatrix_row (see 00041 usage proposed in src/ode/ode_it_solve_ts.cpp) 00042 - array_row 00043 and an array output function 00044 - array_out() 00045 00046 \todo Ensure that array_row works, either here or in 00047 src/ode/ode_it_solve_ts.cpp 00048 */ 00049 00050 #include <iostream> 00051 #include <cmath> 00052 #include <string> 00053 #include <fstream> 00054 #include <sstream> 00055 #include <o2scl/err_hnd.h> 00056 #include <gsl/gsl_ieee_utils.h> 00057 #include <gsl/gsl_sort.h> 00058 00059 #ifndef DOXYGENP 00060 namespace o2scl 00061 { 00062 #endif 00063 00064 /** 00065 \brief A simple class to provide an \c allocate() function 00066 for arrays 00067 00068 The functions here are blank, as fixed-length arrays are 00069 automatically allocated and destroyed by the compiler. This 00070 class is present to provide an analog to \ref pointer_alloc and 00071 \ref ovector_alloc 00072 */ 00073 template<class vec_t> class array_alloc { 00074 public: 00075 /// Allocate \c v for \c i elements 00076 void allocate(vec_t &v, size_t i) {} 00077 /// Free memory 00078 void free(vec_t &v) {} 00079 }; 00080 00081 /** 00082 \brief A simple class to provide an \c allocate() function 00083 for 2-dimensional arrays 00084 00085 The functions here are blank, as fixed-length arrays are 00086 automatically allocated and destroyed by the compiler. This 00087 class is present to provide an analog to \ref pointer_2d_alloc 00088 and \ref omatrix_alloc 00089 */ 00090 template<class mat_t> class array_2d_alloc { 00091 public: 00092 /// Allocate \c v for \c i elements 00093 void allocate(mat_t &v, size_t i, size_t j) {} 00094 /// Free memory 00095 void free(mat_t &v, size_t i) {} 00096 }; 00097 00098 /** 00099 \brief A simple class to provide an \c allocate() function 00100 for pointers 00101 00102 Uses \c new and \c delete. 00103 */ 00104 template<class base_t> class pointer_alloc { 00105 public: 00106 /// Allocate \c v for \c i elements 00107 void allocate(base_t *&v, size_t i) { v=new base_t[i]; } 00108 /// Free memory 00109 void free(base_t *&v) { delete[] v; } 00110 }; 00111 00112 /** 00113 \brief A simple class to provide an \c allocate() function 00114 for pointers 00115 00116 Uses \c new and \c delete. 00117 */ 00118 template<class base_t> class pointer_2d_alloc { 00119 public: 00120 /// Allocate \c v for \c i elements 00121 void allocate(base_t **&v, size_t i, size_t j) { 00122 v=new base_t *[i]; 00123 for(size_t m=0;m<i;m++) v[m]=new base_t[j]; 00124 } 00125 /// Free memory 00126 void free(base_t **&v, size_t i) { 00127 for(size_t m=0;m<i;m++) delete[] v[m]; 00128 delete[] v; 00129 } 00130 }; 00131 00132 /** 00133 \brief A simple class which reverses the order of an array 00134 */ 00135 template<size_t sz> class array_reverse { 00136 00137 #ifndef DOXYGEN_INTERNAL 00138 00139 protected: 00140 00141 /// The array pointer 00142 double *a; 00143 00144 #endif 00145 00146 public: 00147 00148 /// Create a reversed array from \c arr of size \c sz 00149 array_reverse(double *arr) { 00150 a=arr; 00151 } 00152 00153 /** 00154 \brief Array-like indexing 00155 */ 00156 double &operator[](size_t i) { 00157 return a[sz-1-i]; 00158 } 00159 00160 /** 00161 \brief Array-like indexing 00162 */ 00163 const double &operator[](size_t i) const { 00164 return a[sz-1-i]; 00165 } 00166 00167 }; 00168 00169 /** 00170 \brief A simple class which reverses the order of an array 00171 00172 00173 */ 00174 template<size_t sz> class array_const_reverse { 00175 00176 #ifndef DOXYGEN_INTERNAL 00177 00178 protected: 00179 00180 /// The array pointer 00181 double *a; 00182 00183 #endif 00184 00185 public: 00186 00187 /// Create a reversed array from \c arr of size \c sz 00188 array_const_reverse(const double *arr) { 00189 a=(double *)arr; 00190 } 00191 00192 /** 00193 \brief Array-like indexing 00194 */ 00195 const double &operator[](size_t i) const { 00196 return a[sz-1-i]; 00197 } 00198 00199 }; 00200 00201 /** 00202 \brief A simple subvector class for an array (without error checking) 00203 00204 00205 */ 00206 class array_subvector { 00207 00208 #ifndef DOXYGEN_INTERNAL 00209 00210 protected: 00211 00212 /// The array pointer 00213 double *a; 00214 00215 /// The offset 00216 size_t off; 00217 00218 /// The subvector length 00219 size_t len; 00220 00221 #endif 00222 00223 public: 00224 00225 /// Create a reversed array from \c arr of size \c sz 00226 array_subvector(double *arr, size_t offset, size_t n) { 00227 a=arr+offset; 00228 off=offset; 00229 len=n; 00230 } 00231 00232 /** 00233 \brief Array-like indexing 00234 */ 00235 double &operator[](size_t i) { 00236 return a[i]; 00237 } 00238 00239 /** 00240 \brief Array-like indexing 00241 */ 00242 const double &operator[](size_t i) const { 00243 return a[i]; 00244 } 00245 00246 }; 00247 00248 /** \brief A simple subvector class for a const array 00249 (without error checking) 00250 00251 00252 */ 00253 class array_const_subvector { 00254 00255 #ifndef DOXYGEN_INTERNAL 00256 00257 protected: 00258 00259 /// The array pointer 00260 double *a; 00261 00262 /// The offset 00263 size_t off; 00264 00265 /// The subvector length 00266 size_t len; 00267 00268 #endif 00269 00270 public: 00271 00272 /// Create a reversed array from \c arr of size \c sz 00273 array_const_subvector(const double *arr, size_t offset, size_t n) { 00274 a=(double *)arr+offset; 00275 off=offset; 00276 len=n; 00277 } 00278 00279 /** 00280 \brief Array-like indexing 00281 */ 00282 const double &operator[](size_t i) const { 00283 return a[i]; 00284 } 00285 00286 }; 00287 00288 /** 00289 \brief Reverse a subvector of an array 00290 00291 00292 */ 00293 class array_subvector_reverse { 00294 00295 #ifndef DOXYGEN_INTERNAL 00296 00297 protected: 00298 00299 /// The array pointer 00300 double *a; 00301 00302 /// The offset 00303 size_t off; 00304 00305 /// The subvector length 00306 size_t len; 00307 00308 #endif 00309 00310 public: 00311 00312 /// Create a reversed array from \c arr of size \c sz 00313 array_subvector_reverse(double *arr, size_t offset, size_t n) { 00314 a=arr+offset; 00315 off=offset; 00316 len=n; 00317 } 00318 00319 /** 00320 \brief Array-like indexing 00321 */ 00322 double &operator[](size_t i) { 00323 return a[len-1-i]; 00324 } 00325 00326 /** 00327 \brief Array-like indexing 00328 */ 00329 const double &operator[](size_t i) const { 00330 return a[len-1-i]; 00331 } 00332 00333 }; 00334 00335 /** 00336 \brief Reverse a subvector of a const array 00337 00338 00339 */ 00340 class array_const_subvector_reverse { 00341 00342 #ifndef DOXYGEN_INTERNAL 00343 00344 protected: 00345 00346 /// The array pointer 00347 double *a; 00348 00349 /// The offset 00350 size_t off; 00351 00352 /// The subvector length 00353 size_t len; 00354 00355 #endif 00356 00357 public: 00358 00359 /// Create a reversed array from \c arr of size \c sz 00360 array_const_subvector_reverse(const double *arr, size_t offset, size_t n) { 00361 a=(double *)arr+offset; 00362 off=offset; 00363 len=n; 00364 } 00365 00366 /** 00367 \brief Array-like indexing 00368 */ 00369 const double &operator[](size_t i) const { 00370 return a[len-1-i]; 00371 } 00372 00373 }; 00374 00375 /// Extract a row of a C-style 2d-array 00376 template <class data_t, class array_2d_t> class array_row { 00377 public: 00378 data_t *p; 00379 array_row(array_2d_t &a, size_t i) { 00380 p=&(a[i][0]); 00381 } 00382 data_t &operator[](size_t i) { 00383 return p[i]; 00384 } 00385 }; 00386 00387 /** 00388 \brief Output a vector to a stream 00389 00390 No trailing space is output after the last element, and an 00391 endline is output only if \c endline is set to \c true. If the 00392 parameter \c n is zero, this function silently does nothing. 00393 */ 00394 template<class vec_t> 00395 int vector_out(std::ostream &os, size_t n, vec_t &v, 00396 bool endline=false) { 00397 // This next line is important since n-1 is not well-defined if n=0 00398 if (n==0) return 0; 00399 for(size_t i=0;i<n-1;i++) os << v[i] << " "; 00400 os << v[n-1]; 00401 if (endline) os << std::endl; 00402 return 0; 00403 } 00404 00405 /// Provide a downheap() function for vector_sort() 00406 template<class data_t, class vec_t> 00407 void sort_downheap(vec_t &data, const size_t N, size_t k) { 00408 00409 data_t v=data[k]; 00410 00411 while (k<=N/2) { 00412 size_t j=2*k; 00413 00414 if (j<N && data[j] < data[(j+1)]) j++; 00415 if (!(v < data[j])) break; 00416 data[k]=data[j]; 00417 k=j; 00418 } 00419 00420 data[k]=v; 00421 } 00422 00423 /// Sort a vector 00424 template<class data_t, class vec_t> 00425 int vector_sort(const size_t n, vec_t &data) { 00426 00427 size_t N; 00428 size_t k; 00429 00430 if (n==0) return 0; 00431 00432 N=n-1; 00433 k=N/2; 00434 k++; 00435 do { 00436 k--; 00437 sort_downheap<data_t,vec_t>(data,N,k); 00438 } while (k > 0); 00439 00440 while (N > 0) { 00441 double tmp=data[0]; 00442 data[0]=data[N]; 00443 data[N]=tmp; 00444 N--; 00445 sort_downheap<data_t,vec_t>(data,N,0); 00446 } 00447 00448 return 0; 00449 } 00450 00451 /// "Rotate" a vector so that the kth element is now the beginning 00452 template<class data_t, class vec_t> 00453 int vector_rotate(const size_t n, vec_t &data, size_t k) { 00454 00455 double *tmp=new double[n]; 00456 for(size_t i=0;i<n;i++) { 00457 tmp[i]=data[(i+k)%n]; 00458 } 00459 for(size_t i=0;i<n;i++) { 00460 data[i]=tmp[i]; 00461 } 00462 delete[] tmp; 00463 00464 return 0; 00465 } 00466 00467 /// Compute the maximum of the first \c n elements of a vector 00468 template<class data_t, class vec_t> 00469 int vector_max(const size_t n, vec_t &data, data_t &max, size_t &ix) { 00470 00471 if (n==0) { 00472 set_err_ret("Send size=0 to vector_max().",gsl_efailed); 00473 } 00474 ix=0; 00475 max=data[0]; 00476 for(size_t i=1;i<n;i++) { 00477 if (data[i]>max) { 00478 max=data[i]; 00479 ix=i; 00480 } 00481 } 00482 return 0; 00483 } 00484 00485 /// Compute the minimum of the first \c n elements of a vector 00486 template<class data_t, class vec_t> 00487 int vector_min(const size_t n, vec_t &data, data_t &min, size_t &ix) { 00488 00489 if (n==0) { 00490 set_err_ret("Send size=0 to vector_min().",gsl_efailed); 00491 } 00492 ix=0; 00493 min=data[0]; 00494 for(size_t i=1;i<n;i++) { 00495 if (data[i]<min) { 00496 min=data[i]; 00497 ix=i; 00498 } 00499 } 00500 return 0; 00501 } 00502 00503 /** 00504 \brief Compute the sum of the first \c n elements of a vector 00505 00506 If \c n is zero, this will set \c avg to zero and return 00507 \ref gsl_success. 00508 */ 00509 template<class data_t, class vec_t> 00510 int vector_sum(const size_t n, vec_t &data, data_t &sum) { 00511 00512 sum=0; 00513 for(size_t i=0;i<n;i++) { 00514 sum+=data[i]; 00515 } 00516 return 0; 00517 } 00518 00519 /** 00520 \brief Compute the average of the first \c n elements of a vector 00521 00522 If \c n is zero, this will set \c avg to zero and return 00523 \ref gsl_success. 00524 */ 00525 template<class data_t, class vec_t> 00526 int vector_avg(const size_t n, vec_t &data, data_t &avg) { 00527 00528 avg=0; 00529 for(size_t i=0;i<n;i++) { 00530 avg+=data[i]; 00531 } 00532 avg/=((data_t)n); 00533 return 0; 00534 } 00535 00536 /** \brief Compute the variance of the first \c n elements of a vector 00537 given the mean \c mean. 00538 00539 If \c n is zero, this will set \c avg to zero and return 00540 \ref gsl_success. 00541 */ 00542 template<class data_t, class vec_t> 00543 int vector_variance(const size_t n, vec_t &data, data_t &mean, 00544 data_t &var) { 00545 00546 var=0; 00547 for(size_t i=0;i<n;i++) { 00548 data_t delta=(data[i]-mean); 00549 var+=(delta*delta-var)/((double)(i+1)); 00550 } 00551 return 0; 00552 } 00553 00554 /** \brief Reverse a vector 00555 */ 00556 template<class data_t, class vec_t> 00557 int vector_reverse(const size_t n, vec_t &data) { 00558 data_t tmp; 00559 00560 for(size_t i=0;i<n/2;i++) { 00561 tmp=data[n-1-i]; 00562 data[n-1-i]=data[i]; 00563 data[i]=tmp; 00564 } 00565 return 0; 00566 } 00567 00568 /** \brief Compute the standard deviation of the first \c n elements 00569 of a vector 00570 */ 00571 template<class data_t, class vec_t> 00572 int vector_stdev(const size_t n, vec_t &data, data_t &var) { 00573 00574 data_t mean; 00575 vector_mean(n,data,mean); 00576 var=0; 00577 for(size_t i=0;i<n;i++) { 00578 data_t delta=(data[i]-mean); 00579 var+=(delta*delta-var)/((double)(i+1)); 00580 } 00581 var=sqrt(var*((double)n)/((double)(n-1))); 00582 return 0; 00583 } 00584 00585 #ifndef DOXYGENP 00586 } 00587 #endif 00588 00589 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page