00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 #ifdef ROOT_FIT_PARALLEL
00014 
00015 #include "Fit/FitUtilParallel.h"
00016 
00017 #include "Fit/BinData.h"
00018 #include "Fit/UnBinData.h"
00019 #include "Fit/FitUtil.h"
00020 
00021 #include "Math/IParamFunction.h"
00022 
00023 int ncalls = 0; 
00024 
00025 #include <limits>
00026 #include <cmath>
00027 #include <numeric>
00028 
00029 
00030 #ifdef DEBUG
00031 #include <iostream> 
00032 #endif
00033 
00034 #ifdef USE_PTHREAD
00035 
00036 #include <pthread.h>
00037 
00038 
00039 #define NUMBER_OF_THREADS 2
00040 
00041 #else 
00042 #include <omp.h>
00043 
00044 #define MAX_THREAD 8
00045 
00046 #endif
00047 
00048 namespace ROOT { 
00049 
00050    namespace Fit { 
00051 
00052       namespace FitUtilParallel { 
00053 
00054 #ifdef USE_PTHREAD
00055 
00056 class ThreadData { 
00057 
00058 public: 
00059 
00060    ThreadData() : 
00061       fBegin(0),fEnd(0) ,
00062       fData(0),
00063       fFunc(0)
00064    {}
00065    
00066    ThreadData(unsigned int i1, unsigned int i2, const BinData & data, IModelFunction &func) : 
00067       fBegin(i1), fEnd(i2), 
00068       fData(&data),
00069       fFunc(&func)
00070    {}
00071 
00072    void Set(unsigned int nrej, double sum) { 
00073       fNRej = nrej; 
00074       fSum = sum; 
00075    }
00076 
00077    const BinData & Data() const { return *fData; }
00078 
00079    IModelFunction & Func() { return *fFunc; }
00080 
00081    double Sum() const { return fSum; } 
00082 
00083    unsigned int NRej() const { return fNRej; } 
00084 
00085    unsigned int Begin() const { return fBegin; } 
00086    unsigned int End() const { return fEnd; } 
00087    
00088 private: 
00089 
00090    const unsigned int fBegin; 
00091    const unsigned int fEnd; 
00092    const BinData * fData; 
00093    IModelFunction * fFunc; 
00094    double fSum;
00095    unsigned int fNRej; 
00096 };
00097 
00098 
00099 
00100 void *EvaluateResidual(void * ptr) { 
00101 
00102    ThreadData * t = (ThreadData *) ptr; 
00103 
00104    unsigned int istart = t->Begin();
00105    unsigned int iend = t->End(); 
00106    double chi2 = 0;
00107    unsigned int nRejected = 0;
00108    const int nthreads = NUMBER_OF_THREADS; 
00109 
00110    const BinData & data = t->Data();
00111    IModelFunction & func = t->Func();
00112    for (unsigned int i = istart; i < iend; i+=nthreads) { 
00113       const double * x = data.Coords(i); 
00114       double y = data.Value(i);
00115       double invError = data.InvError(i);
00116       double fval = 0; 
00117       fval = func ( x ); 
00118 
00119 
00120 
00121 
00122 
00123 
00124 
00125 
00126       
00127       
00128       if (fval > - std::numeric_limits<double>::max() && fval < std::numeric_limits<double>::max() ) { 
00129          
00130          double tmp = ( y -fval )* invError;      
00131          chi2 += tmp*tmp;
00132       }
00133       else 
00134          nRejected++; 
00135       
00136    }
00137 
00138 #ifdef DEBUG
00139    std::cout << "end loop " << istart << "  " << iend << " chi2 = " << chi2 << " nrej " << nRejected << std::endl; 
00140 #endif
00141    t->Set(nRejected,chi2);
00142    return 0;
00143 }
00144 
00145 double EvaluateChi2(IModelFunction & func, const BinData & data, const double * p, unsigned int & nPoints) {  
00146    
00147    
00148 
00149    const int nthreads = NUMBER_OF_THREADS; 
00150    
00151    unsigned int n = data.Size();
00152 
00153 #ifdef DEBUG
00154    std::cout << "\n\nFit data size = " << n << std::endl;
00155 #endif
00156 
00157    func.SetParameters(p); 
00158 
00159    
00160    pthread_t   thread[nthreads];
00161    ThreadData *  td[nthreads]; 
00162    unsigned int istart = 0; 
00163    for (int ithread = 0; ithread < nthreads; ++ithread) { 
00164 
00165 
00166 
00167       int iend = n; 
00168       istart = ithread;  
00169       td[ithread] = new ThreadData(istart,iend,data,func);
00170       pthread_create(&thread[ithread], NULL, EvaluateResidual, td[ithread]); 
00171       
00172    }
00173 
00174    for (int ithread = 0; ithread < nthreads; ++ithread)  
00175       pthread_join(thread[ithread], NULL); 
00176 
00177    
00178 
00179    double chi2 = 0; 
00180    int nRejected = 0;
00181    for (int ithread = 0; ithread < nthreads; ++ithread) { 
00182       nRejected += td[ithread]->NRej();
00183       chi2 += td[ithread]->Sum();
00184       delete td[ithread]; 
00185    }
00186    
00187 
00188 #ifdef DEBUG
00189    std::cout << "chi2 = " << chi2 << " n = " << nRejected << std::endl;
00190 #endif
00191    
00192    nPoints = n - nRejected;
00193    return chi2;
00194 
00195 }
00196 
00197 
00198 
00199 #else 
00200 
00201 double EvaluateChi2(IModelFunction & func, const BinData & data, const double * p, unsigned int & nPoints) { 
00202    return FitUtil::EvaluateChi2(func,data,p,nPoints); 
00203 }
00204 
00205 
00206 
00207 
00208 inline double EvalLogF(double fval) { 
00209    
00210    
00211    
00212    const static double epsilon = 2.*std::numeric_limits<double>::min();
00213    if(fval<= epsilon) 
00214       return fval/epsilon + std::log(epsilon) - 1; 
00215    else      
00216       return std::log(fval);
00217 }
00218 
00219 
00220 double EvaluateLogL(IModelFunction & func, const UnBinData & data, const double * p, unsigned int &nPoints) {  
00221    
00222 
00223    unsigned int n = data.Size();
00224 
00225 #ifdef DEBUG
00226    std::cout << "\n\nFit data size = " << n << std::endl;
00227    
00228    std::cout << "\tpar = [ " << func.NPar() << " ] =  "; 
00229    for (unsigned int ipar = 0; ipar < func.NPar(); ++ipar) 
00230       std::cout << p[ipar] << ", ";
00231    std::cout <<"---------------------------\n";
00232 #endif
00233 
00234    double logl = 0;
00235    int nRejected = 0; 
00236 
00237 
00238    
00239 
00240 
00241 #pragma omp parallel
00242 #pragma omp for reduction (+:logl,nRejected)
00243 
00244 
00245 
00246    for (unsigned int i = 0; i < n; ++ i) { 
00247 
00248       
00249       
00250 
00251       const double * x = data.Coords(i);
00252       double fval = func ( x, p  ); 
00253 
00254       if (fval < 0) { 
00255          nRejected++; 
00256       }
00257       else 
00258          
00259          logl += EvalLogF( fval); 
00260 
00261 #ifdef DEBUG      
00262 #pragma omp critical 
00263       {     std::cout << " ==== i = " << i << " thread " << omp_get_thread_num() 
00264                       << "fval = " << fval << " logl = " << logl << std::endl;}
00265 
00266 
00267 
00268 #endif
00269    }
00270    
00271    
00272    if (nRejected != 0)  nPoints = n - nRejected;
00273 
00274 #ifdef DEBUG
00275    ncalls++;
00276    int pr = std::cout.precision(15);
00277    std::cout << "ncalls " << ncalls << " Logl = " << logl << " np = " << nPoints << std::endl;
00278    std::cout.precision(pr);
00279    assert(ncalls<3);
00280 #endif
00281 
00282    
00283 
00284    double result = - logl; 
00285    return result;
00286 }
00287 
00288 #endif
00289 
00290 
00291       } 
00292 
00293    } 
00294 
00295 } 
00296 
00297 
00298 
00299 
00300 #endif