![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006-2012, 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_PLANAR_INTP_H 00024 #define O2SCL_PLANAR_INTP_H 00025 00026 #include <iostream> 00027 #include <string> 00028 #include <cmath> 00029 #include <o2scl/err_hnd.h> 00030 #include <o2scl/omatrix_tlate.h> 00031 00032 #ifndef DOXYGENP 00033 namespace o2scl { 00034 #endif 00035 00036 /** \brief Interpolate among two independent variables with planes 00037 00038 This class is experimental. 00039 00040 This is an analog of 1-d linear interpolation for two 00041 dimensions, which is particularly useful with the data points 00042 are not arranged in a specified order (i.e. on a grid). For a 00043 set of data \f${x_i,y_i,f_{j,i}}\f$, the values of \f$f_j\f$ are 00044 predicted given a value of x and y. In contrast to \ref 00045 twod_intp, the data need not be presented in a grid. This 00046 interpolation is performed by finding the plane that goes 00047 through three closest points in the data set. Distances are 00048 determined with 00049 \f[ 00050 d_{ij} = \sqrt{\left(\frac{x_i-x_j}{\Delta x}\right)^2 + 00051 \left(\frac{y_i-y_j}{\Delta y}\right)^2} 00052 \f] 00053 where \f$ \Delta x = x_{\mathrm{max}}-x_{\mathrm{min}} \f$ 00054 and \f$ \Delta y = y_{\mathrm{max}}-y_{\mathrm{min}} \f$ . 00055 00056 This procedure will fail if the three closest points are 00057 co-linear, and interp() will then call the error handler. 00058 00059 There is no caching so the numeric values of the data may be 00060 freely changed between calls to interp(). 00061 00062 The vector and matrix types can be any types which have 00063 suitably defined functions \c operator[]. 00064 00065 For example, with 10 random points in the x-y plane with \f$ 00066 -1<x<1 \f$ and \f$ -1<y<1 \f$, the planes are demarcated 00067 according to 00068 \image html in_planar_intp.png "in_planar_intp.png" 00069 \image latex in_planar_intp.pdf "in_planar_intp.pdf" width=9cm 00070 where each polygonal region represents the set of all points 00071 in the domain which will be mapped to the same plane to approximate 00072 the function. 00073 00074 \future Rewrite so that it never fails unless all the points in 00075 the data set lie on a line. This would probably demand sorting 00076 all of the points by distance from desired location. (11/1/11 - 00077 Maybe not. We might be able to get away with just storing an 00078 vector which records the index j of the point in the data set 00079 which is the nearest neighbor for each point i in the data set.) 00080 00081 \future Instead of comparing <tt>den</tt> with zero, add the 00082 option of comparing it with a small number. 00083 00084 \future Allow the user to specify the scales used for x and 00085 y, instead of computing max-min every time. 00086 00087 \future The generalization to more dimensions might be 00088 straight forward. 00089 */ 00090 template<class vec_t, class mat_t> class planar_intp { 00091 public: 00092 00093 planar_intp() { 00094 data_set=false; 00095 } 00096 00097 /** \brief Initialize the data for the planar interpolation 00098 */ 00099 int set_data(size_t n_points, vec_t &x, vec_t &y, size_t n_dat, 00100 mat_t &dat) { 00101 if (n_points<3) { 00102 O2SCL_ERR2_RET("Must provide at least three points in ", 00103 "planar_intp::set_data()",gsl_efailed); 00104 } 00105 np=n_points; 00106 nd=n_dat; 00107 ux=&x; 00108 uy=&y; 00109 udat=&dat; 00110 data_set=true; 00111 return 0; 00112 } 00113 00114 /** \brief Perform the planar interpolation 00115 00116 It is assumed that \c ip is properly allocated beforehand. 00117 */ 00118 int interp(double x, double y, vec_t &ip) { 00119 double x1, x2, x3, y1, y2, y3; 00120 size_t i1, i2, i3; 00121 return interp(x,y,ip,i1,x1,y1,i2,x2,y2,i3,x3,y3); 00122 } 00123 00124 /** \brief Planar interpolation returning the closest points 00125 00126 This function interpolates \c x and \c y into the data 00127 returning \c ip. It also returns the three closest x- and 00128 y-values used for computing the plane. 00129 00130 It is assumed that \c ip is properly allocated beforehand. 00131 */ 00132 int interp(double x, double y, vec_t &ip, 00133 size_t &i1, double &x1, double &y1, 00134 size_t &i2, double &x2, double &y2, 00135 size_t &i3, double &x3, double &y3) { 00136 00137 size_t i, j; 00138 double c1, c2, c3, c4; 00139 size_t i4; 00140 00141 if (data_set==false) { 00142 O2SCL_ERR_RET("Data not set in planar_intp::interp().",gsl_einval); 00143 } 00144 00145 // Find scaling 00146 double minx=(*ux)[0], maxx=(*ux)[0]; 00147 double miny=(*uy)[0], maxy=(*uy)[0]; 00148 for(i=1;i<np;i++) { 00149 if ((*ux)[i]<minx) minx=(*ux)[i]; 00150 if ((*ux)[i]>maxx) maxx=(*ux)[i]; 00151 if ((*uy)[i]<miny) miny=(*uy)[i]; 00152 if ((*uy)[i]>maxy) maxy=(*uy)[i]; 00153 } 00154 double dx=maxx-minx; 00155 double dy=maxy-miny; 00156 00157 /// Put in initial points 00158 i1=0; i2=1; i3=2; 00159 c1=sqrt(pow((x-(*ux)[0])/dx,2.0)+pow((y-(*uy)[0])/dy,2.0)); 00160 c2=sqrt(pow((x-(*ux)[1])/dx,2.0)+pow((y-(*uy)[1])/dy,2.0)); 00161 c3=sqrt(pow((x-(*ux)[2])/dx,2.0)+pow((y-(*uy)[2])/dy,2.0)); 00162 00163 /// Sort initial points 00164 if (c2<c1) { 00165 if (c3<c2) { 00166 // 321 00167 swap(i1,c1,i3,c3); 00168 } else if (c3<c1) { 00169 // 231 00170 swap(i1,c1,i2,c2); 00171 swap(i2,c2,i3,c3); 00172 } else { 00173 // 213 00174 swap(i1,c1,i2,c2); 00175 } 00176 } else { 00177 if (c3<c1) { 00178 // 312 00179 swap(i1,c1,i3,c3); 00180 swap(i2,c2,i3,c3); 00181 } else if (c3<c2) { 00182 // 132 00183 swap(i3,c3,i2,c2); 00184 } 00185 // 123 00186 } 00187 00188 // Go through remaining points and sort accordingly 00189 for(j=3;j<np;j++) { 00190 i4=j; 00191 c4=sqrt(pow((x-(*ux)[i4])/dx,2.0)+pow((y-(*uy)[i4])/dy,2.0)); 00192 if (c4<c1) { 00193 swap(i4,c4,i3,c3); 00194 swap(i3,c3,i2,c2); 00195 swap(i2,c2,i1,c1); 00196 } else if (c4<c2) { 00197 swap(i4,c4,i3,c3); 00198 swap(i3,c3,i2,c2); 00199 } else if (c4<c3) { 00200 swap(i4,c4,i3,c3); 00201 } 00202 } 00203 00204 // Solve for denominator: 00205 double a, b, c, den, z1, z2, z3; 00206 x1=(*ux)[i1]; 00207 x2=(*ux)[i2]; 00208 x3=(*ux)[i3]; 00209 y1=(*uy)[i1]; 00210 y2=(*uy)[i2]; 00211 y3=(*uy)[i3]; 00212 den=(-x2*y1+x3*y1+x1*y2-x3*y2-x1*y3+x2*y3); 00213 00214 if (den==0.0) { 00215 00216 O2SCL_ERR2_RET("Interpolation failed in.", 00217 "planar_intp::interp().",gsl_efailed); 00218 00219 } else { 00220 00221 // Calculate value of function 00222 for(i=0;i<nd;i++) { 00223 z1=(*(udat[i]))[i1]; 00224 z2=(*(udat[i]))[i2]; 00225 z3=(*(udat[i]))[i3]; 00226 a=-(-y2*z1+y3*z1+y1*z2-y3*z2-y1*z3+y2*z3)/den; 00227 b=-(x2*z1-x3*z1-x1*z2+x3*z2+x1*z3-x2*z3)/den; 00228 c=-(x3*y2*z1-x2*y3*z1-x3*y1*z2+x1*y3*z2+x2*y1*z3-x1*y2*z3)/den; 00229 00230 ip[i]=a*x+b*y+c; 00231 } 00232 } 00233 00234 return 0; 00235 } 00236 00237 #ifndef DOXYGEN_INTERNAL 00238 00239 protected: 00240 00241 /// The number of points 00242 size_t np; 00243 /// The number of functions 00244 size_t nd; 00245 /// The x-values 00246 vec_t *ux; 00247 /// The y-values 00248 vec_t *uy; 00249 /// The data 00250 mat_t *udat; 00251 /// True if the data has been specified 00252 bool data_set; 00253 00254 /// Swap points 1 and 2. 00255 int swap(size_t &i1, double &c1, size_t &i2, double &c2) { 00256 int t; 00257 double tc; 00258 00259 t=i1; tc=c1; 00260 i1=i2; c1=c2; 00261 i2=t; c2=tc; 00262 00263 return 0; 00264 } 00265 00266 #endif 00267 00268 }; 00269 00270 #ifndef DOXYGENP 00271 } 00272 #endif 00273 00274 #endif 00275 00276 00277
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).