00001
00002
00003
00004
00005
00006
00007
00008
00009 #define ECM 14000
00010
00011
00012 #include <iostream>
00013 using namespace std;
00014
00015 #include "TStopwatch.h"
00016 #include "TH1F.h"
00017 #include "TH2F.h"
00018 #include "TF1.h"
00019 #include "TPythia6.h"
00020 #include "TCanvas.h"
00021 #include "TMath.h"
00022 #include "TRandom2.h"
00023 #include "TLorentzVector.h"
00024 #include "TVector3.h"
00025 #include "TVector2.h"
00026 #include "TFile.h"
00027 #include "TROOT.h"
00028
00029
00030
00031
00032 int makePyEvts2(int nEvents,float Higgsmass = 300 )
00033 {
00034
00035
00036 TPythia6* pythia = new TPythia6;
00037
00038
00039
00040 pythia->SetMSEL(0);
00041 pythia->SetMSUB(102,1);
00042 pythia->SetPMAS(25,1,Higgsmass);
00043 pythia->SetMSTP(61,0);
00044
00045
00046
00047 int KF = 23;
00048 int idcmin = pythia->GetMDCY(KF,2);
00049 int idcmax = idcmin + pythia->GetMDCY(KF,3);
00050 cout<<"MyOutput: Z idcmin="<<idcmin<<" idcmax="<<idcmax<<endl;
00051 for (int idc = idcmin; idc < idcmax; idc++)
00052 {
00053 pythia->SetMDME(idc,1,0);
00054 }
00055 pythia->SetMDME(184,1,1);
00056
00057
00058 KF = 25;
00059 idcmin = pythia->GetMDCY(KF,2);
00060 idcmax = idcmin + pythia->GetMDCY(KF,3);
00061 cout<<"MyOutput: Higgs idcmin="<<idcmin<<" idcmax="<<idcmax<<endl;
00062 for (Int_t idc = idcmin; idc < idcmax; idc++)
00063 {
00064 pythia->SetMDME(idc,1,0);
00065 }
00066 pythia->SetMDME(225,1,1);
00067
00068
00069 pythia->Initialize("cms", "p", "p", ECM);
00070
00071 float higgswidth = pythia->GetPMAS(25,2);
00072 cout<<"The higgs width from Pythia is "<<higgswidth<<endl;
00073
00074
00075 TStopwatch watcher;
00076 watcher.Reset();
00077 watcher.Start();
00078
00079
00080 float Higgsmass_max = 3000;
00081
00082
00083
00084
00085 int NbinIVplot = 100;
00086
00087
00088 float Higgsmaxplot = 3000;
00089 float Higgsminplot = 0.0;
00090 int Nbins = 300;
00091
00092 if (Higgsmass < 191)
00093 {
00094 Nbins = 500;
00095 Higgsminplot = Higgsmass - 5*higgswidth;
00096 Higgsmaxplot = Higgsmass + 5*higgswidth;
00097 }
00098
00099
00100
00101 TH1F *Hmass = new TH1F("Hmass","Higgs Mass",Nbins,Higgsminplot,Higgsmaxplot);
00102
00103 TH1F *Z1plot = new TH1F("Z1mass","Z1mass",Nbins,0,150);
00104 TH1F *Z2plot = new TH1F("Z2mass","Z2mass",Nbins,0,150);
00105
00106 TH1F *Cosine1plot = new TH1F("Cosine1graph","Cos(#theta) angle of muon1 in Z1 rest frame",NbinIVplot,-1,1);
00107 TH1F *Cosine2plot = new TH1F("Cosine2graph","Cos(#theta) angle of muon2 in Z2 rest frame",NbinIVplot,-1,1);
00108 TH2F *Correlationplot = new TH2F("Correlationplot","Theta angle muon 1 versus theta angle muon 2",NbinIVplot,-1,1,NbinIVplot,-1,1);
00109
00110 TH1F *Phiplot1 = new TH1F("Phiplot1","Phi angle of muon 1 in Z1 rest frame",NbinIVplot,0,2*TMath::Pi());
00111 TH1F *Phiplot2 = new TH1F("Phiplot2","Phi angle of muon 2 in Z2 rest frame",NbinIVplot,0,2*TMath::Pi());
00112
00113
00114 Int_t KF_muon =13;
00115
00116
00117
00118 TRandom2 R;
00119
00120
00121 for (Int_t ievt = 0; ievt < nEvents; ievt++)
00122 {
00123
00124 if (ievt % 100 == 0)
00125 cout << "Event # " << ievt << endl;
00126
00127
00128 pythia->GenerateEvent();
00129
00130
00131
00132
00133 TLorentzVector muonsum;
00134 TLorentzVector muon1plus, muon2plus, muon1minus , muon2minus;
00135
00136 int Nmuonplus=0;
00137 int Nmuonminus=0;
00138
00139 int Nz=0;
00140
00141
00142
00143 Int_t N = pythia->GetN();
00144
00145
00146
00147
00148
00149 Int_t N_muons =0;
00150 for (Int_t itk =1; itk<=N; itk++)
00151 {
00152
00153 Int_t KS = pythia->GetK(itk,1);
00154
00155 if ( KS!=21) continue;
00156
00157 Int_t KF = pythia->GetK(itk,2);
00158
00159
00160
00161 if ((abs(KF)== KF_muon))
00162 {
00163 TLorentzVector lvMuon = TLorentzVector(pythia->GetP(itk,1), pythia->GetP(itk,2), pythia->GetP(itk,3), pythia->GetP(itk,4));
00164
00165
00166
00167 if(KS==21)
00168 {
00169 N_muons++;
00170
00171 muonsum+=lvMuon;
00172
00173
00174 if (KF<0)
00175 {
00176 if (Nmuonplus==1)
00177 {
00178 muon2plus = lvMuon;
00179 }
00180 if (Nmuonplus==0)
00181 {
00182 muon1plus = lvMuon;
00183 Nmuonplus=1;
00184 }
00185 }
00186 else
00187 {
00188 if (Nmuonminus==1)
00189 {
00190 muon2minus = lvMuon;
00191 }
00192 if (Nmuonminus==0)
00193 {
00194 muon1minus = lvMuon;
00195 Nmuonminus=1;
00196 }
00197 }
00198 }
00199 else
00200 {
00201
00202 }
00203 }
00204
00205
00206 }
00207
00208
00209
00210
00211
00212
00213
00214 TLorentzVector Z1,Z2;
00215 Z1 = muon1plus + muon1minus;
00216 Z2 = muon2plus + muon2minus;
00217 Z1plot->Fill(Z1.M());
00218 Z2plot->Fill(Z2.M());
00219
00220
00221
00222
00223
00224 TVector3 boostZ1,boostZ2;
00225
00226
00227
00228
00229 boostZ1 = - Z1.BoostVector();
00230 boostZ2 = - Z2.BoostVector();
00231
00232
00233
00234
00235
00236 TLorentzVector Z1_star,muon1_star,antimuon1_star;
00237 Z1_star = Z1;
00238 muon1_star = muon1minus;
00239 antimuon1_star = muon1plus;
00240
00241 Z1_star.Boost(boostZ1);
00242 muon1_star.Boost(boostZ1);
00243 antimuon1_star.Boost(boostZ1);
00244
00245
00246
00247 TLorentzVector Z2_star,muon2_star,antimuon2_star;
00248 Z2_star = Z2;
00249 muon2_star = muon2minus;
00250 antimuon2_star = muon2plus;
00251
00252 Z2_star.Boost(boostZ2);
00253 muon2_star.Boost(boostZ2);
00254 antimuon2_star.Boost(boostZ2);
00255
00256
00257 TVector3 uZ1=-boostZ1.Unit();
00258 TVector3 uZ2=-boostZ2.Unit();
00259
00260 float cosine1, cosine2;
00261 TVector3 unit_muon1_star = muon1_star.Vect().Unit();
00262 cosine1 = uZ1.Dot(unit_muon1_star);
00263
00264
00265
00266
00267 TVector3 unit_muon2_star = muon2_star.Vect().Unit();
00268 cosine2 = uZ2.Dot(unit_muon2_star);
00269
00270 Cosine1plot->Fill(cosine1);
00271 Cosine2plot->Fill(cosine2);
00272 Correlationplot->Fill(cosine1,cosine2);
00273
00274 TVector3 x1(1.0,0.0,0.0),y1(0.0,1.0,0.0);
00275 x1.RotateUz(uZ1);
00276 y1.RotateUz(uZ1);
00277
00278
00279 TVector2 muon1_star2d(x1.Dot(unit_muon1_star),y1.Dot(unit_muon1_star));
00280 Phiplot1->Fill(muon1_star2d.Phi());
00281
00282 TVector3 x2(1.0,0.0,0.0),y2(0.0,1.0,0.0);
00283 x2.RotateUz(uZ2);
00284 y2.RotateUz(uZ2);
00285
00286
00287 TVector2 muon2_star2d(x2.Dot(unit_muon2_star),y2.Dot(unit_muon2_star));
00288 Phiplot2->Fill(muon2_star2d.Phi());
00289
00290
00291
00292 Hmass->Fill(muonsum.M());
00293
00294 }
00295
00296
00297
00298
00299
00300
00301
00302 char name[100];
00303
00304
00305
00306 TCanvas * c3 = new TCanvas("c3","**",800,800);
00307 c3->Divide(2,2);
00308 c3->SetFillColor(9);
00309 c3->cd(1);
00310 Cosine1plot->Draw();
00311
00312 Cosine1plot->SetTitle("Cos(#theta) of muon1 in Z1 rest frame");
00313
00314 TF1 *f2 = new TF1("f2","[0]*([1]*(3/4)*(1-x**2)+(3/8)*(1-[1])*(1+x**2))",-1,1);
00315 TF1 *f3 = new TF1("f3","[0]*( (1/(1+[1]))*(3/4)*(1-x**2) + (3/8)*([1]/(1+[1]))*(1+x**2) )",-1,1);
00316 TF1 *f4 = new TF1("f4","[0]*( ((-1-[1])/([1]-3)) * (3/4)*(1-x**2) + (3/8)*((2*[1]-2)/([1]-3))*(1+x**2) )",-1,1);
00317
00318 f2->SetParNames("Normalization","alpha");
00319 f3->SetParNames("Normalization","R");
00320 f4->SetParNames("Normalization","tullyR");
00321
00322 f2->SetParameter(1,0.9);
00323 f3->SetParameter(1,0.1);
00324 f4->SetParameter(1,0.8);
00325
00326 f2->FixParameter(0,Cosine1plot->Integral("bin width"));
00327 f3->FixParameter(0,Cosine1plot->Integral("bin width"));
00328 f4->FixParameter(0,Cosine1plot->Integral("bin width"));
00329
00330 Cosine1plot->Fit("f2","L");
00331 Cosine1plot->Fit("f3","L+");
00332 Cosine1plot->Fit("f4","L+");
00333
00334
00335
00336
00337 c3->cd(2);
00338 Cosine2plot->Draw();
00339 c3->cd(3);
00340 Correlationplot->Draw();
00341
00342
00343
00344
00345
00346 TCanvas *c4 = new TCanvas("c4","**",800,800);
00347 c4->Divide(2,2);
00348 c4->SetFillColor(9);
00349 c4->cd(1);
00350 Phiplot1->Draw();
00351 c4->cd(2);
00352 Phiplot2->Draw();
00353
00354
00355 TCanvas * c5 = new TCanvas("c5","**",800,800);
00356 c5->SetFillColor(9);
00357 Hmass->Draw("e");
00358
00359 TCanvas * c6 = new TCanvas("c6","**",800,800);
00360 c6->cd();
00361 Z1plot->Draw();
00362
00363 TCanvas * c7 = new TCanvas("c6","**",800,800);
00364 c7->cd();
00365 Z2plot->Draw();
00366
00367 sprintf(name,"histograms_correct_HtoZZ_angular__incorrectboost_2_Hm%5.2f_Nevt%d.root",Higgsmass,nEvents);
00368 TFile f(name,"recreate");
00369 f.cd();
00370 int N=gROOT->GetList()->GetSize();
00371 for (int key=0; key<N ; key++)
00372 {
00373
00374 TObject * obj = gROOT->GetList()->At(key);
00375 if (obj->InheritsFrom("TH1F") || obj->InheritsFrom("TH2F"))
00376 {
00377 cout <<"write "<<obj->GetName()<<endl;
00378 obj->Write();
00379
00380 }
00381 }
00382 f.Close();
00383
00384
00385 for (int key=N-1; key>=0 ; key--)
00386 {
00387
00388 TObject * obj = gROOT->GetList()->At(key);
00389 cout <<"delete "<<obj->GetName()<<endl;
00390 if (obj) delete obj;
00391 cout <<"deleted "<<gROOT->GetList()->GetSize()<<" objects remaining"<<endl;
00392 }
00393
00394 watcher.Stop();
00395 cout <<"Time spent for "<<nEvents<<" events :"<< watcher.RealTime()<<" [s]"<<endl;
00396
00397 return -97534;
00398
00399 }
00400