rf704_amplitudefit.cxx

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 "RooTruthModel.h"
00020 #include "RooFormulaVar.h"
00021 #include "RooRealSumPdf.h"
00022 #include "RooPolyVar.h"
00023 #include "RooProduct.h"
00024 #include "TH1.h"
00025 #include "TCanvas.h"
00026 #include "RooPlot.h"
00027 using namespace RooFit ;
00028 
00029 
00030 // Elementary operations on a gaussian PDF
00031 class TestBasic704 : public RooFitTestUnit
00032 {
00033 public: 
00034   TestBasic704(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Amplitude sum operator p.d.f",refFile,writeRef,verbose) {} ;
00035   Bool_t testCode() {
00036 
00037   // 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
00038   // -------------------------------------------------------
00039 
00040   // Observables
00041   RooRealVar t("t","time",-1.,15.);
00042   RooRealVar cosa("cosa","cos(alpha)",-1.,1.);
00043 
00044   // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
00045   RooRealVar tau("tau","#tau",1.5);  
00046   RooRealVar deltaGamma("deltaGamma","deltaGamma", 0.3);
00047   RooTruthModel tm("tm","tm",t) ;
00048   RooFormulaVar coshGBasis("coshGBasis","exp(-@0/ @1)*cosh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
00049   RooFormulaVar sinhGBasis("sinhGBasis","exp(-@0/ @1)*sinh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
00050   RooAbsReal* coshGConv = tm.convolution(&coshGBasis,&t);
00051   RooAbsReal* sinhGConv = tm.convolution(&sinhGBasis,&t);
00052     
00053   // Construct polynomial amplitudes in cos(a) 
00054   RooPolyVar poly1("poly1","poly1",cosa,RooArgList(RooConst(0.5),RooConst(0.2),RooConst(0.2)),0);
00055   RooPolyVar poly2("poly2","poly2",cosa,RooArgList(RooConst(1),RooConst(-0.2),RooConst(3)),0);
00056 
00057   // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
00058   RooProduct  ampl1("ampl1","amplitude 1",RooArgSet(poly1,*coshGConv));
00059   RooProduct  ampl2("ampl2","amplitude 2",RooArgSet(poly2,*sinhGConv));
00060 
00061 
00062 
00063   // C o n s t r u c t   a m p l i t u d e   s u m   p d f 
00064   // -----------------------------------------------------
00065 
00066   // Amplitude strengths
00067   RooRealVar f1("f1","f1",1,0,2) ;
00068   RooRealVar f2("f2","f2",0.5,0,2) ;
00069   
00070   // Construct pdf
00071   RooRealSumPdf pdf("pdf","pdf",RooArgList(ampl1,ampl2),RooArgList(f1,f2)) ;
00072 
00073   // Generate some toy data from pdf
00074   RooDataSet* data = pdf.generate(RooArgSet(t,cosa),10000);
00075 
00076   // Fit pdf to toy data with only amplitude strength floating
00077   pdf.fitTo(*data) ;
00078 
00079 
00080 
00081   // P l o t   a m p l i t u d e   s u m   p d f 
00082   // -------------------------------------------
00083 
00084   // Make 2D plots of amplitudes
00085   TH1* hh_cos = ampl1.createHistogram("hh_cos",t,Binning(50),YVar(cosa,Binning(50))) ;
00086   TH1* hh_sin = ampl2.createHistogram("hh_sin",t,Binning(50),YVar(cosa,Binning(50))) ;
00087   hh_cos->SetLineColor(kBlue) ;
00088   hh_sin->SetLineColor(kBlue) ;
00089 
00090   
00091   // Make projection on t, plot data, pdf and its components
00092   // Note component projections may be larger than sum because amplitudes can be negative
00093   RooPlot* frame1 = t.frame();
00094   data->plotOn(frame1);
00095   pdf.plotOn(frame1);
00096   pdf.plotOn(frame1,Components(ampl1),LineStyle(kDashed));
00097   pdf.plotOn(frame1,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
00098   
00099   // Make projection on cosa, plot data, pdf and its components
00100   // Note that components projection may be larger than sum because amplitudes can be negative
00101   RooPlot* frame2 = cosa.frame();
00102   data->plotOn(frame2);
00103   pdf.plotOn(frame2);
00104   pdf.plotOn(frame2,Components(ampl1),LineStyle(kDashed));
00105   pdf.plotOn(frame2,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
00106   
00107 
00108   regPlot(frame1,"rf704_plot1") ;
00109   regPlot(frame2,"rf704_plot2") ;
00110   regTH(hh_cos,"rf704_hh_cos") ;
00111   regTH(hh_sin,"rf704_hh_sin") ;
00112 
00113   delete data ;
00114   delete coshGConv ;
00115   delete sinhGConv ;
00116 
00117   return kTRUE ;
00118   }
00119 } ;

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