00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, 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_TWOD_INTP_H 00024 #define O2SCL_TWOD_INTP_H 00025 00026 #include <iostream> 00027 #include <string> 00028 #include <vector> 00029 #include <o2scl/smart_interp.h> 00030 #include <o2scl/omatrix_tlate.h> 00031 00032 #ifndef DOXYGENP 00033 namespace o2scl { 00034 #endif 00035 00036 /** 00037 \brief Two-dimensional interpolation class 00038 00039 This class implements two-dimensional interpolation. Derivatives 00040 and integrals along both x- and y-directions can be computed. 00041 00042 The storage of the matrix to be specified in the function 00043 set_data() and this function is designed to follow the format: 00044 \f[ 00045 \begin{array}{cccc} 00046 & x_0 & x_1 & x_2 \\ 00047 y_0 & M_{00} & M_{01} & M_{02} \\ 00048 y_1 & M_{10} & M_{11} & M_{12} \\ 00049 y_2 & M_{20} & M_{21} & M_{22} 00050 \end{array} 00051 \f] 00052 thus the matrix should be \c M[i][j] where \c i is the y index 00053 and \c j is the row index. (See also the discussion in the 00054 User's guide in the section called \ref rowcol_subsect.) 00055 00056 The function set_data(), does not copy the data, it 00057 stores pointers to the data. If the data is modified, then the 00058 function reset_interp() must be called to reset the 00059 interpolation information with the original pointer information. 00060 The storage for the data, including the arrays \c x_grid and \c 00061 y_grid are all managed by the user. 00062 00063 By default, cubic spline interpolation with natural boundary 00064 conditions is used. This can be changed with the set_interp() 00065 function. 00066 00067 There is an example for the usage of this class given 00068 in <tt>examples/ex_twod_intp.cpp</tt>. 00069 00070 \future Could also include mixed second/first derivatives: 00071 deriv_xxy and deriv_xyy. 00072 \future Implement an improved caching system in case, for example 00073 \c xfirst is true and the last interpolation used the same 00074 value of \c x. 00075 00076 \comment 00077 00078 There's a whole list of things to think about here: 00079 00080 1. One could 00081 generalize the matrix and vector types. This would demand 00082 adding subvector types for the interpolation, etc. 00083 00084 2. One could copy the user data instead of storing pointers to it. 00085 00086 3. One could include get() and set() methods, but one would 00087 have to be careful about when to reset the interpolation. 00088 00089 4. xy order vs. yx order for the matrix 00090 00091 5/16/08 - There has been some talk about GSL implementing 2d 00092 interpolation, and I think the table3d class will be much more 00093 useful in the end. Thus, I think it's probably best to delay 00094 working on this for awhile until GSL 2d interpolation or table3d 00095 are more settled. 00096 00097 \endcomment 00098 00099 */ 00100 class twod_intp { 00101 public: 00102 00103 twod_intp(); 00104 00105 virtual ~twod_intp(); 00106 00107 /** 00108 \brief Initialize the data for the 2-dimensional interpolation 00109 00110 The interpolation type (passed directly to int_type) is 00111 specified in \c int_type and the data is specified in 00112 \c data. The data should be arranged so that the first 00113 array index is the y-value (the "row") and the second array 00114 index is the x-value (the "column"). The arrays \c xfun 00115 and \c yfun specify the two independent 00116 variables. \c xfun should be an array of length 00117 \c nx , and should be an array of length \c ny. The 00118 array \c data should be a two-dimensional array of size 00119 [ny][nx]. 00120 00121 If \c x_first is true, then set_data() creates interpolation 00122 objects for each of the rows. Calls to interp() then uses 00123 these to create a column at the specified value of \c x. An 00124 interpolation object is created at this column to find the 00125 value of the function at the specified value \c y. If \c 00126 x_first is false, the opposite strategy is employed. These two 00127 options may give slightly different results. 00128 00129 */ 00130 int set_data(size_t n_x, size_t n_y, ovector &x_grid, 00131 ovector &y_grid, omatrix &data, 00132 bool x_first=true); 00133 00134 /** 00135 \brief Reset the stored interpolation since the data has changed 00136 00137 This will return an error if the set_data() has not been called 00138 */ 00139 int reset_interp(); 00140 00141 /** \brief Perform the 2-d interpolation 00142 */ 00143 double interp(double x, double y); 00144 00145 /** \brief Compute the partial derivative in the x-direction 00146 */ 00147 double deriv_x(double x, double y); 00148 00149 /** \brief Compute the partial second derivative in the x-direction 00150 */ 00151 double deriv2_x(double x, double y); 00152 00153 /** \brief Compute the integral in the x-direction between x=x0 00154 and x=x1 00155 */ 00156 double integ_x(double x0, double x1, double y); 00157 00158 /** \brief Compute the partial derivative in the y-direction 00159 */ 00160 double deriv_y(double x, double y); 00161 00162 /** \brief Compute the partial second derivative in the y-direction 00163 */ 00164 double deriv2_y(double x, double y); 00165 00166 /** \brief Compute the integral in the y-direction between y=y0 00167 and y=y1 00168 */ 00169 double integ_y(double x, double y0, double y1); 00170 00171 /** \brief Compute the mixed partial derivative 00172 \f$ \frac{\partial^2 f}{\partial x \partial y} \f$ 00173 */ 00174 double deriv_xy(double x, double y); 00175 00176 /** 00177 \brief Specify the base interpolation objects to use 00178 00179 This allows the user to provide new interpolation objects for 00180 use in the two-dimensional interpolation. For 00181 example, 00182 00183 \code 00184 twod_intp ti; 00185 ovector x(20), y(40); 00186 omatrix d(40,20); 00187 00188 // Fill x, y, and d with the data, choose linear interpolation 00189 // instead of the default cubic spline 00190 def_interp_mgr<ovector_const_view,linear_interp> dim1; 00191 def_interp_mgr<ovector_const_subvector,linear_interp> dim2; 00192 00193 ti.set_interp(dim1,dim2); 00194 ti.set_data(20,40,x,y,d,true); 00195 \endcode 00196 00197 This function automatically calls reset_interp() and 00198 reset_data() if the data has already been set to reset the 00199 internal interpolation objects. 00200 */ 00201 int set_interp(base_interp_mgr<ovector_const_view> &b1, 00202 base_interp_mgr<ovector_const_subvector> &b2) { 00203 bimp1=&b1; 00204 bimp2=&b2; 00205 return 0; 00206 } 00207 00208 /** 00209 \brief Inform the class the data has been modified or 00210 changed in a way that set_data() will need to be called 00211 again. 00212 */ 00213 int unset_data(); 00214 00215 #ifndef DOXYGENP 00216 00217 protected: 00218 00219 /// \name Default interpolation manager 00220 //@{ 00221 def_interp_mgr<ovector_const_view,cspline_interp> dim1; 00222 def_interp_mgr<ovector_const_subvector,cspline_interp> dim2; 00223 //@} 00224 00225 /// \name Interpolation manager pointers 00226 //@{ 00227 base_interp_mgr<ovector_const_view> *bimp1; 00228 base_interp_mgr<ovector_const_subvector> *bimp2; 00229 //@} 00230 00231 /// The array of sm_interp pointers 00232 sm_interp_vec **si; 00233 /// An array of omatrix_rows 00234 omatrix_row **ari; 00235 /// An array of omatrix_cols 00236 omatrix_col **aci; 00237 00238 /// The number of x grid points 00239 size_t nx; 00240 00241 /// The number of y grid points 00242 size_t ny; 00243 00244 /// If true, the user has specified the base interpolation objects 00245 bool it_set; 00246 00247 /// True if the x interpolation should be done first 00248 bool xfirst; 00249 00250 /// True if the data has been specified by the user 00251 bool data_set; 00252 00253 /// The x grid 00254 ovector *xfun; 00255 00256 /// The y grid 00257 ovector *yfun; 00258 00259 /// The data 00260 omatrix *datap; 00261 00262 /// Allocate the base interpolation objects 00263 int alloc_interp(size_t n); 00264 00265 /// Free the base interpolation objects 00266 int free_interp(); 00267 00268 /// Free the data and interpolation objects 00269 int free_data(); 00270 00271 #endif 00272 00273 }; 00274 00275 #ifndef DOXYGENP 00276 } 00277 #endif 00278 00279 #endif 00280 00281 00282
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