]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskPythiaMpi.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskPythiaMpi.cxx
1 // **************************************
2 // Task used for Pythia MPI studies in ZZ train 
3 // Output is stored in an exchance container
4 // *******************************************
5
6
7 /**************************************************************************
8  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
9  *                                                                        *
10  * Author: The ALICE Off-line Project.                                    *
11  * Contributors are mentioned in the code where appropriate.              *
12  *                                                                        *
13  * Permission to use, copy, modify and distribute this software and its   *
14  * documentation strictly for non-commercial purposes is hereby granted   *
15  * without fee, provided that the above copyright notice appears in all   *
16  * copies and that both the copyright notice and this permission notice   *
17  * appear in the supporting documentation. The authors make no claims     *
18  * about the suitability of this software for any purpose. It is          *
19  * provided "as is" without express or implied warranty.                  *
20  **************************************************************************/
21
22 #include "TChain.h"
23
24 #include "AliAnalysisDataSlot.h"
25 #include "AliAnalysisDataContainer.h"
26 #include "AliAnalysisManager.h"
27 #include "AliMCEvent.h"
28 #include "AliVEvent.h"
29 #include "AliAnalysisTaskPythiaMpi.h"
30 #include "AliGenEventHeader.h"
31 #include "AliVParticle.h"
32 #include "AliPythia.h"
33 #include "AliPythia8.h"
34
35 using namespace std;
36
37 ClassImp(AliAnalysisTaskPythiaMpi)
38
39 //_________________________________________________________| Constructor
40
41 AliAnalysisTaskPythiaMpi::AliAnalysisTaskPythiaMpi(): 
42   fMcEvent(0x0),  
43   fMcHandler(0x0),
44   fOutputList(0),
45   fHistEvents(0),
46   fHistPt(0),
47   fHistEta(0),
48   fHistMpi(0),
49   fHistMultMpi(0),
50   fHistdNdetaMpi(0)
51 {
52   //
53   // Default constructor
54   //
55 }
56
57 ////_________________________________________________________| Specific Constructor
58
59 AliAnalysisTaskPythiaMpi::AliAnalysisTaskPythiaMpi(const char* name):
60   AliAnalysisTaskSE(name),
61   fMcEvent(0x0),  
62   fMcHandler(0x0),
63   fOutputList(0),
64   fHistEvents(0),
65   fHistPt(0),
66   fHistEta(0),
67   fHistMpi(0),
68   fHistMultMpi(0),
69   fHistdNdetaMpi(0)
70 {
71   //
72   // Specific constructor. Initialization of Inputs and Outputs
73   //
74   Info("AliAnalysisMpi","Calling Constructor");
75   // Output slot #1 writes into a TList container 
76   DefineInput(0, TChain::Class());
77   DefineOutput(1, TList::Class());  
78 }
79
80 //_________________________________________________________| Assignment Operator
81
82 AliAnalysisTaskPythiaMpi& AliAnalysisTaskPythiaMpi::operator=(const AliAnalysisTaskPythiaMpi& c)
83 {
84   // Assigment operator
85   if(this!=&c) {
86     AliAnalysisTaskSE::operator=(c);
87     fMcEvent = c.fMcEvent;
88     fMcHandler = c.fMcHandler;
89     fOutputList = c.fOutputList;
90     fHistEvents = c.fHistEvents;
91     fHistPt = c.fHistPt;
92     fHistEta = c.fHistEta;
93     fHistMpi = c.fHistMpi;
94     fHistMultMpi = c.fHistMultMpi;
95     fHistdNdetaMpi = c.fHistdNdetaMpi;
96
97   }
98   return *this;
99
100 }
101 //_________________________________________________________| Copy Constructor
102
103 AliAnalysisTaskPythiaMpi::AliAnalysisTaskPythiaMpi(const AliAnalysisTaskPythiaMpi& c):
104   AliAnalysisTaskSE(c),
105   fMcEvent(c.fMcEvent),
106   fMcHandler(c.fMcHandler), 
107   fOutputList(c.fOutputList),
108   fHistEvents(c.fHistEvents),
109   fHistPt(c.fHistPt),
110   fHistEta(c.fHistEta),
111   fHistMpi(c.fHistMpi),
112   fHistMultMpi(c.fHistMultMpi),
113   fHistdNdetaMpi(c.fHistdNdetaMpi)
114 {
115   //Copy Constructor
116 }
117
118 //_________________________________________________________| Destructor
119
120 AliAnalysisTaskPythiaMpi::~AliAnalysisTaskPythiaMpi()
121 {
122   //
123   // Destructor
124   //
125   Info("~AliAnalysisTaskPythiaMpi","Calling Destructor");
126   if(fOutputList) delete fOutputList;
127   if(fHistEvents) delete fHistEvents;
128   if(fHistPt) delete fHistPt;
129   if(fHistEta) delete fHistEta;
130   if(fHistMpi) delete fHistMpi;
131   if(fHistMultMpi) delete fHistMultMpi;
132   if(fHistdNdetaMpi) delete fHistdNdetaMpi;
133   
134 }
135
136 //_____________________________________________________________________
137
138 void AliAnalysisTaskPythiaMpi::UserCreateOutputObjects()
139 {
140   Info("CreateOutputObjects","CreateOutputObjects of task %s",GetName());
141   
142   fOutputList = new TList();
143   fOutputList->SetOwner(kTRUE);
144  
145   fHistEvents = new TH1I("fHistNEvents","fHistNEvents",2,0,2);
146   fHistEvents->GetXaxis()->SetBinLabel(1,"All events");
147   fHistEvents->GetXaxis()->SetBinLabel(2,"Events with |zVtx|<10cm");
148
149   fHistPt = new TH1F("fHistPt","p_{T} distribution",15,0.2,5.0);
150   fHistPt->GetXaxis()->SetTitle("p_{T} (GeV/c)");
151   fHistPt->GetYaxis()->SetTitle("dN/dp_{T} (GeV/c^{-1})");
152   fHistPt->SetMarkerStyle(kFullCircle);
153
154   fHistEta = new TH1F("fHistEta","eta distribution",20,-1,1);
155   fHistEta->GetXaxis()->SetTitle("#eta");
156   fHistEta->GetYaxis()->SetTitle("dN/d#eta");
157   fHistEta->SetMarkerStyle(kFullCircle);
158
159   fHistMpi = new TH1F("fHistMpi","MPIs distribution",50,0,50);
160   fHistMpi->GetXaxis()->SetTitle("#MPIs");
161   fHistMpi->GetYaxis()->SetTitle("dN/dMPIs");
162   fHistMpi->SetMarkerStyle(kFullCircle);
163
164   fHistMultMpi = new TH2F("fHistMultMpi","Multiplicity vs MPIs",50,0,50,100,0,1e4);
165   fHistMultMpi->GetXaxis()->SetTitle("#MPIs");
166   fHistMultMpi->GetYaxis()->SetTitle("MC particles");
167   fHistMultMpi->SetMarkerStyle(kFullDotSmall);
168
169   fHistdNdetaMpi = new TH2F("fHistdNdetaMpi","dN/d#etadMPIs",50,0,50,100,-5,5);
170   fHistdNdetaMpi->GetXaxis()->SetTitle("#MPIs");
171   fHistdNdetaMpi->GetYaxis()->SetTitle("#eta");
172   //fHistdNdetaMpi->SetMarkerStyle();
173
174
175   fOutputList->Add(fHistEvents);
176   fOutputList->Add(fHistPt);
177   fOutputList->Add(fHistEta);
178   fOutputList->Add(fHistMpi);
179   fOutputList->Add(fHistMultMpi);
180   fOutputList->Add(fHistdNdetaMpi);
181
182   PostData(1,fOutputList);
183
184   return;
185 }
186 //_________________________________________________________________________
187
188 void AliAnalysisTaskPythiaMpi::Init()
189 {
190    // MC handler
191    if(fDebug > 1) printf("AnalysisTaskPythiaMpi::Init() \n");
192    fMcHandler = dynamic_cast<AliInputEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()); 
193 }
194
195 //_________________________________________________________________________
196
197 void AliAnalysisTaskPythiaMpi::LocalInit()
198 {
199    // MC handler
200    //
201    // Initialization
202    //
203
204    if(fDebug > 1) printf("AnalysisTaskPythiaMpi::LocalInit() \n");
205    Init();
206 }
207
208 //____________________________________________| User Exec
209
210 void AliAnalysisTaskPythiaMpi::UserExec(Option_t*)
211 {
212
213    // handle and reset the output jet branch 
214   Info("UserExec","Start of method");
215    if(fDebug > 1) printf("AliAnalysisTaskPythiaMpi::UserExec \n");
216
217    //
218    // Execute analysis for current event
219    //
220    Init();
221    fHistEvents->Fill(0.5);
222    if(fMcHandler){
223       fMcEvent = fMcHandler->MCEvent(); 
224    }else{
225       if(fDebug > 1) printf("AnalysisTaskPythiaMpi::Handler() fMcHandler=NULL\n");
226       return;
227    }
228    if(!fMcEvent){
229       if(fDebug > 1) printf("AnalysisTaskPythiaMpi::Exec()   fMcEvent=NULL \n");
230       return;
231    }
232
233    const AliVVertex *vtxMC = fMcEvent->GetPrimaryVertex();
234    Float_t zVtx = vtxMC->GetZ();
235    if(TMath::Abs(zVtx)<10) fHistEvents->Fill(1.5);
236    else return;
237
238    Printf("MC particles: %d",fMCEvent->GetNumberOfTracks());
239
240
241    //MPIs
242    Double_t Nmpi = (AliPythia8::Instance())->GetNMPI();
243    fHistMpi->Fill(Nmpi);
244
245    fHistMultMpi->Fill(Nmpi,fMCEvent->GetNumberOfTracks());
246
247    //loop over the tracks 
248
249    for(Int_t iTracks = 0; iTracks < fMCEvent->GetNumberOfTracks(); iTracks++){
250      AliVParticle* track = fMCEvent->GetTrack(iTracks);
251      if(!track){
252        Printf("ERROR: Could not receive track %d (mc loop)",iTracks);
253        continue;
254      }
255      fHistPt->Fill(track->Pt());
256      fHistEta->Fill(track->Eta());
257
258      fHistdNdetaMpi->Fill(Nmpi,track->Eta());
259    }
260    
261
262    PostData(1, fOutputList);
263    return; 
264    
265 }
266 //____________________________________________________________________________
267
268 void AliAnalysisTaskPythiaMpi::Terminate(Option_t*){
269   //
270   // Terminate analysis
271   //
272   Info("Terminate","Start and end of Method");
273   AliAnalysisTaskSE::Terminate();
274
275   fOutputList = dynamic_cast<TList*>(GetOutputData(1));
276   if (!fOutputList) {
277     printf("ERROR: fOutputList not available\n");
278     return;
279   }   
280     
281 }
282