Object-oriented Scientific Computing Library: Version 0.910
gsl_monte.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/plain.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_MONTE_H
00043 #define O2SCL_GSL_MONTE_H
00044 
00045 #include <iostream>
00046 #include <o2scl/mcarlo_inte.h>
00047 #include <o2scl/gsl_rnga.h>
00048 #include <o2scl/multi_funct.h>
00049 #include <gsl/gsl_monte.h>
00050 #include <gsl/gsl_monte_plain.h>
00051 
00052 #ifndef DOXYGENP
00053 namespace o2scl {
00054 #endif
00055 
00056   /** \brief Multidimensional integration using plain Monte Carlo (GSL)
00057   */
00058   template<class func_t=multi_funct<>, class rng_t=gsl_rnga, 
00059     class vec_t=ovector_base, class alloc_vec_t=ovector,
00060     class alloc_t=ovector_alloc> class gsl_monte : 
00061   public mcarlo_inte<func_t,rng_t,vec_t> {
00062 
00063   public:
00064 
00065   virtual ~gsl_monte() {}
00066 
00067   /// Integrate function \c func from x=a to x=b.
00068   virtual int minteg_err(func_t &func, size_t ndim, const vec_t &a, 
00069                          const vec_t &b, double &res, double &err) {
00070 
00071     double r;
00072       
00073     alloc_vec_t cc;
00074     alloc_t ao;
00075     ao.allocate(cc,ndim);
00076       
00077     double vol=1.0, m=0.0, q=0.0;
00078     for(size_t i=0;i<ndim;i++) vol*=b[i]-a[i];
00079       
00080     for(size_t n=0;n<this->n_points;n++) {
00081       for(size_t i=0;i<ndim;i++) {
00082         do {
00083           r=this->def_rng.random();
00084         } while (r==0.0);
00085         cc[i]=a[i]+r*(b[i]-a[i]);
00086       }
00087       double fval;
00088       fval=func(ndim,cc);
00089       double d=fval-m;
00090       m+=d/(n+1.0);
00091       q+=d*d*(n/(n+1.0));
00092     }
00093 
00094     res=vol*m;
00095       
00096     if (this->n_points<2) {
00097       err=0.0;
00098     } else {
00099       err=vol*sqrt(q/(this->n_points*(this->n_points-1.0)));
00100     }
00101           
00102     ao.free(cc);
00103 
00104     return 0;
00105 
00106   }
00107       
00108   /** \brief Integrate function \c func over the hypercube from
00109       \f$ x_i=a_i \f$ to \f$ x_i=b_i \f$ for \f$ 0<i< \f$ ndim-1
00110   */
00111   virtual double minteg(func_t &func, size_t ndim, const vec_t &a, 
00112                         const vec_t &b) {
00113     double res;
00114     minteg_err(func,ndim,a,b,res,this->interror);
00115     return res;
00116   }
00117 
00118   /// Return string denoting type ("gsl_monte")
00119   virtual const char *type() { return "gsl_monte"; }
00120       
00121   };
00122   
00123 #ifndef DOXYGENP
00124 }
00125 #endif
00126 
00127 #endif
00128 
 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.