MCMCInterval.cxx

Go to the documentation of this file.
00001 // @(#)root/roostats:$Id: MCMCInterval.cxx 35810 2010-09-27 20:05:22Z moneta $
00002 // Authors: Kevin Belasco        17/06/2009
00003 // Authors: Kyle Cranmer         17/06/2009
00004 /*************************************************************************
00005  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //_________________________________________________
00013 /*
00014 BEGIN_HTML
00015 <p>
00016 MCMCInterval is a concrete implementation of the RooStats::ConfInterval
00017 interface.  It takes as input Markov Chain of data points in the parameter
00018 space generated by Monte Carlo using the Metropolis algorithm.  From the Markov
00019 Chain, the confidence interval can be determined in two ways:
00020 </p>
00021 
00022 <p>
00023 Using a Kernel-Estimated PDF: (not the default method)
00024 </p>
00025 <p>
00026 A RooNDKeysPdf is constructed from the data set using adaptive kernel width.
00027 With this RooNDKeysPdf F, we then integrate over the most likely domain in the
00028 parameter space (tallest points in the posterior RooNDKeysPdf) until the target
00029 confidence level is reached within an acceptable neighborhood as defined by
00030 SetEpsilon(). More specifically: we calculate the following for different
00031 cutoff values C until we reach the target confidence level: \int_{ F >= C } F
00032 d{normset}.
00033 Important note: this is not the default method because of a bug in constructing
00034 the RooNDKeysPdf from a weighted data set.  Configure to use this method by
00035 calling SetUseKeys(true), and the data set will be interpreted without weights.
00036 </p>
00037 
00038 <p>
00039 Using a binned data set: (the default method)
00040 </p>
00041 This is the binned analog of the continuous integrative method that uses the
00042 kernel-estimated PDF.  The points in the Markov Chain are put into a binned
00043 data set and the interval is then calculated by adding the heights of the bins
00044 in decreasing order until the desired level of confidence has been reached.
00045 Note that this means the actual confidence level is >= the confidence level
00046 prescribed by the client (unless the user calls SetHistStrict(kFALSE)).  This
00047 method is the default but may not remain as such in future releases, so you may
00048 wish to explicitly configure to use this method by calling SetUseKeys(false)
00049 </p>
00050 
00051 <p>
00052 These are not the only ways for the confidence interval to be determined, and
00053 other possibilities are being considered being added, especially for the
00054 1-dimensional case.
00055 </p>
00056 
00057 <p>
00058 One can ask an MCMCInterval for the lower and upper limits on a specific
00059 parameter of interest in the interval.  Note that this works better for some
00060 distributions (ones with exactly one local maximum) than others, and sometimes
00061 has little value.
00062 </p>
00063 END_HTML
00064 */
00065 //_________________________________________________
00066 
00067 #ifndef ROOT_Rtypes
00068 #include "Rtypes.h"
00069 #endif
00070 
00071 #ifndef ROOT_TMath
00072 #include "TMath.h"
00073 #endif
00074 
00075 #ifndef RooStats_MCMCInterval
00076 #include "RooStats/MCMCInterval.h"
00077 #endif
00078 #ifndef ROOSTATS_MarkovChain
00079 #include "RooStats/MarkovChain.h"
00080 #endif
00081 #ifndef RooStats_Heaviside
00082 #include "RooStats/Heaviside.h"
00083 #endif
00084 #ifndef ROO_DATA_HIST
00085 #include "RooDataHist.h"
00086 #endif
00087 #ifndef ROO_KEYS_PDF
00088 #include "RooNDKeysPdf.h"
00089 #endif
00090 #ifndef ROO_PRODUCT
00091 #include "RooProduct.h"
00092 #endif
00093 #ifndef RooStats_RooStatsUtils
00094 #include "RooStats/RooStatsUtils.h"
00095 #endif
00096 #ifndef ROO_REAL_VAR
00097 #include "RooRealVar.h"
00098 #endif
00099 #ifndef ROO_ARG_LIST
00100 #include "RooArgList.h"
00101 #endif
00102 #ifndef ROOT_TIterator
00103 #include "TIterator.h"
00104 #endif
00105 #ifndef ROOT_TH1
00106 #include "TH1.h"
00107 #endif
00108 #ifndef ROOT_TH1F
00109 #include "TH1F.h"
00110 #endif
00111 #ifndef ROOT_TH2F
00112 #include "TH2F.h"
00113 #endif
00114 #ifndef ROOT_TH3F
00115 #include "TH3F.h"
00116 #endif
00117 #ifndef ROO_MSG_SERVICE
00118 #include "RooMsgService.h"
00119 #endif
00120 #ifndef ROO_GLOBAL_FUNC
00121 #include "RooGlobalFunc.h"
00122 #endif
00123 #ifndef ROOT_TObject
00124 #include "TObject.h"
00125 #endif
00126 #ifndef ROOT_THnSparse
00127 #include "THnSparse.h"
00128 #endif
00129 #ifndef ROO_NUMBER
00130 #include "RooNumber.h"
00131 #endif
00132 //#ifndef ROOT_TFile
00133 //#include "TFile.h"
00134 //#endif
00135 
00136 #include <cstdlib>
00137 #include <string>
00138 #include <algorithm>
00139 
00140 ClassImp(RooStats::MCMCInterval);
00141 
00142 using namespace RooFit;
00143 using namespace RooStats;
00144 using namespace std;
00145 
00146 static const Double_t DEFAULT_EPSILON = 0.01;
00147 static const Double_t DEFAULT_DELTA   = 10e-6;
00148 
00149 MCMCInterval::MCMCInterval(const char* name)
00150    : ConfInterval(name)
00151 {
00152    fConfidenceLevel = 0.0;
00153    fHistConfLevel = 0.0;
00154    fKeysConfLevel = 0.0;
00155    fTFConfLevel = 0.0;
00156    fFull = 0.0;
00157    fChain = NULL;
00158    fAxes = NULL;
00159    fDataHist = NULL;
00160    fSparseHist = NULL;
00161    fVector.clear();
00162    fKeysPdf = NULL;
00163    fProduct = NULL;
00164    fHeaviside = NULL;
00165    fKeysDataHist = NULL;
00166    fCutoffVar = NULL;
00167    fHist = NULL;
00168    fNumBurnInSteps = 0;
00169    fHistCutoff = -1;
00170    fKeysCutoff = -1;
00171    fDimension = 1;
00172    fUseKeys = kFALSE;
00173    fUseSparseHist = kFALSE;
00174    fIsHistStrict = kTRUE;
00175    fEpsilon = DEFAULT_EPSILON;
00176    fDelta = DEFAULT_DELTA;
00177    fIntervalType = kShortest;
00178    fTFLower = -1.0 * RooNumber::infinity();
00179    fTFUpper = RooNumber::infinity();
00180    fVecWeight = 0;
00181    fLeftSideTF = -1;
00182 }
00183 
00184 MCMCInterval::MCMCInterval(const char* name,
00185         const RooArgSet& parameters, MarkovChain& chain) : ConfInterval(name)
00186 {
00187    fNumBurnInSteps = 0;
00188    fConfidenceLevel = 0.0;
00189    fHistConfLevel = 0.0;
00190    fKeysConfLevel = 0.0;
00191    fTFConfLevel = 0.0;
00192    fFull = 0.0;
00193    fAxes = NULL;
00194    fChain = &chain;
00195    fDataHist = NULL;
00196    fSparseHist = NULL;
00197    fVector.clear();
00198    fKeysPdf = NULL;
00199    fProduct = NULL;
00200    fHeaviside = NULL;
00201    fKeysDataHist = NULL;
00202    fCutoffVar = NULL;
00203    fHist = NULL;
00204    fHistCutoff = -1;
00205    fKeysCutoff = -1;
00206    fUseKeys = kFALSE;
00207    fUseSparseHist = kFALSE;
00208    fIsHistStrict = kTRUE;
00209    fEpsilon = DEFAULT_EPSILON;
00210    SetParameters(parameters);
00211    fDelta = DEFAULT_DELTA;
00212    fIntervalType = kShortest;
00213    fTFLower = -1.0 * RooNumber::infinity();
00214    fTFUpper = RooNumber::infinity();
00215    fVecWeight = 0;
00216    fLeftSideTF = -1;
00217 }
00218 
00219 MCMCInterval::~MCMCInterval()
00220 {
00221    // destructor
00222    delete[] fAxes;
00223    delete fHist;
00224    delete fChain;
00225    // kbelasco: check here for memory management errors
00226    delete fDataHist;
00227    delete fSparseHist;
00228    delete fKeysPdf;
00229    delete fProduct;
00230    delete fHeaviside;
00231    delete fKeysDataHist;
00232    delete fCutoffVar;
00233 }
00234 
00235 struct CompareDataHistBins {
00236    CompareDataHistBins(RooDataHist* hist) : fDataHist(hist) {}
00237    bool operator() (Int_t bin1 , Int_t bin2) { 
00238       fDataHist->get(bin1);
00239       Double_t n1 = fDataHist->weight();
00240       fDataHist->get(bin2);
00241       Double_t n2 = fDataHist->weight();
00242       return (n1 < n2);
00243    }
00244    RooDataHist* fDataHist; 
00245 };
00246 
00247 struct CompareSparseHistBins {
00248    CompareSparseHistBins(THnSparse* hist) : fSparseHist(hist) {}
00249    bool operator() (Long_t bin1, Long_t bin2) { 
00250       Double_t n1 = fSparseHist->GetBinContent(bin1);
00251       Double_t n2 = fSparseHist->GetBinContent(bin2);
00252       return (n1 < n2);
00253    }
00254    THnSparse* fSparseHist; 
00255 };
00256 
00257 struct CompareVectorIndices {
00258    CompareVectorIndices(MarkovChain* chain, RooRealVar* param) :
00259       fChain(chain), fParam(param) {}
00260    bool operator() (Int_t i, Int_t j) {
00261       Double_t xi = fChain->Get(i)->getRealValue(fParam->GetName());
00262       Double_t xj = fChain->Get(j)->getRealValue(fParam->GetName());
00263       return (xi < xj);
00264    }
00265    MarkovChain* fChain;
00266    RooRealVar* fParam;
00267 };
00268 
00269 // kbelasco: for this method, consider running DetermineInterval() if 
00270 // fKeysPdf==NULL, fSparseHist==NULL, fDataHist==NULL, or fVector.size()==0
00271 // rather than just returning false.  Though this should not be an issue
00272 // because nobody should be able to get an MCMCInterval that has their interval
00273 // or posterior representation NULL/empty since they should only get this
00274 // through the MCMCCalculator
00275 Bool_t MCMCInterval::IsInInterval(const RooArgSet& point) const 
00276 {
00277    if (fIntervalType == kShortest) {
00278       if (fUseKeys) {
00279          if (fKeysPdf == NULL)
00280             return false;
00281 
00282          // evaluate keyspdf at point and return whether >= cutoff
00283          RooStats::SetParameters(&point, const_cast<RooArgSet *>(&fParameters));
00284          return fKeysPdf->getVal(&fParameters) >= fKeysCutoff;
00285       } else {
00286          if (fUseSparseHist) {
00287             if (fSparseHist == NULL)
00288                return false;
00289 
00290             // evalute sparse hist at bin where point lies and return
00291             // whether >= cutoff
00292             RooStats::SetParameters(&point,
00293                                     const_cast<RooArgSet*>(&fParameters));
00294             Long_t bin;
00295             // kbelasco: consider making x static
00296             Double_t* x = new Double_t[fDimension];
00297             for (Int_t i = 0; i < fDimension; i++)
00298                x[i] = fAxes[i]->getVal();
00299             bin = fSparseHist->GetBin(x, kFALSE);
00300             Double_t weight = fSparseHist->GetBinContent((Long64_t)bin);
00301             delete[] x;
00302             return (weight >= (Double_t)fHistCutoff);
00303          } else {
00304             if (fDataHist == NULL)
00305                return false;
00306 
00307             // evaluate data hist at bin where point lies and return whether
00308             // >= cutoff
00309             Int_t bin;
00310             bin = fDataHist->getIndex(point);
00311             fDataHist->get(bin);
00312             return (fDataHist->weight() >= (Double_t)fHistCutoff);
00313          }
00314       }
00315    } else if (fIntervalType == kTailFraction) {
00316       if (fVector.size() == 0)
00317          return false;
00318 
00319       // return whether value of point is within the range
00320       Double_t x = point.getRealValue(fAxes[0]->GetName());
00321       if (fTFLower <= x && x <= fTFUpper)
00322          return true;
00323 
00324       return false;
00325    }
00326 
00327    coutE(InputArguments) << "Error in MCMCInterval::IsInInterval: "
00328       << "Interval type not set.  Returning false." << endl;
00329    return false;
00330 }
00331 
00332 void MCMCInterval::SetConfidenceLevel(Double_t cl)
00333 {
00334    fConfidenceLevel = cl;
00335    DetermineInterval();
00336 }
00337 
00338 // kbelasco: update this or just take it out
00339 // kbelasco: consider keeping this around but changing the implementation
00340 // to set the number of bins for each RooRealVar and then reacreating the
00341 // histograms
00342 //void MCMCInterval::SetNumBins(Int_t numBins)
00343 //{
00344 //   if (numBins > 0) {
00345 //      fPreferredNumBins = numBins;
00346 //      for (Int_t d = 0; d < fDimension; d++)
00347 //         fNumBins[d] = numBins;
00348 //   }
00349 //   else {
00350 //      coutE(Eval) << "* Error in MCMCInterval::SetNumBins: " <<
00351 //                     "Negative number of bins given: " << numBins << endl;
00352 //      return;
00353 //   }
00354 //
00355 //   // If the histogram already exists, recreate it with the new bin numbers
00356 //   if (fHist != NULL)
00357 //      CreateHist();
00358 //}
00359 
00360 void MCMCInterval::SetAxes(RooArgList& axes)
00361 {
00362    Int_t size = axes.getSize();
00363    if (size != fDimension) {
00364       coutE(InputArguments) << "* Error in MCMCInterval::SetAxes: " <<
00365                                "number of variables in axes (" << size <<
00366                                ") doesn't match number of parameters ("
00367                                << fDimension << ")" << endl;
00368       return;
00369    }
00370    for (Int_t i = 0; i < size; i++)
00371       fAxes[i] = (RooRealVar*)axes.at(i);
00372 }
00373 
00374 void MCMCInterval::CreateKeysPdf()
00375 {
00376    // kbelasco: check here for memory leak.  does RooNDKeysPdf use
00377    // the RooArgList passed to it or does it make a clone?
00378    // also check for memory leak from chain, does RooNDKeysPdf clone that?
00379    if (fAxes == NULL || fParameters.getSize() == 0) {
00380       coutE(InputArguments) << "Error in MCMCInterval::CreateKeysPdf: "
00381          << "parameters have not been set." << endl;
00382       return;
00383    }
00384 
00385    if (fNumBurnInSteps >= fChain->Size()) {
00386       coutE(InputArguments) <<
00387          "MCMCInterval::CreateKeysPdf: creation of Keys PDF failed: " <<
00388          "Number of burn-in steps (num steps to ignore) >= number of steps " <<
00389          "in Markov chain." << endl;
00390       delete fKeysPdf;
00391       delete fCutoffVar;
00392       delete fHeaviside;
00393       delete fProduct;
00394       fKeysPdf = NULL;
00395       fCutoffVar = NULL;
00396       fHeaviside = NULL;
00397       fProduct = NULL;
00398       return;
00399    }
00400 
00401    RooDataSet* chain = fChain->GetAsDataSet(SelectVars(fParameters),
00402          EventRange(fNumBurnInSteps, fChain->Size()));
00403    RooArgList* paramsList = new RooArgList();
00404    for (Int_t i = 0; i < fDimension; i++)
00405       paramsList->add(*fAxes[i]);
00406 
00407    // kbelasco: check for memory leaks with chain.  who owns it? does
00408    // RooNDKeysPdf take ownership?
00409    fKeysPdf = new RooNDKeysPdf("keysPDF", "Keys PDF", *paramsList, *chain, "a");
00410    fCutoffVar = new RooRealVar("cutoff", "cutoff", 0);
00411    fHeaviside = new Heaviside("heaviside", "Heaviside", *fKeysPdf, *fCutoffVar);
00412    fProduct = new RooProduct("product", "Keys PDF & Heaviside Product",
00413                                         RooArgSet(*fKeysPdf, *fHeaviside));
00414 }
00415 
00416 void MCMCInterval::CreateHist()
00417 {
00418    if (fAxes == NULL || fChain == NULL) {
00419       coutE(Eval) << "* Error in MCMCInterval::CreateHist(): " <<
00420                      "Crucial data member was NULL." << endl;
00421       coutE(Eval) << "Make sure to fully construct/initialize." << endl;
00422       return;
00423    }
00424    if (fHist != NULL)
00425       delete fHist;
00426 
00427    if (fNumBurnInSteps >= fChain->Size()) {
00428       coutE(InputArguments) <<
00429          "MCMCInterval::CreateHist: creation of histogram failed: " <<
00430          "Number of burn-in steps (num steps to ignore) >= number of steps " <<
00431          "in Markov chain." << endl;
00432       fHist = NULL;
00433       return;
00434    }
00435 
00436    if (fDimension == 1)
00437       fHist = new TH1F("posterior", "MCMC Posterior Histogram",
00438             fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax());
00439 
00440    else if (fDimension == 2)
00441       fHist = new TH2F("posterior", "MCMC Posterior Histogram",
00442             fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax(),
00443             fAxes[1]->numBins(), fAxes[1]->getMin(), fAxes[1]->getMax());
00444 
00445    else if (fDimension == 3) 
00446       fHist = new TH3F("posterior", "MCMC Posterior Histogram",
00447             fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax(),
00448             fAxes[1]->numBins(), fAxes[1]->getMin(), fAxes[1]->getMax(),
00449             fAxes[2]->numBins(), fAxes[2]->getMin(), fAxes[2]->getMax());
00450 
00451    else {
00452       coutE(Eval) << "* Error in MCMCInterval::CreateHist() : " <<
00453                      "TH1* couldn't handle dimension: " << fDimension << endl;
00454       return;
00455    }
00456 
00457    // Fill histogram
00458    Int_t size = fChain->Size();
00459    const RooArgSet* entry;
00460    for (Int_t i = fNumBurnInSteps; i < size; i++) {
00461       entry = fChain->Get(i);
00462       if (fDimension == 1)
00463          ((TH1F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
00464                               fChain->Weight());
00465       else if (fDimension == 2)
00466          ((TH2F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
00467                               entry->getRealValue(fAxes[1]->GetName()),
00468                               fChain->Weight());
00469       else
00470          ((TH3F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
00471                               entry->getRealValue(fAxes[1]->GetName()),
00472                               entry->getRealValue(fAxes[2]->GetName()),
00473                               fChain->Weight());
00474    }
00475 
00476    if (fDimension >= 1)
00477       fHist->GetXaxis()->SetTitle(fAxes[0]->GetName());
00478    if (fDimension >= 2)
00479       fHist->GetYaxis()->SetTitle(fAxes[1]->GetName());
00480    if (fDimension >= 3)
00481       fHist->GetZaxis()->SetTitle(fAxes[2]->GetName());
00482 }
00483 
00484 void MCMCInterval::CreateSparseHist()
00485 {
00486    if (fAxes == NULL || fChain == NULL) {
00487       coutE(InputArguments) << "* Error in MCMCInterval::CreateSparseHist(): "
00488                             << "Crucial data member was NULL." << endl;
00489       coutE(InputArguments) << "Make sure to fully construct/initialize."
00490                             << endl;
00491       return;
00492    }
00493    if (fSparseHist != NULL)
00494       delete fSparseHist;
00495 
00496    Double_t* min = new Double_t[fDimension];
00497    Double_t* max = new Double_t[fDimension];
00498    Int_t* bins = new Int_t[fDimension];
00499    for (Int_t i = 0; i < fDimension; i++) {
00500       min[i] = fAxes[i]->getMin();
00501       max[i] = fAxes[i]->getMax();
00502       bins[i] = fAxes[i]->numBins();
00503    }
00504    fSparseHist = new THnSparseF("posterior", "MCMC Posterior Histogram",
00505          fDimension, bins, min, max);
00506 
00507    delete[] min;
00508    delete[] max;
00509    delete[] bins;
00510 
00511    // kbelasco: it appears we need to call Sumw2() just to get the
00512    // histogram to keep a running total of the weight so that Getsumw doesn't
00513    // just return 0
00514    fSparseHist->Sumw2();
00515 
00516    if (fNumBurnInSteps >= fChain->Size()) {
00517       coutE(InputArguments) <<
00518          "MCMCInterval::CreateSparseHist: creation of histogram failed: " <<
00519          "Number of burn-in steps (num steps to ignore) >= number of steps " <<
00520          "in Markov chain." << endl;
00521    }
00522 
00523    // Fill histogram
00524    Int_t size = fChain->Size();
00525    const RooArgSet* entry;
00526    Double_t* x = new Double_t[fDimension];
00527    for (Int_t i = fNumBurnInSteps; i < size; i++) {
00528       entry = fChain->Get(i);
00529       for (Int_t ii = 0; ii < fDimension; ii++)
00530          x[ii] = entry->getRealValue(fAxes[ii]->GetName());
00531       fSparseHist->Fill(x, fChain->Weight());
00532    }
00533    delete[] x;
00534 }
00535 
00536 void MCMCInterval::CreateDataHist()
00537 {
00538    if (fParameters.getSize() == 0 || fChain == NULL) {
00539       coutE(Eval) << "* Error in MCMCInterval::CreateDataHist(): " <<
00540                      "Crucial data member was NULL or empty." << endl;
00541       coutE(Eval) << "Make sure to fully construct/initialize." << endl;
00542       return;
00543    }
00544 
00545    if (fNumBurnInSteps >= fChain->Size()) {
00546       coutE(InputArguments) <<
00547          "MCMCInterval::CreateDataHist: creation of histogram failed: " <<
00548          "Number of burn-in steps (num steps to ignore) >= number of steps " <<
00549          "in Markov chain." << endl;
00550       fDataHist = NULL;
00551       return;
00552    }
00553 
00554    fDataHist = fChain->GetAsDataHist(SelectVars(fParameters),
00555          EventRange(fNumBurnInSteps, fChain->Size()));
00556 }
00557 
00558 void MCMCInterval::CreateVector(RooRealVar* param)
00559 {
00560    fVector.clear();
00561    fVecWeight = 0;
00562 
00563    if (fChain == NULL) {
00564       coutE(InputArguments) << "* Error in MCMCInterval::CreateVector(): " <<
00565                      "Crucial data member (Markov chain) was NULL." << endl;
00566       coutE(InputArguments) << "Make sure to fully construct/initialize."
00567                             << endl;
00568       return;
00569    }
00570 
00571    if (fNumBurnInSteps >= fChain->Size()) {
00572       coutE(InputArguments) <<
00573          "MCMCInterval::CreateVector: creation of vector failed: " <<
00574          "Number of burn-in steps (num steps to ignore) >= number of steps " <<
00575          "in Markov chain." << endl;
00576    }
00577 
00578    // Fill vector
00579    Int_t size = fChain->Size() - fNumBurnInSteps;
00580    fVector.resize(size);
00581    Int_t i;
00582    Int_t chainIndex;
00583    for (i = 0; i < size; i++) {
00584       chainIndex = i + fNumBurnInSteps;
00585       fVector[i] = chainIndex;
00586       fVecWeight += fChain->Weight(chainIndex);
00587    }
00588 
00589    stable_sort(fVector.begin(), fVector.end(),
00590                CompareVectorIndices(fChain, param));
00591 }
00592 
00593 void MCMCInterval::SetParameters(const RooArgSet& parameters)
00594 {
00595    fParameters.removeAll();
00596    fParameters.add(parameters);
00597    fDimension = fParameters.getSize();
00598    if (fAxes != NULL)
00599       delete[] fAxes;
00600    fAxes = new RooRealVar*[fDimension];
00601    TIterator* it = fParameters.createIterator();
00602    Int_t n = 0;
00603    TObject* obj;
00604    while ((obj = it->Next()) != NULL) {
00605       if (dynamic_cast<RooRealVar*>(obj) != NULL)
00606          fAxes[n] = (RooRealVar*)obj;
00607       else
00608          coutE(Eval) << "* Error in MCMCInterval::SetParameters: " <<
00609                      obj->GetName() << " not a RooRealVar*" << std::endl;
00610       n++;
00611    }
00612    delete it;
00613 }
00614 
00615 void MCMCInterval::DetermineInterval()
00616 {
00617    switch (fIntervalType) {
00618       case kShortest:
00619          DetermineShortestInterval();
00620          break;
00621       case kTailFraction:
00622          DetermineTailFractionInterval();
00623          break;
00624       default:
00625          coutE(InputArguments) << "MCMCInterval::DetermineInterval(): " <<
00626             "Error: Interval type not set" << endl;
00627          break;
00628    }
00629 }
00630 
00631 void MCMCInterval::DetermineShortestInterval()
00632 {
00633    if (fUseKeys)
00634       DetermineByKeys();
00635    else
00636       DetermineByHist();
00637 }
00638 
00639 void MCMCInterval::DetermineTailFractionInterval()
00640 {
00641    if (fLeftSideTF < 0 || fLeftSideTF > 1) {
00642       coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval: "
00643          << "Fraction must be in the range [0, 1].  "
00644          << fLeftSideTF << "is not allowed." << endl;
00645       return;
00646    }
00647 
00648    if (fDimension != 1) {
00649       coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval(): "
00650          << "Error: Can only find a tail-fraction interval for 1-D intervals"
00651          << endl;
00652       return;
00653    }
00654 
00655    if (fAxes == NULL) {
00656       coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval(): "
00657                             << "Crucial data member was NULL." << endl;
00658       coutE(InputArguments) << "Make sure to fully construct/initialize."
00659                             << endl;
00660       return;
00661    }
00662 
00663    // kbelasco: fill in code here to find interval
00664    //
00665    // also make changes so that calling GetPosterior...() returns NULL
00666    // when fIntervalType == kTailFraction, since there really
00667    // is no posterior for this type of interval determination
00668    if (fVector.size() == 0)
00669       CreateVector(fAxes[0]);
00670 
00671    if (fVector.size() == 0 || fVecWeight == 0) {
00672       // if size is still 0, then creation failed.
00673       // if fVecWeight == 0, then there are no entries (indicates the same
00674       // error as fVector.size() == 0 because that only happens when
00675       // fNumBurnInSteps >= fChain->Size())
00676       // either way, reset and return
00677       fVector.clear();
00678       fTFLower = -1.0 * RooNumber::infinity();
00679       fTFUpper = RooNumber::infinity();
00680       fTFConfLevel = 0.0;
00681       fVecWeight = 0;
00682       return;
00683    }
00684 
00685    RooRealVar* param = fAxes[0];
00686 
00687    Double_t c = fConfidenceLevel;
00688    Double_t leftTailCutoff  = fVecWeight * (1 - c) * fLeftSideTF;
00689    Double_t rightTailCutoff = fVecWeight * (1 - c) * (1 - fLeftSideTF);
00690    Double_t leftTailSum  = 0;
00691    Double_t rightTailSum = 0;
00692 
00693    // kbelasco: consider changing these values to +infinty and -infinity
00694    Double_t ll = param->getMin();
00695    Double_t ul = param->getMax();
00696 
00697    Double_t x;
00698    Double_t w;
00699 
00700    // save a lot of GetName() calls if compiler does not already optimize this
00701    const char* name = param->GetName();
00702 
00703    // find lower limit
00704    Int_t i;
00705    for (i = 0; i < (Int_t)fVector.size(); i++) {
00706       x = fChain->Get(fVector[i])->getRealValue(name);
00707       w = fChain->Weight();
00708 
00709       if (TMath::Abs(leftTailSum + w - leftTailCutoff) <
00710           TMath::Abs(leftTailSum - leftTailCutoff)) {
00711          // moving the lower limit to x would bring us closer to the desired
00712          // left tail size
00713          ll = x;
00714          leftTailSum += w;
00715       } else
00716          break;
00717    }
00718 
00719    // find upper limit
00720    for (i = (Int_t)fVector.size() - 1; i >= 0; i--) {
00721       x = fChain->Get(fVector[i])->getRealValue(name);
00722       w = fChain->Weight();
00723 
00724       if (TMath::Abs(rightTailSum + w - rightTailCutoff) <
00725           TMath::Abs(rightTailSum - rightTailCutoff)) {
00726          // moving the lower limit to x would bring us closer to the desired
00727          // left tail size
00728          ul = x;
00729          rightTailSum += w;
00730       } else
00731          break;
00732    }
00733 
00734    fTFLower = ll;
00735    fTFUpper = ul;
00736    fTFConfLevel = 1 - (leftTailSum + rightTailSum) / fVecWeight;
00737 }
00738 
00739 void MCMCInterval::DetermineByKeys()
00740 {
00741    if (fKeysPdf == NULL)
00742       CreateKeysPdf();
00743 
00744    if (fKeysPdf == NULL) {
00745       // if fKeysPdf is still NULL, then it means CreateKeysPdf failed
00746       // so clear all the data members this function would normally determine
00747       // and return
00748       fFull = 0.0;
00749       fKeysCutoff = -1;
00750       fKeysConfLevel = 0.0;
00751       return;
00752    }
00753 
00754    // now we have a keys pdf of the posterior
00755 
00756    Double_t cutoff = 0.0;
00757    fCutoffVar->setVal(cutoff);
00758    RooAbsReal* integral = fProduct->createIntegral(fParameters, NormSet(fParameters));
00759    Double_t full = integral->getVal(fParameters);
00760    fFull = full;
00761    delete integral;
00762    if (full < 0.98) {
00763       coutW(Eval) << "Warning: Integral of Keys PDF came out to " << full
00764          << " intead of expected value 1.  Will continue using this "
00765          << "factor to normalize further integrals of this PDF." << endl;
00766    }
00767 
00768    // kbelasco: Is there a better way to set the search range?
00769    // from 0 to max value of Keys
00770    // kbelasco: how to get max value?
00771    //Double_t max = product.maxVal(product.getMaxVal(fParameters));
00772 
00773    Double_t volume = 1.0;
00774    TIterator* it = fParameters.createIterator();
00775    RooRealVar* var;
00776    while ((var = (RooRealVar*)it->Next()) != NULL)
00777       volume *= (var->getMax() - var->getMin());
00778    delete it;
00779 
00780    Double_t topCutoff = full / volume;
00781    Double_t bottomCutoff = topCutoff;
00782    Double_t confLevel = CalcConfLevel(topCutoff, full);
00783    if (AcceptableConfLevel(confLevel)) {
00784       fKeysConfLevel = confLevel;
00785       fKeysCutoff = topCutoff;
00786       return;
00787    }
00788    Bool_t changed = kFALSE;
00789    // find high end of range
00790    while (confLevel > fConfidenceLevel) {
00791       topCutoff *= 2.0;
00792       confLevel = CalcConfLevel(topCutoff, full);
00793       if (AcceptableConfLevel(confLevel)) {
00794          fKeysConfLevel = confLevel;
00795          fKeysCutoff = topCutoff;
00796          return;
00797       }
00798       changed = kTRUE;
00799    }
00800    if (changed) {
00801       bottomCutoff = topCutoff / 2.0;
00802    } else {
00803       changed = kFALSE;
00804       bottomCutoff /= 2.0;
00805       confLevel = CalcConfLevel(bottomCutoff, full);
00806       if (AcceptableConfLevel(confLevel)) {
00807          fKeysConfLevel = confLevel;
00808          fKeysCutoff = bottomCutoff;
00809          return;
00810       }
00811       while (confLevel < fConfidenceLevel) {
00812          bottomCutoff /= 2.0;
00813          confLevel = CalcConfLevel(bottomCutoff, full);
00814          if (AcceptableConfLevel(confLevel)) {
00815             fKeysConfLevel = confLevel;
00816             fKeysCutoff = bottomCutoff;
00817             return;
00818          }
00819          changed = kTRUE;
00820       }
00821       if (changed) {
00822          topCutoff = bottomCutoff * 2.0;
00823       }
00824    }
00825 
00826    coutI(Eval) << "range set: [" << bottomCutoff << ", " << topCutoff << "]"
00827                << endl;
00828 
00829    cutoff = (topCutoff + bottomCutoff) / 2.0;
00830    confLevel = CalcConfLevel(cutoff, full);
00831 
00832    // need to use WithinDeltaFraction() because sometimes the integrating the
00833    // posterior in this binary search seems to not have enough granularity to
00834    // find an acceptable conf level (small no. of strange cases).
00835    // WithinDeltaFraction causes the search to terminate when
00836    // topCutoff is essentially equal to bottomCutoff (compared to the magnitude
00837    // of their mean).
00838    while (!AcceptableConfLevel(confLevel) &&
00839           !WithinDeltaFraction(topCutoff, bottomCutoff)) {
00840       if (confLevel > fConfidenceLevel)
00841          bottomCutoff = cutoff;
00842       else if (confLevel < fConfidenceLevel)
00843          topCutoff = cutoff;
00844       cutoff = (topCutoff + bottomCutoff) / 2.0;
00845       coutI(Eval) << "cutoff range: [" << bottomCutoff << ", "
00846                   << topCutoff << "]" << endl;
00847       confLevel = CalcConfLevel(cutoff, full);
00848    }
00849 
00850    fKeysCutoff = cutoff;
00851    fKeysConfLevel = confLevel;
00852 }
00853 
00854 void MCMCInterval::DetermineByHist()
00855 {
00856    if (fUseSparseHist)
00857       DetermineBySparseHist();
00858    else
00859       DetermineByDataHist();
00860 }
00861 
00862 void MCMCInterval::DetermineBySparseHist()
00863 {
00864    Long_t numBins;
00865    if (fSparseHist == NULL)
00866       CreateSparseHist();
00867 
00868    if (fSparseHist == NULL) {
00869       // if fSparseHist is still NULL, then CreateSparseHist failed
00870       fHistCutoff = -1;
00871       fHistConfLevel = 0.0;
00872       return;
00873    }
00874 
00875    numBins = (Long_t)fSparseHist->GetNbins();
00876 
00877    std::vector<Long_t> bins(numBins);
00878    for (Int_t ibin = 0; ibin < numBins; ibin++)
00879       bins[ibin] = (Long_t)ibin;
00880    std::stable_sort(bins.begin(), bins.end(), CompareSparseHistBins(fSparseHist));
00881 
00882    Double_t nEntries = fSparseHist->GetSumw();
00883    Double_t sum = 0;
00884    Double_t content;
00885    Int_t i;
00886    // see above note on indexing to understand numBins - 3
00887    for (i = numBins - 1; i >= 0; i--) {
00888       content = fSparseHist->GetBinContent(bins[i]);
00889       if ((sum + content) / nEntries >= fConfidenceLevel) {
00890          fHistCutoff = content;
00891          if (fIsHistStrict) {
00892             sum += content;
00893             i--;
00894             break;
00895          } else {
00896             i++;
00897             break;
00898          }
00899       }
00900       sum += content;
00901    }
00902 
00903    if (fIsHistStrict) {
00904       // keep going to find the sum
00905       for ( ; i >= 0; i--) {
00906          content = fSparseHist->GetBinContent(bins[i]);
00907          if (content == fHistCutoff)
00908             sum += content;
00909          else
00910             break; // content must be < fHistCutoff
00911       }
00912    } else {
00913       // backtrack to find the cutoff and sum
00914       for ( ; i < numBins; i++) {
00915          content = fSparseHist->GetBinContent(bins[i]);
00916          if (content > fHistCutoff) {
00917             fHistCutoff = content;
00918             break;
00919          } else // content == fHistCutoff
00920             sum -= content;
00921          if (i == numBins - 1)
00922             // still haven't set fHistCutoff correctly yet, and we have no bins
00923             // left, so set fHistCutoff to something higher than the tallest bin
00924             fHistCutoff = content + 1.0;
00925       }
00926    }
00927 
00928    fHistConfLevel = sum / nEntries;
00929 }
00930 
00931 void MCMCInterval::DetermineByDataHist()
00932 {
00933    Int_t numBins;
00934    if (fDataHist == NULL)
00935       CreateDataHist();
00936    if (fDataHist == NULL) {
00937       // if fDataHist is still NULL, then CreateDataHist failed
00938       fHistCutoff = -1;
00939       fHistConfLevel = 0.0;
00940       return;
00941    }
00942 
00943    numBins = fDataHist->numEntries();
00944 
00945    std::vector<Int_t> bins(numBins);
00946    for (Int_t ibin = 0; ibin < numBins; ibin++)
00947       bins[ibin] = ibin;
00948    std::stable_sort(bins.begin(), bins.end(), CompareDataHistBins(fDataHist));
00949 
00950    Double_t nEntries = fDataHist->sum(kFALSE);
00951    Double_t sum = 0;
00952    Double_t content;
00953    Int_t i;
00954    for (i = numBins - 1; i >= 0; i--) {
00955       fDataHist->get(bins[i]);
00956       content = fDataHist->weight();
00957       if ((sum + content) / nEntries >= fConfidenceLevel) {
00958          fHistCutoff = content;
00959          if (fIsHistStrict) {
00960             sum += content;
00961             i--;
00962             break;
00963          } else {
00964             i++;
00965             break;
00966          }
00967       }
00968       sum += content;
00969    }
00970 
00971    if (fIsHistStrict) {
00972       // keep going to find the sum
00973       for ( ; i >= 0; i--) {
00974          fDataHist->get(bins[i]);
00975          content = fDataHist->weight();
00976          if (content == fHistCutoff)
00977             sum += content;
00978          else
00979             break; // content must be < fHistCutoff
00980       }
00981    } else {
00982       // backtrack to find the cutoff and sum
00983       for ( ; i < numBins; i++) {
00984          fDataHist->get(bins[i]);
00985          content = fDataHist->weight();
00986          if (content > fHistCutoff) {
00987             fHistCutoff = content;
00988             break;
00989          } else // content == fHistCutoff
00990             sum -= content;
00991          if (i == numBins - 1)
00992             // still haven't set fHistCutoff correctly yet, and we have no bins
00993             // left, so set fHistCutoff to something higher than the tallest bin
00994             fHistCutoff = content + 1.0;
00995       }
00996    }
00997 
00998    fHistConfLevel = sum / nEntries;
00999 }
01000 
01001 Double_t MCMCInterval::GetActualConfidenceLevel()
01002 {
01003    if (fIntervalType == kShortest) {
01004       if (fUseKeys)
01005          return fKeysConfLevel;
01006       else
01007          return fHistConfLevel;
01008    } else if (fIntervalType == kTailFraction) {
01009       return fTFConfLevel;
01010    } else {
01011       coutE(InputArguments) << "MCMCInterval::GetActualConfidenceLevel: "
01012          << "not implemented for this type of interval.  Returning 0." << endl;
01013       return 0;
01014    }
01015 }
01016 
01017 Double_t MCMCInterval::LowerLimit(RooRealVar& param)
01018 {
01019    switch (fIntervalType) {
01020       case kShortest:
01021          return LowerLimitShortest(param);
01022       case kTailFraction:
01023          return LowerLimitTailFraction(param);
01024       default:
01025          coutE(InputArguments) << "MCMCInterval::LowerLimit(): " <<
01026             "Error: Interval type not set" << endl;
01027          return RooNumber::infinity();
01028    }
01029 }
01030 
01031 Double_t MCMCInterval::UpperLimit(RooRealVar& param)
01032 {
01033    switch (fIntervalType) {
01034       case kShortest:
01035          return UpperLimitShortest(param);
01036       case kTailFraction:
01037          return UpperLimitTailFraction(param);
01038       default:
01039          coutE(InputArguments) << "MCMCInterval::UpperLimit(): " <<
01040             "Error: Interval type not set" << endl;
01041          return RooNumber::infinity();
01042    }
01043 }
01044 
01045 Double_t MCMCInterval::LowerLimitTailFraction(RooRealVar& /*param*/)
01046 {
01047    if (fTFLower == -1.0 * RooNumber::infinity())
01048       DetermineTailFractionInterval();
01049 
01050    return fTFLower;
01051 }
01052 
01053 Double_t MCMCInterval::UpperLimitTailFraction(RooRealVar& /*param*/)
01054 {
01055    if (fTFUpper == RooNumber::infinity())
01056       DetermineTailFractionInterval();
01057 
01058    return fTFUpper;
01059 }
01060 
01061 Double_t MCMCInterval::LowerLimitShortest(RooRealVar& param)
01062 {
01063    if (fUseKeys)
01064       return LowerLimitByKeys(param);
01065    else
01066       return LowerLimitByHist(param);
01067 }
01068 
01069 Double_t MCMCInterval::UpperLimitShortest(RooRealVar& param)
01070 {
01071    if (fUseKeys)
01072       return UpperLimitByKeys(param);
01073    else
01074       return UpperLimitByHist(param);
01075 }
01076 
01077 // Determine the lower limit for param on this interval
01078 // using the binned data set
01079 Double_t MCMCInterval::LowerLimitByHist(RooRealVar& param)
01080 {
01081    if (fUseSparseHist)
01082       return LowerLimitBySparseHist(param);
01083    else
01084       return LowerLimitByDataHist(param);
01085 }
01086 
01087 // Determine the upper limit for param on this interval
01088 // using the binned data set
01089 Double_t MCMCInterval::UpperLimitByHist(RooRealVar& param)
01090 {
01091    if (fUseSparseHist)
01092       return UpperLimitBySparseHist(param);
01093    else
01094       return UpperLimitByDataHist(param);
01095 }
01096 
01097 // Determine the lower limit for param on this interval
01098 // using the binned data set
01099 Double_t MCMCInterval::LowerLimitBySparseHist(RooRealVar& param)
01100 {
01101    if (fDimension != 1) {
01102       coutE(InputArguments) << "In MCMCInterval::LowerLimitBySparseHist: "
01103          << "Sorry, will not compute lower limit unless dimension == 1" << endl;
01104       return param.getMin();
01105    }
01106    if (fHistCutoff < 0)
01107       DetermineBySparseHist(); // this initializes fSparseHist
01108 
01109    if (fHistCutoff < 0) {
01110       // if fHistCutoff < 0 still, then determination of interval failed
01111       coutE(Eval) << "In MCMCInterval::LowerLimitBySparseHist: "
01112          << "couldn't determine cutoff.  Check that num burn in steps < num "
01113          << "steps in the Markov chain.  Returning param.getMin()." << endl;
01114       return param.getMin();
01115    }
01116 
01117    std::vector<Int_t> coord(fDimension);
01118    for (Int_t d = 0; d < fDimension; d++) {
01119       if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
01120          Long_t numBins = (Long_t)fSparseHist->GetNbins();
01121          Double_t lowerLimit = param.getMax();
01122          Double_t val;
01123          for (Long_t i = 0; i < numBins; i++) {
01124             if (fSparseHist->GetBinContent(i, &coord[0]) >= fHistCutoff) {
01125                val = fSparseHist->GetAxis(d)->GetBinCenter(coord[d]);
01126                if (val < lowerLimit)
01127                   lowerLimit = val;
01128             }
01129          }
01130          return lowerLimit;
01131       }
01132    }
01133    return param.getMin();
01134 }
01135 
01136 // Determine the lower limit for param on this interval
01137 // using the binned data set
01138 Double_t MCMCInterval::LowerLimitByDataHist(RooRealVar& param)
01139 {
01140    if (fHistCutoff < 0)
01141       DetermineByDataHist(); // this initializes fDataHist
01142 
01143    if (fHistCutoff < 0) {
01144       // if fHistCutoff < 0 still, then determination of interval failed
01145       coutE(Eval) << "In MCMCInterval::LowerLimitByDataHist: "
01146          << "couldn't determine cutoff.  Check that num burn in steps < num "
01147          << "steps in the Markov chain.  Returning param.getMin()." << endl;
01148       return param.getMin();
01149    }
01150 
01151    for (Int_t d = 0; d < fDimension; d++) {
01152       if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
01153          Int_t numBins = fDataHist->numEntries();
01154          Double_t lowerLimit = param.getMax();
01155          Double_t val;
01156          for (Int_t i = 0; i < numBins; i++) {
01157             fDataHist->get(i);
01158             if (fDataHist->weight() >= fHistCutoff) {
01159                val = fDataHist->get()->getRealValue(param.GetName());
01160                if (val < lowerLimit)
01161                   lowerLimit = val;
01162             }
01163          }
01164          return lowerLimit;
01165       }
01166    }
01167    return param.getMin();
01168 }
01169 
01170 // Determine the upper limit for param on this interval
01171 // using the binned data set
01172 Double_t MCMCInterval::UpperLimitBySparseHist(RooRealVar& param)
01173 {
01174    if (fDimension != 1) {
01175       coutE(InputArguments) << "In MCMCInterval::UpperLimitBySparseHist: "
01176          << "Sorry, will not compute upper limit unless dimension == 1" << endl;
01177       return param.getMax();
01178    }
01179    if (fHistCutoff < 0)
01180       DetermineBySparseHist(); // this initializes fSparseHist
01181 
01182    if (fHistCutoff < 0) {
01183       // if fHistCutoff < 0 still, then determination of interval failed
01184       coutE(Eval) << "In MCMCInterval::UpperLimitBySparseHist: "
01185          << "couldn't determine cutoff.  Check that num burn in steps < num "
01186          << "steps in the Markov chain.  Returning param.getMax()." << endl;
01187       return param.getMax();
01188    }
01189 
01190    std::vector<Int_t> coord(fDimension);
01191    for (Int_t d = 0; d < fDimension; d++) {
01192       if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
01193          Long_t numBins = (Long_t)fSparseHist->GetNbins();
01194          Double_t upperLimit = param.getMin();
01195          Double_t val;
01196          for (Long_t i = 0; i < numBins; i++) {
01197             if (fSparseHist->GetBinContent(i, &coord[0]) >= fHistCutoff) {
01198                val = fSparseHist->GetAxis(d)->GetBinCenter(coord[d]);
01199                if (val > upperLimit)
01200                   upperLimit = val;
01201             }
01202          }
01203          return upperLimit;
01204       }
01205    }
01206    return param.getMax();
01207 }
01208 
01209 // Determine the upper limit for param on this interval
01210 // using the binned data set
01211 Double_t MCMCInterval::UpperLimitByDataHist(RooRealVar& param)
01212 {
01213    if (fHistCutoff < 0)
01214       DetermineByDataHist(); // this initializes fDataHist
01215 
01216    if (fHistCutoff < 0) {
01217       // if fHistCutoff < 0 still, then determination of interval failed
01218       coutE(Eval) << "In MCMCInterval::UpperLimitByDataHist: "
01219          << "couldn't determine cutoff.  Check that num burn in steps < num "
01220          << "steps in the Markov chain.  Returning param.getMax()." << endl;
01221       return param.getMax();
01222    }
01223 
01224    for (Int_t d = 0; d < fDimension; d++) {
01225       if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
01226          Int_t numBins = fDataHist->numEntries();
01227          Double_t upperLimit = param.getMin();
01228          Double_t val;
01229          for (Int_t i = 0; i < numBins; i++) {
01230             fDataHist->get(i);
01231             if (fDataHist->weight() >= fHistCutoff) {
01232                val = fDataHist->get()->getRealValue(param.GetName());
01233                if (val > upperLimit)
01234                   upperLimit = val;
01235             }
01236          }
01237          return upperLimit;
01238       }
01239    }
01240    return param.getMax();
01241 }
01242 
01243 // Determine the lower limit for param on this interval
01244 // using the keys pdf
01245 Double_t MCMCInterval::LowerLimitByKeys(RooRealVar& param)
01246 {
01247    if (fKeysCutoff < 0)
01248       DetermineByKeys();
01249 
01250    if (fKeysDataHist == NULL)
01251       CreateKeysDataHist();
01252 
01253    if (fKeysCutoff < 0 || fKeysDataHist == NULL) {
01254       // failure in determination of cutoff and/or creation of histogram
01255       coutE(Eval) << "in MCMCInterval::LowerLimitByKeys(): "
01256          << "couldn't find lower limit, check that the number of burn in "
01257          << "steps < number of total steps in the Markov chain.  Returning "
01258          << "param.getMin()" << endl;
01259       return param.getMin();
01260    }
01261 
01262    for (Int_t d = 0; d < fDimension; d++) {
01263       if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
01264          Int_t numBins = fKeysDataHist->numEntries();
01265          Double_t lowerLimit = param.getMax();
01266          Double_t val;
01267          for (Int_t i = 0; i < numBins; i++) {
01268             fKeysDataHist->get(i);
01269             if (fKeysDataHist->weight() >= fKeysCutoff) {
01270                val = fKeysDataHist->get()->getRealValue(param.GetName());
01271                if (val < lowerLimit)
01272                   lowerLimit = val;
01273             }
01274          }
01275          return lowerLimit;
01276       }
01277    }
01278    return param.getMin();
01279 }
01280 
01281 // Determine the upper limit for param on this interval
01282 // using the keys pdf
01283 Double_t MCMCInterval::UpperLimitByKeys(RooRealVar& param)
01284 {
01285    if (fKeysCutoff < 0)
01286       DetermineByKeys();
01287 
01288    if (fKeysDataHist == NULL)
01289       CreateKeysDataHist();
01290 
01291    if (fKeysCutoff < 0 || fKeysDataHist == NULL) {
01292       // failure in determination of cutoff and/or creation of histogram
01293       coutE(Eval) << "in MCMCInterval::UpperLimitByKeys(): "
01294          << "couldn't find upper limit, check that the number of burn in "
01295          << "steps < number of total steps in the Markov chain.  Returning "
01296          << "param.getMax()" << endl;
01297       return param.getMax();
01298    }
01299 
01300    for (Int_t d = 0; d < fDimension; d++) {
01301       if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
01302          Int_t numBins = fKeysDataHist->numEntries();
01303          Double_t upperLimit = param.getMin();
01304          Double_t val;
01305          for (Int_t i = 0; i < numBins; i++) {
01306             fKeysDataHist->get(i);
01307             if (fKeysDataHist->weight() >= fKeysCutoff) {
01308                val = fKeysDataHist->get()->getRealValue(param.GetName());
01309                if (val > upperLimit)
01310                   upperLimit = val;
01311             }
01312          }
01313          return upperLimit;
01314       }
01315    }
01316    return param.getMax();
01317 }
01318 
01319 // Determine the approximate maximum value of the Keys PDF
01320 Double_t MCMCInterval::GetKeysMax()
01321 {
01322    if (fKeysCutoff < 0)
01323       DetermineByKeys();
01324 
01325    if (fKeysDataHist == NULL)
01326       CreateKeysDataHist();
01327 
01328    if (fKeysDataHist == NULL) {
01329       // failure in determination of cutoff and/or creation of histogram
01330       coutE(Eval) << "in MCMCInterval::KeysMax(): "
01331          << "couldn't find Keys max value, check that the number of burn in "
01332          << "steps < number of total steps in the Markov chain.  Returning 0"
01333          << endl;
01334       return 0;
01335    }
01336 
01337    Int_t numBins = fKeysDataHist->numEntries();
01338    Double_t max = 0;
01339    Double_t w;
01340    for (Int_t i = 0; i < numBins; i++) {
01341       fKeysDataHist->get(i);
01342       w = fKeysDataHist->weight();
01343       if (w > max)
01344          max = w;
01345    }
01346 
01347    return max;
01348 }
01349 
01350 Double_t MCMCInterval::GetHistCutoff()
01351 {
01352    if (fHistCutoff < 0)
01353       DetermineByHist();
01354 
01355    return fHistCutoff;
01356 }
01357 
01358 Double_t MCMCInterval::GetKeysPdfCutoff()
01359 {
01360    if (fKeysCutoff < 0)
01361       DetermineByKeys();
01362 
01363    // kbelasco: if fFull hasn't been set (because Keys creation failed because
01364    // fNumBurnInSteps >= fChain->Size()) then this will return infinity, which
01365    // seems ok to me since it will indicate error
01366 
01367    return fKeysCutoff / fFull;
01368 }
01369 
01370 Double_t MCMCInterval::CalcConfLevel(Double_t cutoff, Double_t full)
01371 {
01372    RooAbsReal* integral;
01373    Double_t confLevel;
01374    fCutoffVar->setVal(cutoff);
01375    integral = fProduct->createIntegral(fParameters, NormSet(fParameters));
01376    confLevel = integral->getVal(fParameters) / full;
01377    coutI(Eval) << "cutoff = " << cutoff << ", conf = " << confLevel << endl;
01378    //cout << "tmp: cutoff = " << cutoff << ", conf = " << confLevel << endl;
01379    delete integral;
01380    return confLevel;
01381 }
01382 
01383 TH1* MCMCInterval::GetPosteriorHist()
01384 {
01385   if(fConfidenceLevel == 0)
01386      coutE(InputArguments) << "Error in MCMCInterval::GetPosteriorHist: "
01387                            << "confidence level not set " << endl;
01388   if (fHist == NULL)
01389      CreateHist();
01390 
01391   if (fHist == NULL)
01392      // if fHist is still NULL, then CreateHist failed
01393      return NULL;
01394 
01395   return (TH1*) fHist->Clone("MCMCposterior_hist");
01396 }
01397 
01398 RooNDKeysPdf* MCMCInterval::GetPosteriorKeysPdf()
01399 {
01400    if (fConfidenceLevel == 0)
01401       coutE(InputArguments) << "Error in MCMCInterval::GetPosteriorKeysPdf: "
01402                             << "confidence level not set " << endl;
01403    if (fKeysPdf == NULL)
01404       CreateKeysPdf();
01405 
01406    if (fKeysPdf == NULL)
01407       // if fKeysPdf is still NULL, then it means CreateKeysPdf failed
01408       return NULL;
01409 
01410    return (RooNDKeysPdf*) fKeysPdf->Clone("MCMCPosterior_keys");
01411 }
01412 
01413 RooProduct* MCMCInterval::GetPosteriorKeysProduct()
01414 {
01415    if (fConfidenceLevel == 0)
01416       coutE(InputArguments) << "MCMCInterval::GetPosteriorKeysProduct: "
01417                             << "confidence level not set " << endl;
01418    if (fProduct == NULL) {
01419       CreateKeysPdf();
01420       DetermineByKeys();
01421    }
01422 
01423    if (fProduct == NULL)
01424       // if fProduct is still NULL, then it means CreateKeysPdf failed
01425       return NULL;
01426 
01427    return (RooProduct*) fProduct->Clone("MCMCPosterior_keysproduct");
01428 }
01429 
01430 RooArgSet* MCMCInterval::GetParameters() const
01431 {  
01432    // returns list of parameters
01433    return new RooArgSet(fParameters);
01434 }
01435 
01436 Bool_t MCMCInterval::AcceptableConfLevel(Double_t confLevel)
01437 {
01438    return (TMath::Abs(confLevel - fConfidenceLevel) < fEpsilon);
01439 }
01440 
01441 Bool_t MCMCInterval::WithinDeltaFraction(Double_t a, Double_t b)
01442 {
01443    return (TMath::Abs(a - b) < TMath::Abs(fDelta * (a + b)/2));
01444 }
01445 
01446 void MCMCInterval::CreateKeysDataHist()
01447 {
01448    if (fAxes == NULL)
01449       return;
01450    if (fProduct == NULL)
01451       DetermineByKeys();
01452    if (fProduct == NULL)
01453       // if fProduct still NULL, then creation failed
01454       return;
01455 
01456    //RooAbsBinning** savedBinning = new RooAbsBinning*[fDimension];
01457    Int_t* savedBins = new Int_t[fDimension];
01458    Int_t i;
01459    Double_t numBins;
01460    RooRealVar* var;
01461 
01462    // kbelasco: Note - the accuracy is only increased here if the binning for
01463    // each RooRealVar is uniform
01464 
01465    // kbelasco: look into why saving the binnings and replacing them doesn't
01466    // work (replaces with 1 bin always).
01467    // Note: this code modifies the binning for the parameters (if they are
01468    // uniform) and sets them back to what they were.  If the binnings are not
01469    // uniform, this code does nothing.
01470 
01471    // first scan through fAxes to make sure all binnings are uniform, or else
01472    // we can't change the number of bins because there seems to be an error
01473    // when setting the binning itself rather than just the number of bins
01474    Bool_t tempChangeBinning = true;
01475    for (i = 0; i < fDimension; i++) {
01476       if (!fAxes[i]->getBinning(NULL, false, false).isUniform()) {
01477          tempChangeBinning = false;
01478          break;
01479       }
01480    }
01481 
01482    // kbelasco: for 1 dimension this should be fine, but for more dimensions
01483    // the total nubmer of bins in the histogram increases exponentially with
01484    // the dimension, so don't do this above 1-D for now.
01485    if (fDimension >= 2)
01486       tempChangeBinning = false;
01487 
01488    if (tempChangeBinning) {
01489       // set high nubmer of bins for high accuracy on lower/upper limit by keys
01490       for (i = 0; i < fDimension; i++) {
01491          var = fAxes[i];
01492          //savedBinning[i] = &var->getBinning("__binning_clone", false, true);
01493          savedBins[i] = var->getBinning(NULL, false, false).numBins();
01494          numBins = (var->getMax() - var->getMin()) / fEpsilon;
01495          var->setBins((Int_t)numBins);
01496       }
01497    }
01498 
01499    fKeysDataHist = new RooDataHist("_productDataHist",
01500          "Keys PDF & Heaviside Product Data Hist", fParameters);
01501    fKeysDataHist = fProduct->fillDataHist(fKeysDataHist, &fParameters, 1.);
01502 
01503    if (tempChangeBinning) {
01504       // set the binning back to normal
01505       for (i = 0; i < fDimension; i++)
01506          //fAxes[i]->setBinning(*savedBinning[i], NULL);
01507          //fAxes[i]->setBins(savedBinning[i]->numBins(), NULL);
01508          fAxes[i]->setBins(savedBins[i], NULL);
01509    }
01510 
01511    //delete[] savedBinning;
01512    delete[] savedBins;
01513 }
01514 
01515 Bool_t MCMCInterval::CheckParameters(const RooArgSet& parameterPoint) const
01516 {  
01517    // check that the parameters are correct
01518 
01519    if (parameterPoint.getSize() != fParameters.getSize() ) {
01520      coutE(Eval) << "MCMCInterval: size is wrong, parameters don't match" << std::endl;
01521      return kFALSE;
01522    }
01523    if ( ! parameterPoint.equals( fParameters ) ) {
01524      coutE(Eval) << "MCMCInterval: size is ok, but parameters don't match" << std::endl;
01525      return kFALSE;
01526    }
01527    return kTRUE;
01528 }

Generated on Tue Jul 5 15:14:35 2011 for ROOT_528-00b_version by  doxygen 1.5.1