]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliJetEmbeddingFromGenTask.cxx
Switch axes during filling the event counter
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliJetEmbeddingFromGenTask.cxx
1 // $Id$
2 //
3 // Jet embedding task.
4 //
5 // Author: S.Aiola, C.Loizides
6
7 #include "AliJetEmbeddingFromGenTask.h"
8
9 #include <TClonesArray.h>
10 #include <TFolder.h>
11 #include <TLorentzVector.h>
12 #include <TParticle.h>
13 #include <TParticlePDG.h>
14 #include <TRandom3.h>
15 #include <TProfile.h>
16 #include "AliAnalysisManager.h"
17 #include "AliEMCALDigit.h"
18 #include "AliEMCALGeometry.h"
19 #include "AliEMCALRecPoint.h"
20 #include "AliGenerator.h"
21 #include "AliHeader.h"
22 #include "AliLog.h"
23 #include "AliPicoTrack.h"
24 #include "AliRun.h"
25 #include "AliRunLoader.h"
26 #include "AliStack.h"
27 #include "AliStack.h"
28 #include "AliVCluster.h"
29 #include "AliVEvent.h"
30 #include "AliGenPythiaEventHeader.h"
31
32 ClassImp(AliJetEmbeddingFromGenTask)
33
34 //________________________________________________________________________
35 AliJetEmbeddingFromGenTask::AliJetEmbeddingFromGenTask() : 
36   AliJetModelBaseTask("AliJetEmbeddingFromGenTask"),
37   fGen(0),
38   fMassless(kFALSE),
39   fChargedOnly(kFALSE),
40   fHistPt(0),
41   fHistEtaPhi(0),
42   fHistTrials(0),
43   fHistXsection(0),
44   fHistPtHard(0)
45 {
46   // Default constructor.
47   SetSuffix("EmbeddedFromGen");
48 }
49
50 //________________________________________________________________________
51 AliJetEmbeddingFromGenTask::AliJetEmbeddingFromGenTask(const char *name, Bool_t drawqa) :
52   AliJetModelBaseTask(name,drawqa),
53   fGen(0),
54   fMassless(kFALSE),
55   fChargedOnly(kFALSE),
56   fHistPt(0),
57   fHistEtaPhi(0),
58   fHistTrials(0),
59   fHistXsection(0),
60   fHistPtHard(0)
61 {
62   // Standard constructor.
63   SetSuffix("EmbeddedFromGen");
64 }
65
66 //________________________________________________________________________
67 AliJetEmbeddingFromGenTask::~AliJetEmbeddingFromGenTask()
68 {
69   // Destructor
70 }
71
72 //________________________________________________________________________
73 void AliJetEmbeddingFromGenTask::UserCreateOutputObjects()
74 {
75   // Create user output.
76
77   if (!fQAhistos)
78     return;
79
80   AliJetModelBaseTask::UserCreateOutputObjects();
81
82   fHistPt = new TH1F("fHistpt","fHistPt;#it{p}_{T};N",100,0.,100.);
83   fOutput->Add(fHistPt);
84
85   fHistEtaPhi = new TH2F("fHistEtapHI","fHistEtaPhi;#eta;#varphi",100,-3.,3.,100.,0.,TMath::TwoPi());
86   fOutput->Add(fHistEtaPhi);
87
88   fHistTrials = new TH1F("fHistTrials", "fHistTrials", 1, 0, 1);
89   fHistTrials->GetYaxis()->SetTitle("trials");
90   fOutput->Add(fHistTrials);
91
92   fHistXsection = new TProfile("fHistXsection", "fHistXsection", 1, 0, 1);
93   fHistXsection->GetYaxis()->SetTitle("xsection");
94   fOutput->Add(fHistXsection);
95
96   fHistPtHard = new TH1F("fHistPtHard", "fHistPtHard", 500, 0., 500.);
97   fHistPtHard->GetXaxis()->SetTitle("p_{T,hard} (GeV/c)");
98   fHistPtHard->GetYaxis()->SetTitle("counts");
99   fOutput->Add(fHistPtHard);
100
101   PostData(1, fOutput);
102 }
103
104 //________________________________________________________________________
105 Bool_t AliJetEmbeddingFromGenTask::ExecOnce() 
106 {
107   // Exec only once.
108
109   if (!gAlice) {
110     new AliRun("gAlice","The ALICE Off-line Simulation Framework");
111     delete gRandom;
112     gRandom = new TRandom3(0);
113   }
114
115   TFolder *folder = new TFolder(GetName(),GetName());
116   AliRunLoader *rl = new AliRunLoader(folder);
117   rl->MakeHeader();
118   rl->MakeStack();
119   AliStack *stack = rl->Stack();
120   fGen->SetStack(stack);
121   fGen->Init();
122
123   if (!(InputEvent()->FindListObject(fTracksName))) {
124     fOutTracks = new TClonesArray("AliPicoTrack", 1000);
125     fOutTracks->SetName(fTracksName);
126     InputEvent()->AddObject(fOutTracks);
127     fNTracks = 0;
128   }
129   
130   return kTRUE;
131 }
132
133 //________________________________________________________________________
134 void AliJetEmbeddingFromGenTask::Run() 
135 {
136   // Embed particles.
137
138   if (fCopyArray) 
139     CopyTracks();
140
141   AliStack *stack = fGen->GetStack();
142   stack->Reset();
143   fGen->Generate();
144   const Int_t nprim = stack->GetNprimary();
145   for (Int_t i=0;i<nprim;++i) {
146     if (!stack->IsPhysicalPrimary(i))
147       continue;
148     TParticle *part = stack->Particle(i);
149     TParticlePDG *pdg = part->GetPDG(1);
150     if (!pdg) 
151       continue;
152     Int_t c = (Int_t)(TMath::Abs(pdg->Charge()));
153     if (fChargedOnly && c==0)  continue;
154     Double_t pt = part->Pt();
155     Double_t eta = part->Eta();
156     Double_t phi = part->Phi();
157     if (eta<fEtaMin)
158       continue;
159     if (eta>fEtaMax)
160       continue;
161     if (phi<fPhiMin)
162       continue;
163     if (phi>fPhiMax)
164       continue;
165     if (pt<fPtMin)
166       continue;
167     if (pt>fPtMax)
168       continue;
169     Double_t mass = part->GetMass();
170     if(fMassless) mass = 0.;
171     fHistPt->Fill(pt);
172     fHistEtaPhi->Fill(eta,phi);
173     AddTrack(pt, eta, phi,0,0,0,0,0,0,c,mass);
174   }
175
176   FillPythiaHistograms();
177 }
178
179 //________________________________________________________________________
180 void AliJetEmbeddingFromGenTask::FillPythiaHistograms() {
181   //Get PYTHIA info: pt-hard, x-section, trials
182
183   if (!fQAhistos)
184     return;
185
186   AliRunLoader *rl = AliRunLoader::Instance();
187   AliGenPythiaEventHeader *genPH = dynamic_cast<AliGenPythiaEventHeader*>(rl->GetHeader()->GenEventHeader());
188   if(genPH) {
189     Float_t xsec = genPH->GetXsection();
190     Int_t trials = genPH->Trials();
191     Float_t pthard = genPH->GetPtHard();
192
193     fHistXsection->Fill(0.5,xsec);
194     fHistTrials->Fill(0.5,trials);
195     fHistPtHard->Fill(pthard);
196   }
197 }