twod_intp.h

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_TWOD_INTP_H
00024 #define O2SCL_TWOD_INTP_H
00025 
00026 #include <iostream>
00027 #include <string>
00028 #include <o2scl/smart_interp.h>
00029 #include <o2scl/omatrix_tlate.h>
00030 
00031 #ifndef DOXYGENP
00032 namespace o2scl {
00033 #endif
00034 
00035   /** 
00036       \brief Two-dimensional interpolation class 
00037 
00038       This class implements two-dimensional interpolation. Derivatives
00039       and integrals along both x- and y-directions can be computed.
00040       The function set_data(), does not copy the data but rather
00041       stores pointers to the data. If the data is modified, then the
00042       function reset_interp() can be called to reset the interpolation
00043       information with the original pointer information.
00044 
00045       The storage for the data, including the arrays \c x_fun 
00046       and \c y_fun are all managed by the user. If the data is
00047       changed without calling reset_interp(), then interp() will
00048       return incorrect results.
00049 
00050       By default, cubic spline interpolation with natural boundary
00051       conditions is used. This can be changed with the set_interp()
00052       function.
00053       
00054       \future Could also include mixed second/first derivatives:
00055       deriv_xxy and deriv_xyy.
00056 
00057   */
00058   class twod_intp {
00059   public:
00060     
00061     twod_intp();
00062 
00063     virtual ~twod_intp();
00064     
00065     /** \brief Initialize the data for the 2-dimensional interpolation
00066         
00067         The interpolation type (passed directly to int_type) is
00068         specified in \c int_type and the data is specified in
00069         \c data. The data should be arranged so that the first
00070         array index is the y-value (the "row") and the second array
00071         index is the x-value (the "column"). The arrays \c xfun 
00072         and \c yfun specify the two independent
00073         variables. \c xfun should be an array of length
00074         \c nx , and should be an array of length \c ny. The
00075         array \c data should be a two-dimensional array of size
00076         [ny][nx].
00077 
00078         If \c x_first is true, then set_data() creates interpolation
00079         objects for each of the rows. Calls to interp() then uses
00080         these to create a column at the specified value of \c x. An
00081         interpolation object is created at this column to find the
00082         value of the function at the specified value \c y. If \c
00083         x_first is false, the opposite strategy is employed. These two
00084         options may give slightly different results. In general, if
00085         the data is "more accurate" in the x direction than in the y
00086         direction, it is probably better to choose \c x_first=true.
00087         
00088     */
00089     int set_data(int n_x, int n_y, ovector &x_fun,
00090                  ovector &y_fun, omatrix &u_data, 
00091                  bool x_first=true);
00092     
00093     /** 
00094         \brief Reset the stored interpolation since the data has changed
00095         
00096         This will return an error if the set_data() has not been called
00097     */
00098     int reset_interp();
00099 
00100     /** \brief Perform the 2-d interpolation 
00101      */
00102     double interp(double x, double y);
00103 
00104     /** \brief Compute the partial derivative in the x-direction
00105      */
00106     double deriv_x(double x, double y);
00107 
00108     /** \brief Compute the partial second derivative in the x-direction
00109      */
00110     double deriv2_x(double x, double y);
00111 
00112     /** \brief Compute the integral in the x-direction between x=x0
00113         and x=x1
00114      */
00115     double integ_x(double x0, double x1, double y);
00116 
00117     /** \brief Compute the partial derivative in the y-direction
00118      */
00119     double deriv_y(double x, double y);
00120 
00121     /** \brief Compute the partial second derivative in the y-direction
00122      */
00123     double deriv2_y(double x, double y);
00124 
00125     /** \brief Compute the integral in the y-direction between y=y0
00126         and y=y1
00127      */
00128     double integ_y(double x, double y0, double y1);
00129 
00130     /** \brief Compute the mixed partial derivative 
00131         \f$ \frac{\partial^2 f}{\partial x \partial y} \f$
00132     */
00133     double deriv_xy(double x, double y);
00134 
00135     /** 
00136         \brief Specify the base interpolation objects to use
00137 
00138         \todo Document this function
00139     */
00140     int set_interp(size_t ni, base_interp<ovector_view> *it,
00141                    base_interp<ovector_const_subvector> *it_sub,
00142                    base_interp<ovector_view> &it2,
00143                    base_interp<ovector_const_subvector> &it2_sub) {
00144       nuit=ni;
00145       uit1=it;
00146       uit2=it_sub;
00147       uit3=&it2;
00148       uit4=&it2_sub;
00149       it_set=true;
00150       return 0;
00151     }
00152     
00153 #ifndef DOXYGENP
00154     
00155   protected:
00156 
00157     /// The array of sm_interp pointers 
00158     sm_interp_vec **si;
00159     /// An array of omatrix_rows
00160     omatrix_row **ari;
00161     /// An array of omatrix_cols
00162     omatrix_col **aci;
00163 
00164     /// The number of x grid points
00165     int nx;
00166 
00167     /// The number of y grid points
00168     int ny;
00169 
00170     /// If true, the user has specified the base interpolation objects
00171     bool it_set;
00172 
00173     /// Number of user-specified interpolation objects
00174     int nuit;
00175 
00176     /// \name User-specified interpolation objects
00177     //@{
00178     base_interp<ovector_view> *uit1;
00179     base_interp<ovector_const_subvector> *uit2;
00180     base_interp<ovector_view> *uit3;
00181     base_interp<ovector_const_subvector> *uit4;
00182     //@}
00183 
00184     /// True if the x interpolation should be done first
00185     bool xfirst;
00186 
00187     /// True if the data has been specified by the user
00188     bool data_set;
00189 
00190     /// The x grid
00191     ovector *xfun;
00192 
00193     /// The y grid
00194     ovector *yfun;
00195 
00196     /// The data
00197     omatrix *data;
00198 
00199 #endif
00200 
00201   };
00202   
00203 #ifndef DOXYGENP
00204 }
00205 #endif
00206 
00207 #endif
00208 
00209 
00210 

Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.

Project hosting provided by SourceForge.net Logo, O2scl Sourceforge Project Page