testUnfold3.C

Go to the documentation of this file.
00001 // Test program for the class TUnfoldSys
00002 // Author: Stefan Schmitt
00003 // DESY, 14.10.2008
00004 
00005 //  Version 16, parallel to changes in TUnfold
00006 //
00007 //  History:
00008 //     Version 15, simple example including background subtraction
00009 
00010 #include <TMath.h>
00011 #include <TCanvas.h>
00012 #include <TRandom3.h>
00013 #include <TFitter.h>
00014 #include <TF1.h>
00015 #include <TStyle.h>
00016 #include <TVector.h>
00017 #include <TGraph.h>
00018 
00019 #include <TUnfoldSys.h>
00020 
00021 using namespace std;
00022 
00023 ///////////////////////////////////////////////////////////////////////
00024 // 
00025 // Test program for the class TUnfoldSys
00026 //
00027 //  the goal is to unfold the underlying "true" distribution of a variable Pt
00028 //
00029 //  the reconstructed Pt is measured in 24 bins from 4 to 28 
00030 //  the generator-level Pt is unfolded into 10 bins from 6 to 26
00031 //    plus underflow bin from 0 to 6
00032 //    plus overflow bin above 26 
00033 //  there are two background sources
00034 //       bgr1 and bgr2
00035 //  the signal has a finite trigger efficiency at a threshold of 8 GeV
00036 //
00037 //  two types of systematic error are studied:
00038 //    SYS1: variation of the input shape
00039 //    SYS2: variation of the trigger efficiency at threshold
00040 //
00041 //  Finally, the unfolding is compared to a "bin-by-bin" correction method
00042 //
00043 ///////////////////////////////////////////////////////////////////////
00044 
00045 TRandom *rnd=0;
00046 
00047 Double_t GenerateEvent(const Double_t *parm,
00048                        const Double_t *triggerParm,
00049                        Double_t *intLumi,
00050                        Bool_t *triggerFlag,
00051                        Double_t *ptGen,Int_t *iType)
00052 {
00053    // generate an event
00054    // input:
00055    //      parameters for the event generator
00056    // return value:
00057    //      reconstructed Pt
00058    // output to pointers:
00059    //      integrated luminosity
00060    //      several variables only accessible on generator level
00061    //
00062    // the parm array defines the physical parameters
00063    //  parm[0]: background source 1 fraction
00064    //  parm[1]: background source 2 fraction
00065    //  parm[2]: lower limit of generated Pt distribution
00066    //  parm[3]: upper limit of generated Pt distribution
00067    //  parm[4]: exponent for Pt distribution signal
00068    //  parm[5]: exponent for Pt distribution background 1
00069    //  parm[6]: exponent for Pt distribution background 2
00070    //  parm[7]: resolution parameter a goes with sqrt(Pt)
00071    //  parm[8]: resolution parameter b goes with Pt
00072    //  triggerParm[0]: trigger threshold turn-on position
00073    //  triggerParm[1]: trigger threshold turn-on width
00074    //  triggerParm[2]: trigger efficiency for high pt
00075    //
00076    // intLumi is advanced bu 1 for each *generated* event
00077    // for data, several events may be generated, until one passes the trigger
00078    //
00079    // some generator-level quantities are also returned:
00080    //   triggerFlag: whether the event passed the trigger threshold
00081    //   ptGen: the generated pt
00082    //   iType: which type of process was simulated
00083    //
00084    // the "triggered" also has another meaning:
00085    //   if(triggerFlag==0)   only triggered events are returned
00086    //
00087    // Usage to generate data events
00088    //      ptObs=GenerateEvent(parm,triggerParm,0,0,0)
00089    //
00090    // Usage to generate MC events
00091    //      ptGen=GenerateEvent(parm,triggerParm,&triggerFlag,&ptGen,&iType);
00092    //
00093    Double_t ptObs;
00094    Bool_t isTriggered=kFALSE;
00095    do {
00096       Int_t itype;
00097       Double_t ptgen;
00098       Double_t f=rnd->Rndm();
00099       // decide whether this is background or signal
00100       itype=0;
00101       if(f<parm[0]) itype=1;
00102       else if(f<parm[0]+parm[1]) itype=2;
00103       // generate Pt according to distribution  pow(ptgen,a)
00104       // get exponent
00105       Double_t a=parm[4+itype];
00106       Double_t a1=a+1.0;
00107       Double_t t=rnd->Rndm();
00108       if(a1 == 0.0) {
00109          Double_t x0=TMath::Log(parm[2]);
00110          ptgen=TMath::Exp(t*(TMath::Log(parm[3])-x0)+x0);
00111       } else {
00112          Double_t x0=pow(parm[2],a1);
00113          ptgen=pow(t*(pow(parm[3],a1)-x0)+x0,1./a1);
00114       }
00115       if(iType) *iType=itype;
00116       if(ptGen) *ptGen=ptgen;
00117       
00118       // smearing in Pt with large asymmetric tail
00119       Double_t sigma=
00120          TMath::Sqrt(parm[7]*parm[7]*ptgen+parm[8]*parm[8]*ptgen*ptgen);
00121       ptObs=rnd->BreitWigner(ptgen,sigma);
00122 
00123       // decide whether event was triggered
00124       Double_t triggerProb =
00125          triggerParm[2]/(1.+TMath::Exp((triggerParm[0]-ptObs)/triggerParm[1]));
00126       isTriggered= rnd->Rndm()<triggerProb;
00127       (*intLumi) ++;
00128    } while((!triggerFlag) && (!isTriggered));
00129    // return trigger decision
00130    if(triggerFlag) *triggerFlag=isTriggered;
00131    return ptObs;
00132 }
00133 
00134 //int main(int argc, char *argv[])
00135 void testUnfold3()
00136 {
00137   // switch on histogram errors
00138   TH1::SetDefaultSumw2();
00139 
00140   // show fit result
00141   gStyle->SetOptFit(1111);
00142 
00143   // random generator
00144   rnd=new TRandom3();
00145 
00146   // data and MC luminosities
00147   Double_t const lumiData= 30000;
00148   Double_t const lumiMC  =1000000;
00149 
00150   // Binning
00151   // reconstructed mass
00152   Int_t const nDet=24;
00153   Double_t const xminDet=4.0;
00154   Double_t const xmaxDet=28.0;
00155   // signal contributions
00156   Int_t const nGen=10;
00157   Double_t const xminGen= 6.0;
00158   Double_t const xmaxGen=26.0;
00159 
00160   //========================
00161   // Step 1: fill histograms
00162 
00163   // =================== data true parameters ============================
00164   //
00165   // in real life these are unknown
00166   //
00167   static Double_t parm_data[]={
00168      0.05, // fraction of background 1 (on generator level)
00169      0.05, // fraction of background 2 (on generator level)
00170      3.5, // lower Pt cut (generator level)
00171      100.,// upper Pt cut (generator level)
00172      -3.0,// signal exponent
00173      0.1, // background 1 exponent
00174      -0.5, // background 2 exponent
00175      0.2, // energy resolution a term
00176      0.01, // energy resolution b term
00177   };
00178   static Double_t triggerParm_data[]={8.0, // threshold is 8 GeV
00179                                4.0, // width is 4 GeV
00180                                0.95 // high Pt efficiency os 95%
00181   };
00182 
00183   //============================================
00184   // generate data distribution
00185   TH1D *histDetDATA=new TH1D("PT(data)",";Pt(obs)",nDet,xminDet,xmaxDet);
00186   // second data distribution, with generator-level bins
00187   // for bin-by-bin correction study
00188   TH1D *histDetDATAbbb=new TH1D("PT(data,bbb)",";Pt(obs,bbb)",
00189                                 nGen,xminGen,xmaxGen);
00190   // in real life we do not have this distribution...
00191   TH1D *histGenDATA=new TH1D("PT(MC,signal,gen)",";Pt(gen)",
00192                              nGen,xminGen,xmaxGen);
00193 
00194   Double_t intLumi=0.0;
00195   while(intLumi<lumiData) {
00196      Int_t iTypeGen;
00197      Bool_t isTriggered;
00198      Double_t ptGen;
00199      Double_t ptObs=GenerateEvent(parm_data,triggerParm_data,&intLumi,
00200                                   &isTriggered,&ptGen,&iTypeGen);
00201      if(isTriggered) {
00202         histDetDATA->Fill(ptObs);
00203         histDetDATAbbb->Fill(ptObs);
00204      }
00205      // fill generator-level signal histogram for data
00206      // for real data, this does not exist
00207      if(iTypeGen==0) histGenDATA->Fill(ptGen);
00208   }
00209 
00210 
00211   // =================== MC parameters ============================
00212   // default settings
00213   static Double_t parm_MC[]={
00214      0.05, // fraction of background 1 (on generator level)
00215      0.05, // fraction of background 2 (on generator level)
00216      3.5, // lower Pt cut (generator level)
00217      100.,// upper Pt cut (generator level)
00218      -4.0,// signal exponent !!! UNKNOWN !!! steeper than in data
00219           //                                 to illustrate bin-by-bin bias
00220      0.1, // background 1 exponent
00221      -0.5, // background 2 exponent
00222      0.2, // energy resolution a term
00223      0.01, // energy resolution b term
00224   };
00225   static Double_t triggerParm_MC[]={8.0, // threshold is 8 GeV
00226                              4.0, // width is 4 GeV
00227                              0.95 // high Pt efficiency os 95%
00228   };
00229 
00230   
00231   // study trigger systematics: change parameters for trigger threshold
00232   static Double_t triggerParm_MCSYS1[]={7.7, // threshold is 7 GeV
00233                                         3.7, // width is 3 GeV
00234                                         0.95 // high Pt efficiency is 9%
00235   };
00236   // study dependence on initial signal shape
00237   static Double_t parm_MC_SYS2[]={
00238      0.0, // fraction of background: not needed
00239      0.0, // fraction of background: not needed
00240      3.5, // lower Pt cut (generator level)
00241      100.,// upper Pt cut (generator level)
00242      -2.0, // signal exponent changed
00243      0.1, // background 1 exponent
00244      -0.5, // background 2 exponent
00245      0.2, // energy resolution a term
00246      0.01, // energy resolution b term
00247   };
00248 
00249   //============================================
00250   // generate MC distributions
00251   // detector-level histograms
00252   TH1D *histDetMC[3];
00253   histDetMC[0]=new TH1D("PT(MC,signal)",";Pt(obs)",nDet,xminDet,xmaxDet);
00254   histDetMC[1]=new TH1D("PT(MC,bgr1)",";Pt(obs)",nDet,xminDet,xmaxDet);
00255   histDetMC[2]=new TH1D("PT(MC,bgr2)",";Pt(obs)",nDet,xminDet,xmaxDet);
00256   TH1D *histDetMCall=new TH1D("PT(MC,all)",";Pt(obs)",nDet,xminDet,xmaxDet);
00257   TH1D *histDetMCbgr=new TH1D("PT(MC,bgr)",";Pt(obs)",nDet,xminDet,xmaxDet);
00258   // another set of detector level histograms, this time with
00259   // generator-level binning, to try the bin-by-bin correction
00260   TH1D *histDetMCbbb[3];
00261   histDetMCbbb[0]=new TH1D("PT(MC,sig)",";Pt(obs,bbb)",nGen,xminGen,xmaxGen);
00262   histDetMCbbb[1]=new TH1D("PT(MC,bgr1)bbb",";Pt(obs,bbb)",nGen,xminGen,xmaxGen);
00263   histDetMCbbb[2]=new TH1D("PT(MC,bgr2)bbb",";Pt(obs,bbb)",nGen,xminGen,xmaxGen);
00264   // true signal distribution, two different binnings
00265   TH1D *histGenMC=new TH1D("PT(MC,signal,gen)bbb",";Pt(gen)",
00266                          nGen,xminGen,xmaxGen);
00267   TH2D *histGenDet=new TH2D("PT(MC,signal,gen,det)",";Pt(gen);Pt(det)",
00268                             nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00269 
00270   // weighting factor to accomodate for the different luminosity in data and MC
00271   Double_t lumiWeight = lumiData/lumiMC;
00272   intLumi=0.0;
00273   while(intLumi<lumiMC) {
00274      Int_t iTypeGen;
00275      Bool_t isTriggered;
00276      Double_t ptGen;
00277      Double_t ptObs=GenerateEvent(parm_MC,triggerParm_MC,&intLumi,&isTriggered,
00278                                   &ptGen,&iTypeGen);
00279      // detector-level distributions
00280      // require trigger
00281      if(isTriggered) {
00282         // fill signal, bgr1,bgr2 into separate histograms
00283         histDetMC[iTypeGen]->Fill(ptObs,lumiWeight);
00284         histDetMCbbb[iTypeGen]->Fill(ptObs,lumiWeight);
00285         // fill total MC prediction (signal+bgr1+bgr2)
00286         histDetMCall->Fill(ptObs,lumiWeight);
00287         // fill bgr only MC prediction (bgr1+bgr2)
00288         if(iTypeGen>0) {
00289            histDetMCbgr->Fill(ptObs,lumiWeight);
00290         }
00291      }
00292      // generator level distributions (signal only)
00293      if(iTypeGen==0) {
00294         histGenMC->Fill(ptGen,lumiWeight);
00295      }
00296      // matrix dectribing the migrations (signal only)
00297      if(iTypeGen==0) {
00298         // if event was triggered, fill histogram with ptObs
00299         if(isTriggered) {
00300            histGenDet->Fill(ptGen,ptObs,lumiWeight);
00301         } else {
00302            // not triggered: fill detector-level underflow bin
00303            histGenDet->Fill(ptGen,-999.,lumiWeight);
00304         }
00305      }
00306   }
00307 
00308   // generate another MC, this time with changed trigger settings
00309   // fill 2D histogram and histograms for bin-by-bin study
00310   TH2D *histGenDetSYS1=new TH2D("PT(MC,signal,gen,det,sys1)",
00311                                 ";Pt(gen);Pt(obs)",
00312                                 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00313   TH1D *histGenSYS1=new TH1D("PT(MC,signal,gen,sys1)",
00314                              ";Pt(obs)",
00315                              nGen,xminGen,xmaxGen);
00316   TH1D *histDetSYS1bbb=new TH1D("PT(MC,signal,det,sys1,bbb)",
00317                                 ";Pt(obs)",
00318                                 nGen,xminGen,xmaxGen);
00319   intLumi=0.;
00320   while(intLumi<lumiMC) {
00321      Int_t iTypeGen;
00322      Bool_t isTriggered;
00323      Double_t ptGen;
00324      Double_t ptObs=GenerateEvent(parm_MC,triggerParm_MCSYS1,&intLumi,
00325                                   &isTriggered,&ptGen,&iTypeGen);
00326      // matrix describing the migrations (signal only)
00327      if(iTypeGen==0) {
00328         // if event was triggered, fill histogram with ptObs
00329         if(isTriggered) {
00330            histGenDetSYS1->Fill(ptGen,ptObs,lumiWeight);
00331         } else {
00332            // not triggered: fill detector-level underflow bin
00333            histGenDetSYS1->Fill(ptGen,-999.,lumiWeight);
00334         }
00335      }
00336      // generator level distribution for bin-by-bin study
00337      if(iTypeGen==0) {
00338         histGenSYS1->Fill(ptGen,lumiWeight);
00339         // detector level for bin-by-bin study
00340         if(isTriggered) {
00341            histDetSYS1bbb->Fill(ptObs,lumiWeight);
00342         }
00343      }
00344   }
00345   // generate a third MC, this time with changed initial distribution
00346   // only the 2D histogram is filled
00347   TH2D *histGenDetSYS2=new TH2D("PT(MC,signal,gen,det,sys2)",
00348                                 ";Pt(gen);Pt(det)",
00349                                 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00350   TH1D *histGenSYS2=new TH1D("PT(MC,signal,gen,sys2)",
00351                              ";Pt(obs)",
00352                              nGen,xminGen,xmaxGen);
00353   TH1D *histDetSYS2bbb=new TH1D("PT(MC,signal,det,sys2,bbb)",
00354                                 ";Pt(obs)",
00355                                 nGen,xminGen,xmaxGen);
00356   intLumi=0.0;
00357   while(intLumi<lumiMC) {
00358      Int_t iTypeGen;
00359      Bool_t isTriggered;
00360      Double_t ptGen;
00361      Double_t ptObs=GenerateEvent(parm_MC_SYS2,triggerParm_MC,&intLumi,
00362                                   &isTriggered,&ptGen,&iTypeGen);
00363      // matrix dectribing the migrations (signal only)
00364      if(iTypeGen==0) {
00365         // if event was triggered, fill histogram with ptObs
00366         if(isTriggered) {
00367            histGenDetSYS2->Fill(ptGen,ptObs,lumiWeight);
00368         } else {
00369            // not triggered: fill detector-level underflow bin
00370            histGenDetSYS2->Fill(ptGen,-999.,lumiWeight);
00371         }
00372      }
00373      if(iTypeGen==0) {
00374         // generator level distribution for bin-by-bin study
00375         histGenSYS2->Fill(ptGen,lumiWeight);
00376         // detector level for bin-by-bin study
00377         if(isTriggered) {
00378            histDetSYS2bbb->Fill(ptObs,lumiWeight);
00379         }
00380      }
00381   }
00382 
00383   //========================
00384   // Step 2: unfolding
00385   //
00386   // this is based on a matrix (TH2D) which describes the connection
00387   // between
00388   //   nDet    Detector bins
00389   //   nGen+2  Generator level bins
00390   //
00391   // why are there nGen+2 Generator level bins?
00392   //   the two extra bins are the underflow and overflow bin
00393   //   These are unfolded as well. Such bins are needed to
00394   //   accomodate for migrations from outside the phase-space into the
00395   //   observed phase-space.
00396   //
00397   // In addition there are overflow/underflow bins
00398   //   for the reconstructed variable. These are used to count events
00399   //   which are -NOT- reconstructed. Either because they migrate
00400   //   out of the reconstructed phase-space. Or because there are detector
00401   //   inefficiencies.
00402   //
00403   TUnfoldSys unfold(histGenDet,TUnfold::kHistMapOutputHoriz,
00404                     TUnfold::kRegModeSize,TUnfold::kEConstraintNone);
00405 
00406   // define the input vector (the measured data distribution)
00407   unfold.SetInput(histDetDATA);
00408 
00409   // subtract background, normalized to data luminosity
00410   //  with 10% scale error each
00411   Double_t scale_bgr=1.0;
00412   Double_t dscale_bgr=0.2;
00413   unfold.SubtractBackground(histDetMC[1],"background1",scale_bgr,dscale_bgr);
00414   unfold.SubtractBackground(histDetMC[2],"background2",scale_bgr,dscale_bgr);
00415 
00416   // add systematic errors
00417   // trigger threshold
00418   unfold.AddSysError(histGenDetSYS1,"trigger_SYS1",
00419                      TUnfold::kHistMapOutputHoriz,
00420                      TUnfoldSys::kSysErrModeMatrix);
00421   // calorimeter response
00422   unfold.AddSysError(histGenDetSYS2,"signalshape_SYS2",
00423                      TUnfold::kHistMapOutputHoriz,
00424                      TUnfoldSys::kSysErrModeMatrix);
00425 
00426   // run the unfolding
00427   Int_t nScan=30;
00428   TSpline *logTauX,*logTauY;
00429   TGraph *lCurve;
00430   // this method scans the parameter tau and finds the kink in the L curve
00431   // finally, the unfolding is done for the best choice of tau
00432   Int_t iBest=unfold.ScanLcurve(nScan,0.,0.,&lCurve,&logTauX,&logTauY);
00433 
00434   // save graphs with one point to visualize best choice of tau
00435   Double_t t[1],x[1],y[1];
00436   logTauX->GetKnot(iBest,t[0],x[0]);
00437   logTauY->GetKnot(iBest,t[0],y[0]);
00438   TGraph *bestLcurve=new TGraph(1,x,y);
00439   TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
00440 
00441   // get unfolding output
00442   // Reminder: there are
00443   //   nGen+2  Generator level bins
00444   // But we are not interested in the additional bins.
00445   // The mapping is required to get the proper errors and correlations
00446 
00447   // Here a mapping is defined from the nGen+2 bins to nGen bins.
00448   Int_t *binMap=new Int_t[nGen+2];
00449   binMap[0] = -1;    // underflow bin should be ignored
00450   binMap[nGen+1]=-1; // overflow bin should be ignored
00451   for(Int_t i=1;i<=nGen;i++) {
00452      binMap[i]=i;  // map inner bins to output histogram bins
00453   }
00454 
00455   // this returns the unfolding output
00456   // The errors included at this point are:
00457   //    * statistical errros
00458   //    * background subtraction errors
00459   TH1D *histUnfoldStatBgr=new TH1D("PT(unfold,signal,stat+bgr)",";Pt(gen)",
00460                             nGen,xminGen,xmaxGen);
00461   unfold.GetOutput(histUnfoldStatBgr,binMap);
00462 
00463   // retreive error matrix of statistical errors only
00464   TH2D *histEmatStat=new TH2D("ErrMat(stat)",";Pt(gen);Pt(gen)",
00465                               nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00466   unfold.GetEmatrixInput(histEmatStat,binMap);
00467 
00468   // retreive full error matrix
00469   // This includes also all systematic errors defined above
00470   TH2D *histEmatTotal=new TH2D("ErrMat(total)",";Pt(gen);Pt(gen)",
00471                                nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00472   unfold.GetEmatrixTotal(histEmatTotal,binMap);
00473 
00474   // create two copies of the unfolded data, one with statistical errors
00475   // the other with total errors
00476   TH1D *histUnfoldStat=new TH1D("PT(unfold,signal,stat)",";Pt(gen)",
00477                                 nGen,xminGen,xmaxGen);
00478   TH1D *histUnfoldTotal=new TH1D("PT(unfold,signal,total)",";Pt(gen)",
00479                                  nGen,xminGen,xmaxGen);
00480   for(Int_t i=0;i<nGen+2;i++) {
00481      Double_t c=histUnfoldStatBgr->GetBinContent(i);
00482 
00483      // histogram with unfolded data and stat errors
00484      histUnfoldStat->SetBinContent(i,c);
00485      histUnfoldStat->SetBinError
00486         (i,TMath::Sqrt(histEmatStat->GetBinContent(i,i)));
00487 
00488      // histogram with unfolded data and total errors
00489      histUnfoldTotal->SetBinContent(i,c); 
00490      histUnfoldTotal->SetBinError
00491         (i,TMath::Sqrt(histEmatTotal->GetBinContent(i,i)));
00492   }
00493 
00494   // create histogram with correlation matrix
00495   TH2D *histCorr=new TH2D("Corr(total)",";Pt(gen);Pt(gen)",
00496                           nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
00497   for(Int_t i=0;i<nGen+2;i++) {
00498      Double_t ei,ej;
00499      ei=TMath::Sqrt(histEmatTotal->GetBinContent(i,i));
00500      if(ei<=0.0) continue;
00501      for(Int_t j=0;j<nGen+2;j++) {
00502         ej=TMath::Sqrt(histEmatTotal->GetBinContent(j,j));
00503         if(ej<=0.0) continue;
00504         histCorr->SetBinContent(i,j,histEmatTotal->GetBinContent(i,j)/ei/ej);
00505      }
00506   }
00507 
00508   delete [] binMap;
00509 
00510   //========================
00511   // Step 3: plots
00512 
00513   TCanvas output;
00514   output.Divide(3,2);
00515   output.cd(1);
00516   histDetDATA->SetMinimum(0.0);
00517   histDetDATA->Draw("E");
00518   histDetMCall->SetMinimum(0.0);
00519   histDetMCall->SetLineColor(kBlue);  
00520   histDetMCbgr->SetLineColor(kRed);  
00521   histDetMC[1]->SetLineColor(kCyan);  
00522   histDetMCall->Draw("SAME HIST");
00523   histDetMCbgr->Draw("SAME HIST");
00524   histDetMC[1]->Draw("SAME HIST");
00525 
00526   output.cd(2);
00527   histUnfoldTotal->SetMinimum(0.0);
00528   histUnfoldTotal->SetMaximum(histUnfoldTotal->GetMaximum()*1.5);
00529   // outer error: total error
00530   histUnfoldTotal->Draw("E");
00531   // middle error: stat+bgr
00532   histUnfoldStatBgr->Draw("SAME E1");
00533   // inner error: stat only
00534   histUnfoldStat->Draw("SAME E1");
00535   histGenDATA->Draw("SAME HIST");
00536   histGenMC->SetLineColor(kBlue);
00537   histGenMC->Draw("SAME HIST");
00538 
00539   output.cd(3);
00540   histGenDet->SetLineColor(kBlue);
00541   histGenDet->Draw("BOX");
00542 
00543   // show tau as a function of chi**2
00544   output.cd(4);
00545   logTauX->Draw();
00546   bestLogTauLogChi2->SetMarkerColor(kRed);
00547   bestLogTauLogChi2->Draw("*");
00548 
00549   // show the L curve
00550   output.cd(5);
00551   lCurve->Draw("AL");
00552   bestLcurve->SetMarkerColor(kRed);
00553   bestLcurve->Draw("*");
00554 
00555   // show correlation matrix
00556   output.cd(6);
00557   histCorr->Draw("BOX");
00558 
00559   output.SaveAs("testUnfold3.ps");
00560 
00561 
00562   // step 4: compare results to
00563   // the so-called bin-by-bin "correction"
00564   for(Int_t i=1;i<=nGen;i++) {
00565      // data contribution in this bin
00566      Double_t data=histDetDATAbbb->GetBinContent(i);
00567      Double_t errData=histDetDATAbbb->GetBinError(i);
00568 
00569      // subtract background contribution
00570      Double_t data_bgr=data;
00571      Double_t errData_bgr=errData;
00572      for(Int_t j=1;j<=2;j++) {
00573         data_bgr -= scale_bgr*histDetMCbbb[j]->GetBinContent(i);
00574         errData_bgr = TMath::Sqrt(errData_bgr*errData_bgr+
00575                                   scale_bgr*histDetMCbbb[j]->GetBinError(i)*
00576                                   scale_bgr*histDetMCbbb[j]->GetBinError(i)+
00577                                   dscale_bgr*histDetMCbbb[j]->GetBinContent(i)*
00578                                   dscale_bgr*histDetMCbbb[j]->GetBinContent(i)
00579                                   );
00580      }
00581      // "correct" the data, using the Monte Carlo and neglecting off-diagonals
00582      Double_t fCorr=(histGenMC->GetBinContent(i)/
00583                       histDetMCbbb[0]->GetBinContent(i));
00584      Double_t data_bbb= data_bgr *fCorr;
00585      // stat only error
00586      Double_t errData_stat_bbb = errData*fCorr;
00587      // stat plus background subtraction
00588      Double_t errData_statbgr_bbb = errData_bgr*fCorr;
00589      // estimate systematic error by repeating the exercise
00590      // using the MC with systematic shifts applied
00591      Double_t fCorr_1=(histGenSYS1->GetBinContent(i)/
00592                        histDetSYS1bbb->GetBinContent(i));
00593      Double_t shift_sys1= data_bgr*fCorr_1 - data_bbb;
00594      Double_t fCorr_2=(histGenSYS2->GetBinContent(i)/
00595                         histDetSYS2bbb->GetBinContent(i));
00596      Double_t shift_sys2= data_bgr*fCorr_2 - data_bbb;
00597 
00598      cout<<data_bbb<<" "<<shift_sys1<<" "<<shift_sys2<<"\n";
00599 
00600      // add systematic shifts quadratically and get total error
00601      Double_t errData_total_bbb=
00602         TMath::Sqrt(errData_statbgr_bbb*errData_statbgr_bbb
00603                     +shift_sys1*shift_sys1
00604                     +shift_sys2*shift_sys2);
00605 
00606      // get results from real unfolding
00607      Double_t data_unfold= histUnfoldStat->GetBinContent(i);
00608      Double_t errData_stat_unfold=histUnfoldStat->GetBinError(i);
00609      Double_t errData_statbgr_unfold=histUnfoldStatBgr->GetBinError(i);
00610      Double_t errData_total_unfold=histUnfoldTotal->GetBinError(i);
00611 
00612      // compare
00613      std::cout<<"Bin "<<i<<": true "<<histGenDATA->GetBinContent(i)
00614               <<" unfold: "<<data_unfold
00615               <<" +/- "<<errData_stat_unfold<<" (stat)"
00616               <<" +/- "<<TMath::Sqrt(errData_statbgr_unfold*
00617                                      errData_statbgr_unfold-
00618                                      errData_stat_unfold*
00619                                      errData_stat_unfold)<<" (bgr)"
00620               <<" +/- "<<TMath::Sqrt(errData_total_unfold*
00621                                      errData_total_unfold-
00622                                      errData_statbgr_unfold*
00623                                      errData_statbgr_unfold)<<" (sys)"<<"\n";
00624      std::cout<<"Bin "<<i<<": true "<<histGenDATA->GetBinContent(i)
00625               <<" binbybin: "<<data_bbb
00626               <<" +/- "<<errData_stat_bbb<<" (stat)"
00627               <<" +/- "<<TMath::Sqrt(errData_statbgr_bbb*
00628                                      errData_statbgr_bbb-
00629                                      errData_stat_bbb*
00630                                      errData_stat_bbb)<<" (bgr)"
00631               <<" +/- "<<TMath::Sqrt(errData_total_bbb*
00632                                      errData_total_bbb-
00633                                      errData_statbgr_bbb*
00634                                      errData_statbgr_bbb)<<" (sys)"
00635               <<"\n";
00636   }
00637 }

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