]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FASTSIM/FASTSmear.C
Fixes for not filled histograms and calculation of Dijet binning
[u/mrichter/AliRoot.git] / FASTSIM / FASTSmear.C
CommitLineData
797082f2 1TH1F *hpt;
2void FASTSmear(char *finname="galice_upsi.root",char *foutname="FAST_upsi.root",Float_t bkg=0,Float_t nevmax=100000000){
3
4 printf ("processing file %s\n",finname);
5 if (gClassTable->GetID("AliRun") < 0) {
6 gROOT->LoadMacro("loadlibs.C");
7 loadlibs();
8 }
9 // create the structure to hold the variables for the branch
10
11 struct fast_t {
12 Float_t accp;
13 Float_t effp;
14 Float_t trigeffp;
15 Float_t pgenp;
16 Float_t psmearp;
17 Float_t thetagenp;
18 Float_t thetasmearp;
19 Float_t phigenp;
20 Float_t phismearp;
21 Float_t accm;
22 Float_t effm;
23 Float_t trigeffm;
24 Float_t pgenm;
25 Float_t psmearm;
26 Float_t thetagenm;
27 Float_t thetasmearm;
28 Float_t phigenm;
29 Float_t phismearm;
30 };
31
32 fast_t fast;
33
34 // open the input file
35 TFile *fin=new TFile(finname);
36 gAlice=(AliRun*)fin->Get("gAlice");
37
38 // create the output file, its tree and branch
39 TFile *fout = new TFile(foutname,"RECREATE");
40 TTree *FASTtrack = new TTree("FASTtrack","mu tracks smeared with fastsim");
41 FASTtrack->Branch("fast",&fast.accp,"accp:effp:trigeffp:pgenp:psmearp:thetagenp:thetasmearp:phigenp:phismearp:accm:effm:trigeffm:pgenm:psmearm:thetagenm:thetasmearm:phigenm:phismearm");
42
43 // create the fast tracker
44
45 AliFastMuonTriggerEff *trigeff = new AliFastMuonTriggerEff();
46 AliFastMuonTrackingAcc *acc = new AliFastMuonTrackingAcc();
47 AliFastMuonTrackingEff *eff = new AliFastMuonTrackingEff();
48 AliFastMuonTrackingRes *res = new AliFastMuonTrackingRes();
49 acc->SetBackground(bkg);
50 eff->SetBackground(bkg);
51 res->SetBackground(bkg);
52 acc->Init();
53 eff->Init();
54 res->Init();
55 AliFastDetector* tracker = new AliFastDetector("Tracker", "Muon Tracker");
56 tracker->AddResponse(trigeff);
57 tracker->AddResponse(acc);
58 tracker->AddResponse(eff);
59 tracker->AddResponse(res);
60 tracker->Init();
61 tracker->Dump();
62
63 // loop over the events
64
880d6abe 65 Int_t nev=AliRunLoader::GetNumberOfEvents();
797082f2 66 TParticle *mup, *mum;
67
68 if (nev>nevmax) nev = nevmax;
69
70 table = AliMUONFastTracking::Instance();
71 printf ("background set to %g \n",table->GetBackground());
72 for (Int_t iev=0; iev<nev; iev++) {
73 Int_t npart=gAlice->GetEvent(iev);
74 // npart=500;
75 for (Int_t ipart=0; ipart<npart; ipart+=3) {
76 if (!(ipart%30)) printf ("Event #%d upsilon #%d \n",iev,ipart/3);
77 part = gAlice->Particle(ipart+1);
78 if (part->GetPdgCode() == 13) mum = part;
79 else if (part->GetPdgCode() == -13) mup = part;
80 part = gAlice->Particle(ipart+2);
81 if (part->GetPdgCode() == 13) mum = part;
82 else if (part->GetPdgCode() == -13) mup = part;
83
84 // the mu+
85 printf ("background set to %g \n",table->GetBackground());
86 res->SetCharge(1);
87 eff->SetCharge(1);
88 acc->SetCharge(1);
89 fast.pgenp = mup->P();
90 Double_t ptp = fast.pgenp * TMath::Sin(mup->Theta());
91 fast.thetagenp = 180./TMath::Pi() * mup->Theta();
92 fast.phigenp = 180./TMath::Pi() * mup->Phi();
93 if (fast.phigenp>180) fast.phigenp -=360;
94 res->Evaluate(fast.pgenp,fast.thetagenp,fast.phigenp,
95 fast.psmearp,fast.thetasmearp,fast.phismearp);
96 fast.effp = eff->Evaluate(ptp,fast.thetagenp,fast.phigenp);
97 fast.accp = acc->Evaluate(ptp,fast.thetagenp,fast.phigenp);
98 fast.trigeffp = trigeff->Evaluate(1,ptp,fast.thetagenp,fast.phigenp);
99
100 // the mu-
101 res->SetCharge(-1);
102 acc->SetCharge(-1);
103 eff->SetCharge(-1);
104 fast.pgenm = mum->P();
105 Double_t ptm = fast.pgenm * TMath::Sin(mum->Theta());
106 fast.thetagenm = 180./TMath::Pi() * mum->Theta();
107 fast.phigenm = 180./TMath::Pi() * mum->Phi();
108 if (fast.phigenm>180) fast.phigenm -=360;
109 res->Evaluate(fast.pgenm,fast.thetagenm,fast.phigenm,
110 fast.psmearm,fast.thetasmearm,fast.phismearm);
111 fast.effm = eff->Evaluate(ptm,fast.thetagenm,fast.phigenm);
112 fast.accm = acc->Evaluate(ptm,fast.thetagenm,fast.phigenm);
113 fast.trigeffm = trigeff->Evaluate(-1,ptm,fast.thetagenm,fast.phigenm);
114
115 // fill the tree
116 FASTtrack->Fill();
117 }
118 }
119 fout->Write();
120 fout->Close();
121 fin->Close();
122}
123
124
125