00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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$, values 00041 for \f$f_j\f$ are predicted given a value of x and y. (In 00042 contrast to 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. 00045 00046 This procedure will fail if the three closest points are 00047 co-linear, and interp() will then call set_err() and 00048 return zero. 00049 00050 There is no caching so the numeric values of the data may be 00051 freely changed between calls to interp(). 00052 00053 The vector and matrix types can be any types which have 00054 suitably defined functions \c operator[]. 00055 00056 \future Rewrite so that it never fails unless all the points in 00057 the data set lie on a line. This would probably demand sorting 00058 all of the points by distance from desired location. 00059 */ 00060 template<class vec_t, class mat_t> class planar_intp { 00061 public: 00062 00063 planar_intp() { 00064 data_set=false; 00065 } 00066 00067 /** \brief Initialize the data for the planar interpolation 00068 */ 00069 int set_data(size_t n_points, vec_t &x, vec_t &y, size_t n_dat, 00070 mat_t &dat) { 00071 if (n_points<3) { 00072 set_err_ret("Must provide at least three points in set_data()", 00073 gsl_efailed); 00074 } 00075 if (x.size()<n_points) { 00076 set_err_ret("x too small.",gsl_einval); 00077 } else if (y.size()<n_points) { 00078 set_err_ret("y too small.",gsl_einval); 00079 } 00080 np=n_points; 00081 nd=n_dat; 00082 ux=&x; 00083 uy=&y; 00084 udat=&dat; 00085 data_set=true; 00086 return 0; 00087 } 00088 00089 /** 00090 \brief Perform the planar interpolation 00091 00092 It is assumed that \c ip is properly allocated beforehand. 00093 */ 00094 int interp(double x, double y, vec_t &ip) { 00095 double x1, x2, x3, y1, y2, y3; 00096 return interp(x,y,ip,x1,y1,x2,y2,x3,y3); 00097 } 00098 00099 /** 00100 \brief Planar interpolation returning the closest points 00101 00102 This function interpolates \c x and \c y into the data 00103 returning \c ip. It also returns the three closest x- and 00104 y-values used for computing the plane. 00105 00106 It is assumed that \c ip is properly allocated beforehand. 00107 */ 00108 int interp(double x, double y, vec_t &ip, double &x1, double &y1, 00109 double &x2, double &y2, double &x3, double &y3) { 00110 size_t i, j; 00111 double c1, c2, c3, c4; 00112 int i1, i2, i3, i4; 00113 00114 if (data_set==false) { 00115 set_err_ret("Data not set in planar_intp::interp().",gsl_einval); 00116 } 00117 00118 /// Put in initial points 00119 i1=0; i2=1; i3=2; 00120 c1=sqrt(pow(x-(*ux)[0],2.0)+pow(y-(*uy)[0],2.0)); 00121 c2=sqrt(pow(x-(*ux)[1],2.0)+pow(y-(*uy)[1],2.0)); 00122 c3=sqrt(pow(x-(*ux)[2],2.0)+pow(y-(*uy)[2],2.0)); 00123 00124 /// Sort initial points 00125 if (c2<c1) { 00126 if (c3<c2) { 00127 // 321 00128 swap(i1,c1,i3,c3); 00129 } else if (c3<c1) { 00130 // 231 00131 swap(i1,c1,i2,c2); 00132 swap(i2,c2,i3,c3); 00133 } else { 00134 // 213 00135 swap(i1,c1,i2,c2); 00136 } 00137 } else { 00138 if (c3<c1) { 00139 // 312 00140 swap(i1,c1,i3,c3); 00141 swap(i2,c2,i3,c3); 00142 } else if (c3<c2) { 00143 // 132 00144 swap(i3,c3,i2,c2); 00145 } 00146 // 123 00147 } 00148 00149 // Go through remaining points and sort accordingly 00150 for(j=3;j<np;j++) { 00151 i4=j; 00152 c4=sqrt(pow(x-(*ux)[i4],2.0)+pow(y-(*uy)[i4],2.0)); 00153 if (c4<c1) { 00154 swap(i4,c4,i3,c3); 00155 swap(i3,c3,i2,c2); 00156 swap(i2,c2,i1,c1); 00157 } else if (c4<c2) { 00158 swap(i4,c4,i3,c3); 00159 swap(i3,c3,i2,c2); 00160 } else if (c4<c3) { 00161 swap(i4,c4,i3,c3); 00162 } 00163 } 00164 00165 // Solve for denominator: 00166 double a, b, c, den, z1, z2, z3; 00167 x1=(*ux)[i1]; 00168 x2=(*ux)[i2]; 00169 x3=(*ux)[i3]; 00170 y1=(*uy)[i1]; 00171 y2=(*uy)[i2]; 00172 y3=(*uy)[i3]; 00173 den=(-x2*y1+x3*y1+x1*y2-x3*y2-x1*y3+x2*y3); 00174 00175 if (den==0.0) { 00176 00177 set_err_ret("Interpolation failed.",gsl_efailed); 00178 00179 } else { 00180 00181 // Calculate value of function 00182 for(i=0;i<nd;i++) { 00183 z1=(*(udat[i]))[i1]; 00184 z2=(*(udat[i]))[i2]; 00185 z3=(*(udat[i]))[i3]; 00186 a=-(-y2*z1+y3*z1+y1*z2-y3*z2-y1*z3+y2*z3)/den; 00187 b=-(x2*z1-x3*z1-x1*z2+x3*z2+x1*z3-x2*z3)/den; 00188 c=-(x3*y2*z1-x2*y3*z1-x3*y1*z2+x1*y3*z2+x2*y1*z3-x1*y2*z3)/den; 00189 00190 ip[i]=a*x+b*y+c; 00191 } 00192 } 00193 00194 return 0; 00195 } 00196 00197 #ifndef DOXYGEN_INTERNAL 00198 00199 protected: 00200 00201 /// The number of points 00202 size_t np; 00203 /// The number of functions 00204 size_t nd; 00205 /// The x-values 00206 vec_t *ux; 00207 /// The y-values 00208 vec_t *uy; 00209 /// The data 00210 mat_t *udat; 00211 /// True if the data has been specified 00212 bool data_set; 00213 00214 /// Swap points 1 and 2. 00215 int swap(int &i1, double &c1, int &i2, double &c2) { 00216 int t; 00217 double tc; 00218 00219 t=i1; tc=c1; 00220 i1=i2; c1=c2; 00221 i2=t; c2=tc; 00222 00223 return 0; 00224 } 00225 00226 #endif 00227 00228 }; 00229 00230 #ifndef DOXYGENP 00231 } 00232 #endif 00233 00234 #endif 00235 00236 00237
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