00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020 #ifndef __CINT__
00021 #include "RooGlobalFunc.h"
00022 #endif
00023 #include "RooRealVar.h"
00024 #include "RooDataSet.h"
00025 #include "RooGaussian.h"
00026 #include "RooGenericPdf.h"
00027 #include "RooFormulaVar.h"
00028 #include "RooFFTConvPdf.h"
00029 #include "RooPlot.h"
00030 #include "TCanvas.h"
00031 #include "TAxis.h"
00032 #include "TH1.h"
00033 using namespace RooFit ;
00034
00035
00036
00037 void rf210_angularconv()
00038 {
00039
00040
00041
00042
00043 RooRealVar psi("psi","psi",0,3.14159268) ;
00044
00045
00046 RooGenericPdf Tpsi("Tpsi","1+sin(2*@0)",psi) ;
00047
00048
00049 RooRealVar gbias("gbias","gbias",0.2,0.,1) ;
00050 RooRealVar greso("greso","greso",0.3,0.1,1.0) ;
00051 RooGaussian Rpsi("Rpsi","Rpsi",psi,gbias,greso) ;
00052
00053
00054 RooRealVar cpsi("cpsi","cos(psi)",-1,1) ;
00055 RooFormulaVar psif("psif","acos(cpsi)",cpsi) ;
00056
00057
00058 RooGenericPdf Tcpsi("T","1+sin(2*@0)",psif) ;
00059
00060
00061
00062
00063
00064
00065
00066 RooFFTConvPdf Mpsi("Mf","Mf",psi,Tpsi,Rpsi) ;
00067
00068
00069 Mpsi.setBufferFraction(0) ;
00070
00071
00072
00073
00074
00075
00076
00077 RooDataSet* data_psi = Mpsi.generate(psi,10000) ;
00078
00079
00080 Mpsi.fitTo(*data_psi) ;
00081
00082
00083 RooPlot* frame1 = psi.frame(Title("Cyclical convolution in angle psi")) ;
00084 data_psi->plotOn(frame1) ;
00085 Mpsi.plotOn(frame1) ;
00086
00087
00088 Tpsi.plotOn(frame1,LineColor(kRed)) ;
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 RooFFTConvPdf Mcpsi("Mf","Mf",psif,psi,Tpsi,Rpsi) ;
00101
00102
00103 Mcpsi.setBufferFraction(0) ;
00104
00105
00106
00107
00108
00109
00110 RooDataSet* data_cpsi = Mcpsi.generate(cpsi,10000) ;
00111
00112
00113 Mcpsi.fitTo(*data_cpsi) ;
00114
00115
00116 RooPlot* frame2 = cpsi.frame(Title("Same convolution in psi, expressed in cos(psi)")) ;
00117 data_cpsi->plotOn(frame2) ;
00118 Mcpsi.plotOn(frame2) ;
00119
00120
00121 Tcpsi.plotOn(frame2,LineColor(kRed)) ;
00122
00123
00124
00125
00126 TCanvas* c = new TCanvas("rf210_angularconv","rf210_angularconv",800,400) ;
00127 c->Divide(2) ;
00128 c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
00129 c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
00130
00131 }