Object-oriented Scientific Computing Library: Version 0.910
gsl_miser.h
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 /* monte/miser.c
00024  * 
00025  * Copyright (C) 1996, 1997, 1998, 1999, 2000 Michael Booth
00026  * 
00027  * This program is free software; you can redistribute it and/or modify
00028  * it under the terms of the GNU General Public License as published by
00029  * the Free Software Foundation; either version 3 of the License, or (at
00030  * your option) any later version.
00031  * 
00032  * This program is distributed in the hope that it will be useful, but
00033  * WITHOUT ANY WARRANTY; without even the implied warranty of
00034  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00035  * General Public License for more details.
00036  * 
00037  * You should have received a copy of the GNU General Public License
00038  * along with this program; if not, write to the Free Software
00039  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 
00040  * 02110-1301, USA.
00041  */
00042 #ifndef O2SCL_GSL_MISER_H
00043 #define O2SCL_GSL_MISER_H
00044 
00045 #include <iostream>
00046 #include <o2scl/misc.h>
00047 #include <o2scl/mcarlo_inte.h>
00048 #include <o2scl/uvector_tlate.h>
00049 #include <o2scl/gsl_rnga.h>
00050 #include <gsl/gsl_math.h>
00051 #include <gsl/gsl_monte.h>
00052 #include <gsl/gsl_machine.h>
00053 #include <gsl/gsl_monte_miser.h>
00054 
00055 #ifndef DOXYGENP
00056 namespace o2scl {
00057 #endif
00058   
00059   /** \brief Multidimensional integration using the MISER Monte Carlo 
00060       algorithm (GSL)
00061 
00062       This class uses recursive stratified sampling to estimate the
00063       value of an integral over a hypercubic region.
00064 
00065       By default the minimum number of calls to estimate the variance
00066       is 16 times the number of dimensions. This ratio is
00067       stored in \ref calls_per_dim. By default the minimum number of
00068       calls per bisection is 32 times \ref calls_per_dim times
00069       the number of dimensions. This ratio is stored in 
00070       \ref bisection_ratio. These ratios are employed by
00071       \ref minteg_err().
00072 
00073       Alternatively, the user can directly set these minimums by \ref
00074       set_min_calls() and \ref set_min_calls_per_bisection() and use
00075       \ref miser_minteg_err() which ignores \ref calls_per_dim and
00076       \ref bisection_ratio. 
00077 
00078       If \ref mcarlo_inte::verbose is greater than zero, then the
00079       status of the integration is output at every level of bisection
00080       less than \ref n_levels_out. If it is greater than 1, then the
00081       boundaries of the current region are also output. Finally, if it
00082       is greater than 2, a keypress is required after each output.
00083 
00084       Based on \ref Press90 .
00085   */
00086   template<class func_t=multi_funct<>, class rng_t=gsl_rnga, 
00087     class vec_t=ovector_base, class alloc_vec_t=ovector,
00088     class alloc_t=ovector_alloc> class gsl_miser : 
00089   public mcarlo_inte<func_t,rng_t,vec_t> {
00090 
00091   public:
00092   
00093   /** \brief Number of calls per dimension (default 16)
00094    */
00095   size_t calls_per_dim;
00096   
00097   /** \brief Factor to set \ref min_calls_per_bisection (default 32)
00098    */
00099   size_t bisection_ratio;
00100   
00101   /** \brief Introduce random variation into bisection (default 0.0)
00102       
00103       From GSL documentation:
00104       \verbatim
00105       This parameter introduces a random fractional variation of size
00106       DITHER into each bisection, which can be used to break the
00107       symmetry of integrands which are concentrated near the exact
00108       center of the hypercubic integration region.  The default value of
00109       dither is zero, so no variation is introduced. If needed, a
00110       typical value of DITHER is 0.1.
00111       \endverbatim
00112   */
00113   double dither;
00114 
00115   /** \brief Specify fraction of function calls for estimating variance
00116       (default 0.1)
00117         
00118       From GSL documentation:
00119       \verbatim
00120       This parameter specifies the fraction of the currently available
00121       number of function calls which are allocated to estimating the
00122       variance at each recursive step. The default value is 0.1.
00123       \endverbatim
00124   */
00125   double estimate_frac;
00126 
00127   /** \brief How estimated variances for two sub-regions are combined
00128       (default 2.0)
00129 
00130       The error handler will be called if this is less than zero.
00131 
00132       From GSL documentation:
00133       \verbatim
00134       This parameter controls how the estimated variances for the two
00135       sub-regions of a bisection are combined when allocating points.
00136       With recursive sampling the overall variance should scale better
00137       than 1/N, since the values from the sub-regions will be obtained
00138       using a procedure which explicitly minimizes their variance.  To
00139       accommodate this behavior the MISER algorithm allows the total
00140       variance to depend on a scaling parameter \alpha,
00141         
00142       \Var(f) = {\sigma_a \over N_a^\alpha} + {\sigma_b \over N_b^\alpha}.
00143         
00144       The authors of the original paper describing MISER recommend the
00145       value \alpha = 2 as a good choice, obtained from numerical
00146       experiments, and this is used as the default value in this
00147       implementation.
00148       \endverbatim
00149   */
00150   double alpha;
00151 
00152   /** \brief Minimum number of calls to estimate the variance
00153 
00154       This is set by minteg() and minteg_err() to be \ref
00155       calls_per_dim times the number of dimensions in the problem. The
00156       default value of calls_per_dim is 16 (which is the GSL default).
00157 
00158       From GSL documentation:
00159       \verbatim
00160       This parameter specifies the minimum number of function calls
00161       required for each estimate of the variance. If the number of
00162       function calls allocated to the estimate using ESTIMATE_FRAC falls
00163       below MIN_CALLS then MIN_CALLS are used instead.  This ensures
00164       that each estimate maintains a reasonable level of accuracy. 
00165       \endverbatim
00166   */
00167   int set_min_calls(size_t &mc) {
00168     min_calls=mc;
00169     return 0;
00170   }
00171 
00172   /** \brief Minimum number of calls required to proceed with bisection
00173         
00174       This is set by minteg() and minteg_err() to be \ref
00175       calls_per_dim times \ref bisection_ratio times the number of
00176       dimensions in the problem. The default values give 512 times the
00177       number of dimensions (also the GSL default).
00178 
00179       From GSL documentation:
00180       \verbatim
00181       This parameter specifies the minimum number of function calls
00182       required to proceed with a bisection step.  When a recursive step
00183       has fewer calls available than MIN_CALLS_PER_BISECTION it performs
00184       a plain Monte Carlo estimate of the current sub-region and
00185       terminates its branch of the recursion. 
00186       \endverbatim
00187   */
00188   int set_min_calls_per_bisection(size_t &mcb) {
00189     min_calls_per_bisection=mcb;
00190     return 0;
00191   }
00192 
00193   /** \brief The number of recursive levels to output when verbose is 
00194       greater than zero (default 5)
00195    */
00196   size_t n_levels_out;
00197   
00198 #ifndef DOXYGEN_INTERNAL
00199 
00200   protected:
00201 
00202   /// Minimum number of calls to estimate the variance
00203   size_t min_calls;
00204 
00205   /// Minimum number of calls required to proceed with bisection
00206   size_t min_calls_per_bisection;
00207 
00208   /// The number of dimensions of memory allocated
00209   size_t dim;
00210 
00211   /// \name Arrays which contain a value for each dimension 
00212   //@{
00213   /// The current midpoint 
00214   uvector xmid;
00215   /// The left variance
00216   uvector sigma_l;
00217   /// The right variance
00218   uvector sigma_r;
00219   /// The maximum function value in the left half 
00220   uvector fmax_l;
00221   /// The maximum function value in the right half
00222   uvector fmax_r;
00223   /// The minimum function value in the left half
00224   uvector fmin_l;
00225   /// The minimum function value in the right half
00226   uvector fmin_r;
00227   /// The sum in the left half
00228   uvector fsum_l;
00229   /// The sum in the right half 
00230   uvector fsum_r;
00231   /// The sum of the squares in the left half
00232   uvector fsum2_l;
00233   /// The sum of the squares in the right half 
00234   uvector fsum2_r;
00235   /// The number of evaluation points in the left half
00236   uvector_size_t hits_l;
00237   /// The number of evaluation points in the right half
00238   uvector_size_t hits_r;
00239   //@}
00240     
00241   /** \brief Estimate the variance
00242 
00243       \future Remove the reference to GSL_POSINF and replace with a
00244       function parameter.
00245   */
00246   virtual int estimate_corrmc(func_t &func, size_t ndim,
00247                               const vec_t &xl, const vec_t &xu,
00248                               size_t calls, double &res,
00249                               double &err, const uvector &lxmid,
00250                               uvector &lsigma_l, uvector &lsigma_r) {
00251     size_t i, n;
00252       
00253     double m=0.0, q=0.0;
00254     double vol=1.0;
00255 
00256     for (i=0;i<dim;i++) {
00257       vol*=xu[i]-xl[i];
00258       hits_l[i]=hits_r[i]=0;
00259       fsum_l[i]=fsum_r[i]=0.0;
00260       fsum2_l[i]=fsum2_r[i]=0.0;
00261       lsigma_l[i]=lsigma_r[i]=-1;
00262     }
00263 
00264     for (n=0;n<calls;n++) {
00265       double fval;
00266 
00267       unsigned int j=(n/2) % dim;
00268       unsigned int side=(n % 2);
00269 
00270       for (i=0;i<dim;i++) {
00271 
00272         // The equivalent of gsl_rng_uniform_pos()
00273         double z;
00274         do { z=this->def_rng.random(); } while (z==0);
00275           
00276         if (i != j) {
00277           x[i]=xl[i]+z*(xu[i]-xl[i]);
00278         } else {
00279           if (side == 0) {
00280             x[i]=lxmid[i]+z*(xu[i]-lxmid[i]);
00281           } else {
00282             x[i]=xl[i]+z*(lxmid[i]-xl[i]);
00283           }
00284         }
00285       }
00286         
00287       fval=func(ndim,x);
00288         
00289       /* recurrence for mean and variance */
00290       {
00291         double d=fval-m;
00292         m+=d/(n+1.0);
00293         q+=d*d*(n/(n+1.0));
00294       }
00295 
00296       /* compute the variances on each side of the bisection */
00297       for (i=0;i<dim;i++) {
00298         if (x[i] <= lxmid[i]) {
00299           fsum_l[i]+=fval;
00300           fsum2_l[i]+=fval*fval;
00301           hits_l[i]++;
00302         } else {
00303           fsum_r[i]+=fval;
00304           fsum2_r[i]+=fval*fval;
00305           hits_r[i]++;
00306         }
00307       }
00308     }
00309 
00310     for (i=0;i<dim;i++) {
00311       double fraction_l=(lxmid[i]-xl[i])/(xu[i]-xl[i]);
00312 
00313       if (hits_l[i] > 0) {
00314         fsum_l[i] /= hits_l[i];
00315         lsigma_l[i]=sqrt(fsum2_l[i]-fsum_l[i]*fsum_l[i]/hits_l[i]);
00316         lsigma_l[i]*=fraction_l*vol/hits_l[i];
00317       }
00318         
00319       if (hits_r[i] > 0) {
00320         fsum_r[i] /= hits_r[i];
00321         lsigma_r[i]=sqrt(fsum2_r[i]-fsum_r[i]*fsum_r[i]/hits_r[i]);
00322         lsigma_r[i]*=(1-fraction_l)*vol/hits_r[i];
00323       }
00324     }
00325       
00326     res=vol*m;
00327       
00328     if (calls<2) {
00329       err=GSL_POSINF;
00330     } else {
00331       err=vol*sqrt(q/(calls*(calls-1.0)));
00332     }
00333       
00334     return gsl_success;
00335   }
00336 
00337   /// Memory allocator
00338   alloc_t ao;
00339 
00340   /// The most recent integration point
00341   alloc_vec_t x;
00342 
00343 #endif
00344     
00345   public:
00346     
00347   gsl_miser() {
00348     estimate_frac=0.1;
00349     alpha=2.0;
00350     dither=0.0;
00351     min_calls=100;
00352     min_calls_per_bisection=3000;
00353     dim=0;
00354     n_levels_out=5;
00355     calls_per_dim=16;
00356     bisection_ratio=32;
00357   }
00358 
00359   virtual ~gsl_miser() {
00360     free();
00361   }
00362   
00363   /// Allocate memory
00364   virtual int allocate(size_t ldim) {
00365 
00366     if (dim>0) {
00367       free();
00368     }
00369 
00370     if (ldim==0) {
00371       O2SCL_ERR2_RET("Can't allocate zero memory in ",
00372                      "gsl_miser::allocate().",gsl_efailed);
00373     }
00374 
00375     ao.allocate(x,ldim);
00376 
00377     xmid.allocate(ldim);
00378     sigma_l.allocate(ldim);
00379     sigma_r.allocate(ldim);
00380     fmax_l.allocate(ldim);
00381     fmax_r.allocate(ldim);
00382     fmin_l.allocate(ldim);
00383     fmin_r.allocate(ldim);
00384     fsum_l.allocate(ldim);
00385     fsum_r.allocate(ldim);
00386     fsum2_l.allocate(ldim);
00387     fsum2_r.allocate(ldim);
00388     hits_l.allocate(ldim);
00389     hits_r.allocate(ldim);
00390     
00391     dim=ldim;
00392 
00393     return 0;
00394   }
00395 
00396   /// Free allocated memory
00397   virtual int free() {
00398     if (dim>0) {
00399       hits_r.free();
00400       hits_l.free();
00401       fsum2_r.free();
00402       fsum2_l.free();
00403       fsum_r.free();
00404       fsum_l.free();
00405       fmin_r.free();
00406       fmin_l.free();
00407       fmax_r.free();
00408       fmax_l.free();
00409       sigma_r.free();
00410       sigma_l.free();
00411       xmid.free();
00412       ao.free(x);
00413       dim=0;
00414     }
00415     return 0;
00416   }
00417     
00418   /** \brief Integrate function \c func over the hypercube from
00419       \f$ x_i=\mathrm{xl}_i \f$ to \f$ x_i=\mathrm{xu}_i \f$ for
00420       \f$ 0<i< \f$ ndim-1
00421 
00422       \note The values of \ref min_calls and \ref
00423       min_calls_per_bisection should be set before calling this
00424       function. The default values if not set are 100 and 3000
00425       respectively, which correspond to the GSL default setting
00426       for a 6 dimensional problem. 
00427   */
00428   virtual int miser_minteg_err(func_t &func, size_t ndim, const vec_t &xl, 
00429                                const vec_t &xu, size_t calls, size_t level,
00430                                double &res, double &err) {
00431 
00432     if (min_calls==0 || min_calls_per_bisection==0) {
00433       O2SCL_ERR2_RET("Variables min_calls or min_calls_per_bisection ",
00434                      "are zero in gsl_miser::miser_minteg_err().",
00435                      gsl_einval);
00436     }
00437       
00438     size_t n, estimate_calls, calls_l, calls_r;
00439     size_t i;
00440     size_t i_bisect;
00441     int found_best;
00442 
00443     double res_est=0, err_est=0;
00444     double res_r=0, err_r=0, res_l=0, err_l=0;
00445     double xbi_l, xbi_m, xbi_r, s;
00446 
00447     double vol;
00448     double weight_l, weight_r;
00449 
00450     for (i=0;i<dim;i++) {
00451       if (xu[i] <= xl[i]) {
00452         std::string str="Upper limit, "+dtos(xu[i])+", must be greater "+
00453           "than lower limit, "+dtos(xl[i])+", in gsl_miser::"+
00454           "miser_minteg_err().";
00455         O2SCL_ERR_RET(str.c_str(),gsl_einval);
00456       }
00457       if (xu[i]-xl[i]>GSL_DBL_MAX) {
00458         O2SCL_ERR2_RET("Range of integration is too large ",
00459                        "in gsl_miser::miser_minteg_err().",gsl_einval);
00460       }
00461     }
00462 
00463     if (alpha<0) {
00464       std::string str="Parameter alpha, "+dtos(alpha)+", must be non-"+
00465       "negative in gsl_miser::mister_minteg_err().";
00466       O2SCL_ERR_RET(str.c_str(),gsl_einval);
00467     }
00468 
00469     /* [GSL] Compute volume */
00470 
00471     vol=1;
00472 
00473     for (i=0;i<dim;i++) {
00474       vol*=xu[i]-xl[i];
00475     }
00476 
00477     if (calls<min_calls_per_bisection) {
00478       double m=0.0, q=0.0;
00479 
00480       if (calls<2) {
00481         O2SCL_ERR2_RET("Insufficient calls for subvolume ", 
00482                        "in gsl_miser::miser_minteg_err().",gsl_einval);
00483       }
00484 
00485       for (n=0;n<calls;n++) {
00486         /* Choose a random point in the integration region */
00487 
00488         for (i=0;i<dim;i++) {
00489             
00490           // The equivalent of gsl_rng_uniform_pos()
00491           double rdn;
00492           do { rdn=this->def_rng.random(); } while (rdn==0);
00493 
00494           x[i]=xl[i]+rdn*(xu[i]-xl[i]);
00495         }
00496 
00497         {
00498           double fval;
00499           fval=func(ndim,x);
00500             
00501           /* [GSL] recurrence for mean and variance */
00502 
00503           double d=fval-m;
00504           m+=d/(n+1.0);
00505           q+=d*d*(n/(n+1.0));
00506         }
00507       }
00508 
00509       res=vol*m;
00510 
00511       err=vol*sqrt(q/(calls*(calls-1.0)));
00512 
00513       return gsl_success;
00514     }
00515     
00516     // The equivalent of what's in GSL but rewritten with explicit
00517     // typecasts
00518     size_t prod=(size_t)(((double)calls)*estimate_frac);
00519     estimate_calls=(min_calls > prod ? min_calls : prod);
00520       
00521     if (estimate_calls<4*dim) {
00522       O2SCL_ERR2_RET("Insufficient calls to sample all halfspaces ", 
00523                      "in gsl_miser::miser_minteg_err().",gsl_esanity);
00524     }
00525 
00526     // [GSL] Flip coins to bisect the integration region with some fuzz 
00527 
00528     for (i=0;i<dim;i++) {
00529       s=(this->def_rng.random()-0.5) >= 0.0 ? dither : -dither;
00530       xmid[i]=(0.5+s)*xl[i]+(0.5-s)*xu[i];
00531     }
00532 
00533     /* 
00534        [GSL] The idea is to chose the direction to bisect based on
00535        which will give the smallest total variance. We could (and may
00536        do so later) use MC to compute these variances. But the NR guys
00537        simply estimate the variances by finding the min and max
00538        function values for each half-region for each bisection.
00539     */
00540     estimate_corrmc(func,dim,xl,xu,estimate_calls,
00541                     res_est,err_est,xmid,sigma_l,sigma_r);
00542 
00543     // [GSL] We have now used up some calls for the estimation 
00544 
00545     calls -= estimate_calls;
00546 
00547     // [GSL] Now find direction with the smallest total "variance"
00548 
00549     {
00550       double best_var=GSL_DBL_MAX;
00551       double beta=2.0/(1.0+alpha);
00552       found_best=0;
00553       i_bisect=0;
00554       weight_l=weight_r=1.0;
00555 
00556       for (i=0;i<dim;i++) {
00557         if (sigma_l[i] >= 0 && sigma_r[i] >= 0) {
00558           // [GSL] Estimates are okay 
00559           double var=pow (sigma_l[i], beta)+pow (sigma_r[i], beta);
00560             
00561           if (var <= best_var) {
00562             found_best=1;
00563             best_var=var;
00564             i_bisect=i;
00565             weight_l=pow (sigma_l[i], beta);
00566             weight_r=pow (sigma_r[i], beta);
00567             if (weight_l==0 && weight_r==0) {
00568               weight_l=1;
00569               weight_r=1;
00570             }
00571           }
00572         } else {
00573           if (sigma_l[i]<0) {
00574             O2SCL_ERR2_RET("No points in left-half space ", 
00575                            "in gsl_miser::miser_minteg_err().",gsl_esanity);
00576           }
00577           if (sigma_r[i]<0) {
00578             O2SCL_ERR2_RET("No points in right-half space ", 
00579                            "in gsl_miser::miser_minteg_err().",gsl_esanity);
00580           }
00581         }
00582       }
00583     }
00584       
00585     if (!found_best) {
00586       // [GSL] All estimates were the same, so chose a direction at
00587       // random
00588       i_bisect=this->def_rng.random_int(dim);
00589     }
00590 
00591     xbi_l=xl[i_bisect];
00592     xbi_m=xmid[i_bisect];
00593     xbi_r=xu[i_bisect];
00594 
00595     // [GSL] Get the actual fractional sizes of the two "halves", and
00596     // distribute the remaining calls among them 
00597     {
00598       double fraction_l=fabs ((xbi_m-xbi_l)/(xbi_r-xbi_l));
00599       double fraction_r=1-fraction_l;
00600 
00601       double a=fraction_l*weight_l;
00602       double b=fraction_r*weight_r;
00603 
00604       calls_l=(size_t)(min_calls+(calls-2*min_calls)*a/(a+b));
00605       calls_r=(size_t)(min_calls+(calls-2*min_calls)*b/(a+b));
00606     }
00607 
00608     if (this->verbose>0 && level<n_levels_out) {
00609       std::cout << "gsl_miser: level,calls_l,calls_r,calls,min_calls_pb: " 
00610       << level << " " << calls_l << " " << calls_r << " " << calls << " "
00611       << min_calls_per_bisection << std::endl;
00612       std::cout << "\tres,err: " << res_est << " " << err_est << std::endl;
00613       if (this->verbose>1) {
00614         std::cout << "\ti,left,mid,right: " << i_bisect << " "
00615                   << xbi_l << " " << xbi_m << " " << xbi_r << std::endl;
00616         for(size_t j=0;j<dim;j++) {
00617           std::cout << "\t\ti,low,high: " << j << " " << xl[j] << " " << xu[j] 
00618                     << std::endl;
00619         }
00620       }
00621       if (this->verbose>2) {
00622         char ch;
00623         std::cin >> ch;
00624       }
00625     }
00626 
00627     /* [GSL] Compute the integral for the left hand side of the
00628        bisection. Due to the recursive nature of the algorithm we must
00629        allocate some new memory for the integration limits for each
00630        recursive call
00631     */
00632     {
00633       int status;
00634         
00635       alloc_vec_t xu_tmp;
00636       ao.allocate(xu_tmp,dim);
00637         
00638       for (i=0;i<dim;i++) {
00639         xu_tmp[i]=xu[i];
00640       }
00641         
00642       xu_tmp[i_bisect]=xbi_m;
00643         
00644       status=miser_minteg_err(func,dim,xl,xu_tmp,calls_l,level+1,
00645                               res_l,err_l);
00646 
00647       ao.free(xu_tmp);
00648 
00649       if (status != gsl_success) {
00650         return status;
00651       }
00652     }
00653 
00654     // [GSL] Compute the integral for the right hand side of the
00655     // bisection 
00656     {
00657       int status;
00658 
00659       alloc_vec_t xl_tmp;
00660       ao.allocate(xl_tmp,dim);
00661 
00662       for (i=0;i<dim;i++) {
00663         xl_tmp[i]=xl[i];
00664       }
00665 
00666       xl_tmp[i_bisect]=xbi_m;
00667 
00668       status=miser_minteg_err(func,dim,xl_tmp,xu,calls_r,level+1,
00669                               res_r,err_r);
00670         
00671       ao.free(xl_tmp);
00672 
00673       if (status != gsl_success) {
00674         return status;
00675       }
00676     }
00677 
00678     res=res_l+res_r;
00679     err=sqrt(err_l*err_l+err_r*err_r);
00680 
00681     return 0;
00682   }
00683     
00684   virtual ~gsl_miser() {}
00685   
00686   /** \brief Integrate function \c func from x=a to x=b.
00687 
00688       This function is just a wrapper to miser_minteg_err() which
00689       allocates the memory if necessary, sets \c min_calls and \c
00690       min_calls_per_bisection, calls \ref miser_minteg_err(), and then
00691       frees the previously allocated memory.
00692   */
00693   virtual int minteg_err(func_t &func, size_t ndim, const vec_t &a, 
00694                          const vec_t &b, double &res, double &err) {
00695     if (ndim!=dim) allocate(ndim);
00696     min_calls=calls_per_dim*ndim;
00697     min_calls_per_bisection=bisection_ratio*min_calls;
00698     int ret=miser_minteg_err(func,ndim,a,b,this->n_points,0,res,err);
00699     // Set these to back to zero to ensure that if 
00700     // miser_minteg_err() is called directly, the user sets 
00701     // min_calls and min_calls_per_bisection first.
00702     min_calls=0;
00703     min_calls_per_bisection=0;
00704     return ret;
00705   }
00706     
00707   /** \brief Integrate function \c func over the hypercube from
00708       \f$ x_i=a_i \f$ to \f$ x_i=b_i \f$ for
00709       \f$ 0<i< \f$ ndim-1
00710         
00711       This function is just a wrapper to minteg_err() which allocates
00712       the memory, sets min_calls and min_calls_per_bisection, calls
00713       miser_minteg_err(), and then frees the previously allocated
00714       memory.
00715   */
00716   virtual double minteg(func_t &func, size_t ndim, const vec_t &a, 
00717                         const vec_t &b) {
00718     double res;
00719     minteg_err(func,ndim,a,b,res,this->interror);
00720     return res;
00721   }
00722 
00723   /// Return string denoting type ("gsl_miser")
00724   virtual const char *type() { return "gsl_miser"; }
00725     
00726   };
00727   
00728 #ifndef DOXYGENP
00729 }
00730 #endif
00731 
00732 #endif
00733 
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.