ROOT logo

From $ROOTSYS/tutorials/roofit/rf704_amplitudefit.C

//////////////////////////////////////////////////////////////////////////
//
// 'SPECIAL PDFS' RooFit tutorial macro #704
// 
// Using a p.d.f defined by a sum of real-valued amplitude components
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooTruthModel.h"
#include "RooFormulaVar.h"
#include "RooRealSumPdf.h"
#include "RooPolyVar.h"
#include "RooProduct.h"
#include "TH1.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf704_amplitudefit()
{
  // 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
  // -------------------------------------------------------

  // Observables
  RooRealVar t("t","time",-1.,15.);
  RooRealVar cosa("cosa","cos(alpha)",-1.,1.);

  // Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
  RooRealVar tau("tau","#tau",1.5);  
  RooRealVar deltaGamma("deltaGamma","deltaGamma", 0.3);
  RooTruthModel tm("tm","tm",t) ;
  RooFormulaVar coshGBasis("coshGBasis","exp(-@0/ @1)*cosh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
  RooFormulaVar sinhGBasis("sinhGBasis","exp(-@0/ @1)*sinh(@0*@2/2)",RooArgList(t,tau,deltaGamma));
  RooAbsReal* coshGConv = tm.convolution(&coshGBasis,&t);
  RooAbsReal* sinhGConv = tm.convolution(&sinhGBasis,&t);
    
  // Construct polynomial amplitudes in cos(a) 
  RooPolyVar poly1("poly1","poly1",cosa,RooArgList(RooConst(0.5),RooConst(0.2),RooConst(0.2)),0);
  RooPolyVar poly2("poly2","poly2",cosa,RooArgList(RooConst(1),RooConst(-0.2),RooConst(3)),0);

  // Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
  RooProduct  ampl1("ampl1","amplitude 1",RooArgSet(poly1,*coshGConv));
  RooProduct  ampl2("ampl2","amplitude 2",RooArgSet(poly2,*sinhGConv));



  // C o n s t r u c t   a m p l i t u d e   s u m   p d f 
  // -----------------------------------------------------

  // Amplitude strengths
  RooRealVar f1("f1","f1",1,0,2) ;
  RooRealVar f2("f2","f2",0.5,0,2) ;
  
  // Construct pdf
  RooRealSumPdf pdf("pdf","pdf",RooArgList(ampl1,ampl2),RooArgList(f1,f2)) ;

  // Generate some toy data from pdf
  RooDataSet* data = pdf.generate(RooArgSet(t,cosa),10000);

  // Fit pdf to toy data with only amplitude strength floating
  pdf.fitTo(*data) ;



  // P l o t   a m p l i t u d e   s u m   p d f 
  // -------------------------------------------

  // Make 2D plots of amplitudes
  TH1* hh_cos = ampl1.createHistogram("hh_cos",t,Binning(50),YVar(cosa,Binning(50))) ;
  TH1* hh_sin = ampl2.createHistogram("hh_sin",t,Binning(50),YVar(cosa,Binning(50))) ;
  hh_cos->SetLineColor(kBlue) ;
  hh_sin->SetLineColor(kRed) ;

  
  // Make projection on t, plot data, pdf and its components
  // Note component projections may be larger than sum because amplitudes can be negative
  RooPlot* frame1 = t.frame();
  data->plotOn(frame1);
  pdf.plotOn(frame1);
  pdf.plotOn(frame1,Components(ampl1),LineStyle(kDashed));
  pdf.plotOn(frame1,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
  
  // Make projection on cosa, plot data, pdf and its components
  // Note that components projection may be larger than sum because amplitudes can be negative
  RooPlot* frame2 = cosa.frame();
  data->plotOn(frame2);
  pdf.plotOn(frame2);
  pdf.plotOn(frame2,Components(ampl1),LineStyle(kDashed));
  pdf.plotOn(frame2,Components(ampl2),LineStyle(kDashed),LineColor(kRed));
  

  
  TCanvas* c = new TCanvas("rf704_amplitudefit","rf704_amplitudefit",800,800) ;
  c->Divide(2,2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw();
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw();
  c->cd(3) ; gPad->SetLeftMargin(0.20) ; hh_cos->GetZaxis()->SetTitleOffset(2.3) ; hh_cos->Draw("surf") ;
  c->cd(4) ; gPad->SetLeftMargin(0.20) ; hh_sin->GetZaxis()->SetTitleOffset(2.3) ; hh_sin->Draw("surf") ;
  
  
}
 rf704_amplitudefit.C:1
 rf704_amplitudefit.C:2
 rf704_amplitudefit.C:3
 rf704_amplitudefit.C:4
 rf704_amplitudefit.C:5
 rf704_amplitudefit.C:6
 rf704_amplitudefit.C:7
 rf704_amplitudefit.C:8
 rf704_amplitudefit.C:9
 rf704_amplitudefit.C:10
 rf704_amplitudefit.C:11
 rf704_amplitudefit.C:12
 rf704_amplitudefit.C:13
 rf704_amplitudefit.C:14
 rf704_amplitudefit.C:15
 rf704_amplitudefit.C:16
 rf704_amplitudefit.C:17
 rf704_amplitudefit.C:18
 rf704_amplitudefit.C:19
 rf704_amplitudefit.C:20
 rf704_amplitudefit.C:21
 rf704_amplitudefit.C:22
 rf704_amplitudefit.C:23
 rf704_amplitudefit.C:24
 rf704_amplitudefit.C:25
 rf704_amplitudefit.C:26
 rf704_amplitudefit.C:27
 rf704_amplitudefit.C:28
 rf704_amplitudefit.C:29
 rf704_amplitudefit.C:30
 rf704_amplitudefit.C:31
 rf704_amplitudefit.C:32
 rf704_amplitudefit.C:33
 rf704_amplitudefit.C:34
 rf704_amplitudefit.C:35
 rf704_amplitudefit.C:36
 rf704_amplitudefit.C:37
 rf704_amplitudefit.C:38
 rf704_amplitudefit.C:39
 rf704_amplitudefit.C:40
 rf704_amplitudefit.C:41
 rf704_amplitudefit.C:42
 rf704_amplitudefit.C:43
 rf704_amplitudefit.C:44
 rf704_amplitudefit.C:45
 rf704_amplitudefit.C:46
 rf704_amplitudefit.C:47
 rf704_amplitudefit.C:48
 rf704_amplitudefit.C:49
 rf704_amplitudefit.C:50
 rf704_amplitudefit.C:51
 rf704_amplitudefit.C:52
 rf704_amplitudefit.C:53
 rf704_amplitudefit.C:54
 rf704_amplitudefit.C:55
 rf704_amplitudefit.C:56
 rf704_amplitudefit.C:57
 rf704_amplitudefit.C:58
 rf704_amplitudefit.C:59
 rf704_amplitudefit.C:60
 rf704_amplitudefit.C:61
 rf704_amplitudefit.C:62
 rf704_amplitudefit.C:63
 rf704_amplitudefit.C:64
 rf704_amplitudefit.C:65
 rf704_amplitudefit.C:66
 rf704_amplitudefit.C:67
 rf704_amplitudefit.C:68
 rf704_amplitudefit.C:69
 rf704_amplitudefit.C:70
 rf704_amplitudefit.C:71
 rf704_amplitudefit.C:72
 rf704_amplitudefit.C:73
 rf704_amplitudefit.C:74
 rf704_amplitudefit.C:75
 rf704_amplitudefit.C:76
 rf704_amplitudefit.C:77
 rf704_amplitudefit.C:78
 rf704_amplitudefit.C:79
 rf704_amplitudefit.C:80
 rf704_amplitudefit.C:81
 rf704_amplitudefit.C:82
 rf704_amplitudefit.C:83
 rf704_amplitudefit.C:84
 rf704_amplitudefit.C:85
 rf704_amplitudefit.C:86
 rf704_amplitudefit.C:87
 rf704_amplitudefit.C:88
 rf704_amplitudefit.C:89
 rf704_amplitudefit.C:90
 rf704_amplitudefit.C:91
 rf704_amplitudefit.C:92
 rf704_amplitudefit.C:93
 rf704_amplitudefit.C:94
 rf704_amplitudefit.C:95
 rf704_amplitudefit.C:96
 rf704_amplitudefit.C:97
 rf704_amplitudefit.C:98
 rf704_amplitudefit.C:99
 rf704_amplitudefit.C:100
 rf704_amplitudefit.C:101
 rf704_amplitudefit.C:102
 rf704_amplitudefit.C:103
 rf704_amplitudefit.C:104
 rf704_amplitudefit.C:105
 rf704_amplitudefit.C:106
 rf704_amplitudefit.C:107
 rf704_amplitudefit.C:108
 rf704_amplitudefit.C:109
 rf704_amplitudefit.C:110
 rf704_amplitudefit.C:111
 rf704_amplitudefit.C:112
 rf704_amplitudefit.C:113
 rf704_amplitudefit.C:114
 rf704_amplitudefit.C:115