]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliCFMuonSingleTask1.cxx
Update for AliAODHeader. (R. Arnaldi)
[u/mrichter/AliRoot.git] / PWG3 / muon / AliCFMuonSingleTask1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 //-----------------------------------------------------------------------
17 // Example of task (running locally, on AliEn and CAF),
18 // which provides standard way of calculating acceptance and efficiency
19 // between different steps of the procedure.
20 // The ouptut of the task is a AliCFContainer from which the efficiencies
21 // can be calculated
22 //-----------------------------------------------------------------------
23 // Author : R. Vernet, Consorzio Cometa - Catania (it)
24 //-----------------------------------------------------------------------
25 // Modification done by X. Lopez - LPC Clermont (fr)
26 //-----------------------------------------------------------------------
27
28 #ifndef ALICFMUONSINGLETASK1_CXX
29 #define ALICFMUONSINGLETASK1_CXX
30
31 #include "AliCFMuonSingleTask1.h"
32 #include "AliHeader.h"
33 #include "AliESDHeader.h"
34 #include "AliStack.h"
35 #include "TParticle.h"
36 #include "TLorentzVector.h"
37 #include "TH1I.h"
38 #include "AliMCEvent.h"
39 #include "AliAnalysisManager.h"
40 #include "AliESDEvent.h"
41 #include "AliAODEvent.h"
42 #include "AliCFManager.h"
43 #include "AliCFCutBase.h"
44 #include "AliCFContainer.h"
45 #include "TChain.h"
46 #include "AliESDtrack.h"
47 #include "AliLog.h"
48 #include "AliESDMuonTrack.h"
49 #include "AliESDtrack.h"
50 #include "AliESDInputHandler.h"
51 #include "TCanvas.h"
52
53 ClassImp(AliCFMuonSingleTask1)
54
55 //__________________________________________________________________________
56 AliCFMuonSingleTask1::AliCFMuonSingleTask1() :
57   fReadAODData(0),
58   fCFManager(0x0),
59   fQAHistList(0x0),
60   fHistEventsProcessed(0x0),
61   fNevt(0)
62 {
63   //
64   //Default ctor
65   //
66 }
67 //___________________________________________________________________________
68 AliCFMuonSingleTask1::AliCFMuonSingleTask1(const Char_t* name) :
69   AliAnalysisTaskSE(name),
70   fReadAODData(0),
71   fCFManager(0x0),
72   fQAHistList(0x0),
73   fHistEventsProcessed(0x0),
74   fNevt(0)
75 {
76   //
77   // Constructor. Initialization of Inputs and Outputs
78   //
79   Info("AliCFMuonSingleTask1","Calling Constructor");
80
81   fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;
82
83   DefineOutput(1,TH1I::Class());
84   DefineOutput(2,AliCFContainer::Class());
85
86 }
87
88 //___________________________________________________________________________
89 AliCFMuonSingleTask1& AliCFMuonSingleTask1::operator=(const AliCFMuonSingleTask1& c) 
90 {
91   //
92   // Assignment operator
93   //
94   if (this!=&c) {
95     AliAnalysisTaskSE::operator=(c) ;
96     fReadAODData = c.fReadAODData ;
97     fCFManager  = c.fCFManager;
98     fQAHistList = c.fQAHistList ;
99     fHistEventsProcessed = c.fHistEventsProcessed;
100     fNevt = c.fNevt ;
101   }
102   return *this;
103 }
104
105 //___________________________________________________________________________
106 AliCFMuonSingleTask1::AliCFMuonSingleTask1(const AliCFMuonSingleTask1& c) :
107   AliAnalysisTaskSE(c),
108   fReadAODData(c.fReadAODData),
109   fCFManager(c.fCFManager),
110   fQAHistList(c.fQAHistList),
111   fHistEventsProcessed(c.fHistEventsProcessed),
112   fNevt(c.fNevt)
113 {
114   //
115   // Copy Constructor
116   //
117 }
118
119 //___________________________________________________________________________
120 AliCFMuonSingleTask1::~AliCFMuonSingleTask1() {
121   //
122   //destructor
123   //
124   Info("~AliCFMuonSingleTask1","Calling Destructor");
125   if (fCFManager)           delete fCFManager ;
126   if (fHistEventsProcessed) delete fHistEventsProcessed ;
127   if (fQAHistList) {fQAHistList->Clear(); delete fQAHistList;}
128 }
129
130 //_________________________________________________
131 void AliCFMuonSingleTask1::UserExec(Option_t *)
132 {
133   //
134   // Main loop function
135   //  
136
137   Info("UserExec","") ;
138   if (!fMCEvent) {
139     Error("UserExec","NO MC EVENT FOUND!");
140     return;
141   }
142
143   fNevt++;
144   fCFManager->SetEventInfo(fMCEvent);
145   Double_t containerInput[15] ;
146
147 ////////
148 //// MC
149 ////////
150
151 // loop on the MC part
152   for (Int_t ipart=0; ipart<fMCEvent->GetNumberOfTracks(); ipart++) { 
153
154       AliMCParticle *mcPart  = (AliMCParticle*) fMCEvent->GetTrack(ipart);
155       TParticle *part = mcPart->Particle();
156
157     // Selection of mu-
158       if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
159     // rapidity and Pt cuts
160       if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mcPart)) continue;
161
162       Float_t emc = part->Energy();
163       Float_t pzmc = part->Pz();           
164       Float_t rapmc = Rap(emc,pzmc);
165       Float_t phimc = part->Phi(); 
166       phimc = Phideg(phimc);
167       Float_t ptmc = part->Pt();
168       Float_t pmc = part->P();    
169       Float_t zmc = part->Vz();
170
171       containerInput[0] = fNevt ;   
172       containerInput[1] = rapmc ;   
173       containerInput[2] = phimc ;   
174       containerInput[3] = ptmc ;   
175       containerInput[4] = pmc ;
176       containerInput[14] = zmc ;
177
178     // 9 var calculated only for ESD
179       for(Int_t i=5; i<=13; i++){
180           containerInput[i] = 1;
181       }
182
183     // fill the container at the first step
184       fCFManager->GetParticleContainer()->Fill(containerInput,0);
185
186   }
187
188 ////////
189 //// ESD
190 ////////
191
192   Int_t trig1MuL = 1 ,trig1MuH = 2, trigUSL = 4, trigUSH = 8, trigLSL = 16, trigLSH = 32;
193   Int_t trig = 0;
194  
195   AliESDEvent *fESD; 
196   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*>
197       (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
198   fESD = esdH->GetEvent();
199
200 // trigger 
201
202   if (fESD->GetTriggerMask() & trig1MuL) trig = 1;
203   if (fESD->GetTriggerMask() & trig1MuH) trig = 2;
204   if (fESD->GetTriggerMask() & trigUSL)  trig = 3;
205   if (fESD->GetTriggerMask() & trigUSH)  trig = 4;
206   if (fESD->GetTriggerMask() & trigLSL)  trig = 5;
207   if (fESD->GetTriggerMask() & trigLSH)  trig = 6;
208  
209 // vertex 
210
211   Float_t vx = -200 , vy = -200 , vz = -200;
212   
213   AliESDVertex* vertex = (AliESDVertex*) fESD->GetVertex();
214   if (vertex->GetNContributors()) {
215       vz = vertex->GetZv();
216       vy = vertex->GetYv();
217       vx = vertex->GetXv();
218   }
219
220 // tracks
221
222   Int_t mult1 = fESD->GetNumberOfMuonTracks() ;
223
224     for (Int_t j = 0; j<mult1; j++) {
225         AliESDMuonTrack* mu1 = new AliESDMuonTrack(*(fESD->GetMuonTrack(j)));
226         Float_t zr1 = mu1->Charge();
227
228 // Select mu-
229         if(zr1<0){
230             Float_t er = mu1->E();
231             Float_t pzr = mu1->Pz();
232             Float_t rapr=Rap(er,pzr);
233             Float_t phir = mu1->Phi(); 
234             phir = Phideg(phir);
235             Float_t ptr = mu1->Pt();
236             Float_t pr = mu1->P();
237             Float_t hit = mu1->GetNHit();       
238             Float_t chi2 = mu1->GetChi2() / (2.0 * hit - 5);
239             Int_t matchtrig = mu1->GetMatchTrigger();
240             Float_t chi2match = mu1->GetChi2MatchTrigger();
241             Float_t dcar = mu1->GetDCA();
242             Float_t zr = mu1->GetZ();
243
244 // rapidity and Pt cuts
245             if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mu1)) continue;
246
247             if(dcar<50){
248
249                 containerInput[0] = fNevt ;   
250                 containerInput[1] = rapr ;   
251                 containerInput[2] = phir ;   
252                 containerInput[3] = ptr ;
253                 containerInput[4] = pr ;
254                 containerInput[5] = hit ;   
255                 containerInput[6] = chi2 ;   
256                 containerInput[7] = matchtrig ;   
257                 containerInput[8] = chi2match ;
258                 containerInput[9] = vx ;   
259                 containerInput[10] = vy ;   
260                 containerInput[11] = vz ;   
261                 containerInput[12] = trig ;
262                 containerInput[13] = dcar ;
263                 containerInput[14] = zr ;
264
265 // fill the container at the second step
266                 fCFManager->GetParticleContainer()->Fill(containerInput,1);
267
268             } // Limit of histos
269             
270         }  // mu- Selection
271     }      // mu Loop
272
273 //  ----------
274   fHistEventsProcessed->Fill(0);
275   PostData(1,fHistEventsProcessed) ;
276   PostData(2,fCFManager->GetParticleContainer()) ;
277
278 }
279 //________________________________________________________________________
280 Float_t AliCFMuonSingleTask1::Rap(Float_t e, Float_t pz) 
281 {
282 // calculate rapidity
283     Float_t rap;
284     if(e!=pz){
285         rap = 0.5*TMath::Log((e+pz)/(e-pz));
286         return rap;
287     }
288     else{
289         rap = -200;
290         return rap;
291     }
292 }
293 //________________________________________________________________________
294 Float_t AliCFMuonSingleTask1::Phideg(Float_t phi) 
295 {
296 // calculate Phi in range [-180,180] 
297     Float_t phideg;
298     
299         phideg = phi-TMath::Pi();
300         phideg = phideg*57.296;
301         return phideg;
302 }
303 //________________________________________________________________________
304 void AliCFMuonSingleTask1::Terminate(Option_t *) 
305 {
306   // project pt (var[3]) from the two steps MC(0) and ESD(1)
307
308     AliCFContainer *cont = dynamic_cast<AliCFContainer*> (GetOutputData(2));
309
310     TH1D *kpt = cont->ShowProjection(3,0);
311     TH1D *rpt = cont->ShowProjection(3,1);
312
313     TCanvas *c1 = new TCanvas("AliCFMuonSingleTask1"," MC & ESD",10,10,510,510);
314     c1->Divide(1,2);
315     c1->cd(1);
316     kpt->Draw("HIST");
317     c1->cd(2);
318     rpt->Draw("HIST");
319
320 }
321 //________________________________________________________________________
322
323 #endif