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