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