rf704_amplitudefit.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'SPECIAL PDFS' RooFit tutorial macro #704
00004 // 
00005 // Using a p.d.f defined by a sum of real-valued amplitude components
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooGaussian.h"
00019 #include "RooConstVar.h"
00020 #include "RooTruthModel.h"
00021 #include "RooFormulaVar.h"
00022 #include "RooRealSumPdf.h"
00023 #include "RooPolyVar.h"
00024 #include "RooProduct.h"
00025 #include "TH1.h"
00026 #include "TCanvas.h"
00027 #include "TAxis.h"
00028 #include "RooPlot.h"
00029 using namespace RooFit ;
00030 
00031 
00032 void rf704_amplitudefit()
00033 {
00034   // S e t u p   2 D   a m p l i t u d e   f u n c t i o n s
00035   // -------------------------------------------------------
00036 
00037   // Observables
00038   RooRealVar t("t","time",-1.,15.);
00039   RooRealVar cosa("cosa","cos(alpha)",-1.,1.);
00040 
00041   // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
00042   RooRealVar tau("tau","#tau",1.5);  
00043   RooRealVar deltaGamma("deltaGamma","deltaGamma", 0.3);
00044   RooTruthModel tm("tm","tm",t) ;
00045   RooFormulaVar coshGBasis("coshGBasis","exp(-@0/ @1)*cosh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
00046   RooFormulaVar sinhGBasis("sinhGBasis","exp(-@0/ @1)*sinh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
00047   RooAbsReal* coshGConv = tm.convolution(&coshGBasis,&t);
00048   RooAbsReal* sinhGConv = tm.convolution(&sinhGBasis,&t);
00049     
00050   // Construct polynomial amplitudes in cos(a) 
00051   RooPolyVar poly1("poly1","poly1",cosa,RooArgList(RooConst(0.5),RooConst(0.2),RooConst(0.2)),0);
00052   RooPolyVar poly2("poly2","poly2",cosa,RooArgList(RooConst(1),RooConst(-0.2),RooConst(3)),0);
00053 
00054   // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
00055   RooProduct  ampl1("ampl1","amplitude 1",RooArgSet(poly1,*coshGConv));
00056   RooProduct  ampl2("ampl2","amplitude 2",RooArgSet(poly2,*sinhGConv));
00057 
00058 
00059 
00060   // C o n s t r u c t   a m p l i t u d e   s u m   p d f 
00061   // -----------------------------------------------------
00062 
00063   // Amplitude strengths
00064   RooRealVar f1("f1","f1",1,0,2) ;
00065   RooRealVar f2("f2","f2",0.5,0,2) ;
00066   
00067   // Construct pdf
00068   RooRealSumPdf pdf("pdf","pdf",RooArgList(ampl1,ampl2),RooArgList(f1,f2)) ;
00069 
00070   // Generate some toy data from pdf
00071   RooDataSet* data = pdf.generate(RooArgSet(t,cosa),10000);
00072 
00073   // Fit pdf to toy data with only amplitude strength floating
00074   pdf.fitTo(*data) ;
00075 
00076 
00077 
00078   // P l o t   a m p l i t u d e   s u m   p d f 
00079   // -------------------------------------------
00080 
00081   // Make 2D plots of amplitudes
00082   TH1* hh_cos = ampl1.createHistogram("hh_cos",t,Binning(50),YVar(cosa,Binning(50))) ;
00083   TH1* hh_sin = ampl2.createHistogram("hh_sin",t,Binning(50),YVar(cosa,Binning(50))) ;
00084   hh_cos->SetLineColor(kBlue) ;
00085   hh_sin->SetLineColor(kRed) ;
00086 
00087   
00088   // Make projection on t, plot data, pdf and its components
00089   // Note component projections may be larger than sum because amplitudes can be negative
00090   RooPlot* frame1 = t.frame();
00091   data->plotOn(frame1);
00092   pdf.plotOn(frame1);
00093   pdf.plotOn(frame1,Components(ampl1),LineStyle(kDashed));
00094   pdf.plotOn(frame1,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
00095   
00096   // Make projection on cosa, plot data, pdf and its components
00097   // Note that components projection may be larger than sum because amplitudes can be negative
00098   RooPlot* frame2 = cosa.frame();
00099   data->plotOn(frame2);
00100   pdf.plotOn(frame2);
00101   pdf.plotOn(frame2,Components(ampl1),LineStyle(kDashed));
00102   pdf.plotOn(frame2,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
00103   
00104 
00105   
00106   TCanvas* c = new TCanvas("rf704_amplitudefit","rf704_amplitudefit",800,800) ;
00107   c->Divide(2,2) ;
00108   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw();
00109   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw();
00110   c->cd(3) ; gPad->SetLeftMargin(0.20) ; hh_cos->GetZaxis()->SetTitleOffset(2.3) ; hh_cos->Draw("surf") ;
00111   c->cd(4) ; gPad->SetLeftMargin(0.20) ; hh_sin->GetZaxis()->SetTitleOffset(2.3) ; hh_sin->Draw("surf") ;
00112   
00113   
00114 }

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