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 makePyEvts(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 pythia->SetMSTP(25,2);
00045
00046
00047
00048
00049 int KF = 23;
00050 int idcmin = pythia->GetMDCY(KF,2);
00051 int idcmax = idcmin + pythia->GetMDCY(KF,3);
00052 cout<<"MyOutput: Z idcmin="<<idcmin<<" idcmax="<<idcmax<<endl;
00053 for (int idc = idcmin; idc < idcmax; idc++)
00054 {
00055 pythia->SetMDME(idc,1,0);
00056 }
00057 pythia->SetMDME(184,1,1);
00058
00059
00060 KF = 25;
00061 idcmin = pythia->GetMDCY(KF,2);
00062 idcmax = idcmin + pythia->GetMDCY(KF,3);
00063 cout<<"MyOutput: Higgs idcmin="<<idcmin<<" idcmax="<<idcmax<<endl;
00064 for (Int_t idc = idcmin; idc < idcmax; idc++)
00065 {
00066 pythia->SetMDME(idc,1,0);
00067 }
00068 pythia->SetMDME(225,1,1);
00069
00070
00071 pythia->Initialize("cms", "p", "p", ECM);
00072
00073 float higgswidth = pythia->GetPMAS(25,2);
00074 cout<<"The higgs width from Pythia is "<<higgswidth<<endl;
00075
00076
00077 TStopwatch watcher;
00078 watcher.Reset();
00079 watcher.Start();
00080
00081
00082 float Higgsmass_max = 3000;
00083
00084
00085
00086
00087 int NbinIVplot = 100;
00088
00089
00090 float Higgsmaxplot = 3000;
00091 float Higgsminplot = 0.0;
00092 int Nbins = 300;
00093
00094 if (Higgsmass < 191)
00095 {
00096 Nbins = 500;
00097 Higgsminplot = Higgsmass - 5*higgswidth;
00098 Higgsmaxplot = Higgsmass + 5*higgswidth;
00099 }
00100
00101
00102
00103 TH1F *Hmass = new TH1F("Hmass","Higgs Mass",Nbins,Higgsminplot,Higgsmaxplot);
00104
00105 TH1F *Z1plot = new TH1F("Z1mass","Z1mass",Nbins,0,150);
00106 TH1F *Z2plot = new TH1F("Z2mass","Z2mass",Nbins,0,150);
00107
00108 TH1F *Cosine1plot = new TH1F("Cosine1graph","Cos(#theta) angle of muon1 in Z1 rest frame",NbinIVplot,-1,1);
00109 TH1F *Cosine2plot = new TH1F("Cosine2graph","Cos(#theta) angle of muon2 in Z2 rest frame",NbinIVplot,-1,1);
00110 TH2F *Correlationplot = new TH2F("Correlationplot","Theta angle muon 1 versus theta angle muon 2",NbinIVplot,-1,1,NbinIVplot,-1,1);
00111
00112 TH1F *Phiplot1 = new TH1F("Phiplot1","Phi angle of muon 1 in Z1 rest frame",NbinIVplot,0,2*TMath::Pi());
00113 TH1F *Phiplot2 = new TH1F("Phiplot2","Phi angle of muon 2 in Z2 rest frame",NbinIVplot,0,2*TMath::Pi());
00114
00115
00116 Int_t KF_muon =13;
00117
00118
00119
00120 TRandom2 R;
00121
00122
00123 for (Int_t ievt = 0; ievt < nEvents; ievt++)
00124 {
00125
00126 if (ievt % 100 == 0)
00127 cout << "Event # " << ievt << endl;
00128
00129
00130 pythia->GenerateEvent();
00131
00132
00133
00134
00135 TLorentzVector muonsum;
00136 TLorentzVector muon1plus, muon2plus, muon1minus , muon2minus;
00137
00138 int Nmuonplus=0;
00139 int Nmuonminus=0;
00140
00141 int Nz=0;
00142
00143
00144
00145
00146 Int_t N = pythia->GetN();
00147
00148
00149
00150 Int_t N_muons =0;
00151
00152 for (Int_t itk =1; itk<=N; itk++)
00153 {
00154
00155 Int_t KS = pythia->GetK(itk,1);
00156
00157 if ( KS!=21) continue;
00158
00159 Int_t KF = pythia->GetK(itk,2);
00160
00161 if ((abs(KF)== KF_muon))
00162 {
00163
00164 TLorentzVector lvMuon = TLorentzVector(pythia->GetP(itk,1), pythia->GetP(itk,2), pythia->GetP(itk,3), pythia->GetP(itk,4));
00165
00166 if(KS==21)
00167 {
00168 N_muons++;
00169
00170 muonsum+=lvMuon;
00171
00172
00173 if (KF<0)
00174 {
00175 if (Nmuonplus==1)
00176 {
00177 muon2plus = lvMuon;
00178 }
00179 if (Nmuonplus==0)
00180 {
00181 muon1plus = lvMuon;
00182 Nmuonplus=1;
00183 }
00184 }
00185 else
00186 {
00187 if (Nmuonminus==1)
00188 {
00189 muon2minus = lvMuon;
00190 }
00191 if (Nmuonminus==0)
00192 {
00193 muon1minus = lvMuon;
00194 Nmuonminus=1;
00195 }
00196 }
00197 }
00198 else
00199 {
00200
00201 }
00202 }
00203
00204
00205 }
00206
00207
00208
00209
00210
00211
00212
00213 TLorentzVector Z1,Z2;
00214 Z1 = muon1plus + muon1minus;
00215 Z2 = muon2plus + muon2minus;
00216 Z1plot->Fill(Z1.M());
00217 Z2plot->Fill(Z2.M());
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253 TVector3 boostH,boostZ1_hrest,boostZ2_hrest;
00254 TLorentzVector H,H_hrest,Z1_hrest,Z2_hrest,muon1_hrest,muon2_hrest,muon1_star,muon2_star,Z1_star,Z2_star;
00255 H = Z1+Z2;
00256 boostH = -H.BoostVector();
00257 H_hrest = H;
00258 Z1_hrest = Z1;
00259 Z2_hrest = Z2;
00260 muon1_hrest = muon1minus;
00261 muon2_hrest = muon2minus;
00262 H_hrest.Boost(boostH);
00263 Z1_hrest.Boost(boostH);
00264 Z2_hrest.Boost(boostH);
00265 muon1_hrest.Boost(boostH);
00266 muon2_hrest.Boost(boostH);
00267
00268 boostZ1_hrest = -Z1_hrest.BoostVector();
00269 boostZ2_hrest = -Z2_hrest.BoostVector();
00270
00271 Z1_star = Z1_hrest;
00272 Z2_star = Z2_hrest;
00273
00274 muon1_star = muon1_hrest;
00275 muon2_star = muon2_hrest;
00276
00277 muon1_star.Boost(boostZ1_hrest);
00278 muon2_star.Boost(boostZ2_hrest);
00279
00280 TVector3 uZ1 = Z1_hrest.Vect().Unit();
00281 TVector3 uZ2 = Z2_hrest.Vect().Unit();
00282
00283 float cosine1, cosine2;
00284 TVector3 unit_muon1_star = muon1_star.Vect().Unit();
00285 cosine1 = uZ1.Dot(unit_muon1_star);
00286
00287
00288
00289
00290 TVector3 unit_muon2_star = muon2_star.Vect().Unit();
00291 cosine2 = uZ2.Dot(unit_muon2_star);
00292
00293 Cosine1plot->Fill(cosine1);
00294 Cosine2plot->Fill(cosine2);
00295 Correlationplot->Fill(cosine1,cosine2);
00296
00297 TVector3 x1(1.0,0.0,0.0),y1(0.0,1.0,0.0);
00298 x1.RotateUz(uZ1);
00299 y1.RotateUz(uZ1);
00300
00301
00302 TVector2 muon1_star2d(x1.Dot(unit_muon1_star),y1.Dot(unit_muon1_star));
00303 Phiplot1->Fill(muon1_star2d.Phi());
00304
00305 TVector3 x2(1.0,0.0,0.0),y2(0.0,1.0,0.0);
00306 x2.RotateUz(uZ2);
00307 y2.RotateUz(uZ2);
00308
00309
00310 TVector2 muon2_star2d(x2.Dot(unit_muon2_star),y2.Dot(unit_muon2_star));
00311 Phiplot2->Fill(muon2_star2d.Phi());
00312
00313
00314 bool didImessmuonplus=false;
00315 bool didImessmuonminus=false;
00316
00317
00318 if(R.Rndm()>0.5)
00319 {
00320 TLorentzVector temp;
00321 temp=muon1plus;
00322 muon1plus=muon2plus;
00323 muon2plus=temp;
00324 didImessmuonplus=true;
00325 }
00326
00327 if(R.Rndm()>0.5)
00328 {
00329 TLorentzVector temp;
00330 temp=muon1minus;
00331 muon1minus=muon2minus;
00332 muon2minus=temp;
00333 didImessmuonminus=true;
00334 }
00335
00336 int MM=didImessmuonminus;
00337 int MP=didImessmuonplus;
00338
00339
00340 bool condition=(MM || MP) && (!MP || !MM);
00341 bool auxcond=((MP+MM)&1);
00342 if (auxcond!=condition)
00343 cout <<condition<<" "<<((MP+MM)&1)<<endl;
00344
00345
00346
00347 Hmass->Fill(muonsum.M());
00348
00349 }
00350
00351
00352
00353
00354
00355
00356
00357 char name[100];
00358
00359
00360
00361 TCanvas * c3 = new TCanvas("c3","**",800,800);
00362 c3->Divide(2,2);
00363 c3->SetFillColor(9);
00364 c3->cd(1);
00365 Cosine1plot->Draw();
00366
00367 Cosine1plot->SetTitle("Cos(#theta) of muon1 in Z1 rest frame");
00368
00369 TF1 *f2 = new TF1("f2","[0]*([1]*(3/4)*(1-x**2)+(3/8)*(1-[1])*(1+x**2))",-1,1);
00370 TF1 *f3 = new TF1("f3","[0]*( (1/(1+[1]))*(3/4)*(1-x**2) + (3/8)*([1]/(1+[1]))*(1+x**2) )",-1,1);
00371 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);
00372
00373 f2->SetParNames("Normalization","alpha");
00374 f3->SetParNames("Normalization","R");
00375 f4->SetParNames("Normalization","tullyR");
00376
00377 f2->SetParameter(1,0.9);
00378 f3->SetParameter(1,0.1);
00379 f4->SetParameter(1,0.8);
00380
00381 f2->FixParameter(0,Cosine1plot->Integral("bin width"));
00382 f3->FixParameter(0,Cosine1plot->Integral("bin width"));
00383 f4->FixParameter(0,Cosine1plot->Integral("bin width"));
00384
00385 Cosine1plot->Fit("f2","L");
00386 Cosine1plot->Fit("f3","L+");
00387 Cosine1plot->Fit("f4","L+");
00388
00389
00390
00391
00392 c3->cd(2);
00393 Cosine2plot->Draw();
00394 c3->cd(3);
00395 Correlationplot->Draw();
00396
00397
00398
00399
00400
00401 TCanvas *c4 = new TCanvas("c4","**",800,800);
00402 c4->Divide(2,2);
00403 c4->SetFillColor(9);
00404 c4->cd(1);
00405 Phiplot1->Draw();
00406 c4->cd(2);
00407 Phiplot2->Draw();
00408
00409
00410 TCanvas * c5 = new TCanvas("c5","**",800,800);
00411 c5->SetFillColor(9);
00412 Hmass->Draw("e");
00413
00414 TCanvas * c6 = new TCanvas("c6","**",800,800);
00415 c6->cd();
00416 Z1plot->Draw();
00417
00418 TCanvas * c7 = new TCanvas("c6","**",800,800);
00419 c7->cd();
00420 Z2plot->Draw();
00421
00422 sprintf(name,"histograms_correct_Hpseudo_toZZ_angular_Hm%5.2f_Nevt%d.root",Higgsmass,nEvents);
00423 TFile f(name,"recreate");
00424 f.cd();
00425 int N=gROOT->GetList()->GetSize();
00426 for (int key=0; key<N ; key++)
00427 {
00428
00429 TObject * obj = gROOT->GetList()->At(key);
00430 if (obj->InheritsFrom("TH1F") || obj->InheritsFrom("TH2F"))
00431 {
00432 cout <<"write "<<obj->GetName()<<endl;
00433 obj->Write();
00434
00435 }
00436 }
00437 f.Close();
00438
00439
00440 for (int key=N-1; key>=0 ; key--)
00441 {
00442
00443 TObject * obj = gROOT->GetList()->At(key);
00444 cout <<"delete "<<obj->GetName()<<endl;
00445 if (obj) delete obj;
00446 cout <<"deleted "<<gROOT->GetList()->GetSize()<<" objects remaining"<<endl;
00447 }
00448
00449 watcher.Stop();
00450 cout <<"Time spent for "<<nEvents<<" events :"<< watcher.RealTime()<<" [s]"<<endl;
00451
00452 return -97534;
00453
00454 }
00455