]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDAnalysisTask.cxx
New task to check TPC-ITS track prolongation eff with cosmics
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDAnalysisTask.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 #ifndef AliHMPIDAnalysisTASK_CXX
17 #define AliHMPIDAnalysisTASK_CXX
18
19 #include "AliHMPIDAnalysisTask.h"
20 #include "TCanvas.h"
21 #include "AliStack.h"
22 #include "TParticle.h"
23 #include "Riostream.h"
24 #include "TH1I.h"
25 #include "TH2F.h"
26 #include "TH1F.h"
27 #include "AliMCEvent.h"
28 #include "AliAnalysisManager.h"
29 #include "AliESDEvent.h"
30 #include "TChain.h"
31 #include "AliESDtrack.h"
32 #include "AliLog.h"
33
34 ClassImp(AliHMPIDAnalysisTask)
35
36 //__________________________________________________________________________
37 AliHMPIDAnalysisTask::AliHMPIDAnalysisTask() :
38   fHistList(0x0),
39   fHistEventsProcessed(0x0),
40   fNevts(0),
41   fTrigNevts(0),
42   fTrigger(0)
43 {
44   //
45   //Default ctor
46   //
47 }
48 //___________________________________________________________________________
49 AliHMPIDAnalysisTask::AliHMPIDAnalysisTask(const Char_t* name) :
50   AliAnalysisTaskSE(name),
51   fHistList(0x0),
52   fHistEventsProcessed(0x0),
53   fNevts(0),
54   fTrigNevts(0),
55   fTrigger(0)
56 {
57   //
58   // Constructor. Initialization of Inputs and Outputs
59   //
60   Info("AliHMPIDAnalysisTask","Calling Constructor");
61
62   /*
63     DefineInput(0) and DefineOutput(0)
64     are taken care of by AliAnalysisTaskSE constructor
65   */
66   DefineOutput(1,TH1I::Class());
67   DefineOutput(2,TList::Class());
68 }
69
70 //___________________________________________________________________________
71 AliHMPIDAnalysisTask& AliHMPIDAnalysisTask::operator=(const AliHMPIDAnalysisTask& c) 
72 {
73   //
74   // Assignment operator
75   //
76   if (this!=&c) {
77     AliAnalysisTaskSE::operator=(c) ;
78     fHistList = c.fHistList ;
79     fHistEventsProcessed = c.fHistEventsProcessed;
80     fNevts = c.fNevts;
81     fTrigNevts = c.fTrigNevts;
82     fTrigger = c.fTrigger;
83   }
84   return *this;
85 }
86
87 //___________________________________________________________________________
88 AliHMPIDAnalysisTask::AliHMPIDAnalysisTask(const AliHMPIDAnalysisTask& c) :
89   AliAnalysisTaskSE(c),
90   fHistList(c.fHistList),
91   fHistEventsProcessed(c.fHistEventsProcessed),
92   fNevts(c.fNevts),
93   fTrigNevts(c.fTrigNevts),
94   fTrigger(c.fTrigger)
95 {
96   //
97   // Copy Constructor
98   //
99 }
100
101 //___________________________________________________________________________
102 AliHMPIDAnalysisTask::~AliHMPIDAnalysisTask() {
103   //
104   //destructor
105   //
106   Info("~AliHMPIDAnalysisask","Calling Destructor");
107   if (fHistEventsProcessed) delete fHistEventsProcessed ;
108   if (fHistList) {fHistList->Clear(); delete fHistList;}
109 }
110
111 //_________________________________________________
112 void AliHMPIDAnalysisTask::UserExec(Option_t *)
113 {
114   //
115   // Main loop function, executed on Event basis
116   //
117   Info("UserExec","") ;
118
119   AliESDEvent* fESD = dynamic_cast<AliESDEvent*>(fInputEvent);
120   if (!fESD) {
121     Error("UserExec","NO ESD FOUND!");
122     return;
123   }
124   fNevts++;
125  
126   if(fESD->GetTriggerMask()&fTrigger == fTrigger) fTrigNevts++;
127
128
129   for (Int_t iTrack = 0; iTrack<fESD->GetNumberOfTracks(); iTrack++) {
130
131     AliESDtrack* track = fESD->GetTrack(iTrack);
132     if(!track) continue;
133     Double_t rin[3], rout[3];
134     track->GetInnerXYZ(rin);
135     track->GetOuterXYZ(rout);
136   
137      ((TH2F *)fHistList->At(22))->Fill(rin[0],rin[1]);
138      ((TH2F *)fHistList->At(23))->Fill(rout[0],rout[1]);
139
140
141      TH1F *h = ((TH1F*)(fHistList->At(0)));
142      if(track->GetHMPIDsignal() == -20) h->Fill(0);
143      else if(track->GetHMPIDsignal() == -9) h->Fill(1);
144      else if(track->GetHMPIDsignal() == -5) h->Fill(2);
145      else if(track->GetHMPIDsignal() == -11) h->Fill(3);
146      else if(track->GetHMPIDsignal() == -22) h->Fill(4);
147      else h->Fill(4);
148     
149      if(track->GetHMPIDcluIdx() == 99099999) continue;
150      if(track->GetHMPIDcluIdx() < 0) continue;
151
152     Int_t ch = track->GetHMPIDcluIdx()/1000000; 
153
154     Float_t x, y; Int_t q, nph; 
155     track->GetHMPIDmip(x,y,q,nph);
156
157     if(x==0 && y==0) continue;
158     Float_t xpc, ypc, th, ph;
159     track->GetHMPIDtrk(xpc,ypc,th,ph); //special setting in local cosmic rec!!!
160                                        //do not use with standard rec
161
162     Double_t dist = TMath::Sqrt( (xpc-x)*(xpc-x) + (ypc - y)*(ypc - y));
163
164      if(ch >=0 && ch < 7) {
165      if(dist < 3) ((TH2F *)fHistList->At(ch+1))->Fill(x,y);
166      ((TH1F *)fHistList->At(ch+7+1))->Fill(q);
167      ((TH2F *)fHistList->At(ch+14+1))->Fill(dist,q);
168      }
169   
170     }//track loop
171
172   fHistEventsProcessed->Fill(0);
173   
174   /* PostData(0) is taken care of by AliAnalysisTaskSE */
175   PostData(1,fHistEventsProcessed) ;
176   PostData(2,fHistList) ;
177 }
178
179
180 //___________________________________________________________________________
181 void AliHMPIDAnalysisTask::Terminate(Option_t*)
182 {
183   // The Terminate() function is the last function to be called during
184   // a query. It always runs on the client, it can be used to present
185   // the results graphically or save the results to file.
186
187   Info("Terminate","");
188   AliAnalysisTaskSE::Terminate();
189
190 }
191
192
193 //___________________________________________________________________________
194 void AliHMPIDAnalysisTask::UserCreateOutputObjects() {
195   //
196   //HERE ONE CAN CREATE OUTPUT OBJECTS
197   //TO BE SET BEFORE THE EXECUTION OF THE TASK
198   //
199   Info("CreateOutputObjects","CreateOutputObjects of task %s", GetName());
200
201   //slot #1
202   OpenFile(1);
203   fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;
204
205    OpenFile(2);
206    fHistList = new TList();
207    
208    //0
209    TH1F *trkH = new TH1F("trkH","signal flags in HMPID",6,0,6);
210    TString summary[6] =  {"NotPerformed","MipDistCut", "MipQdcCut", "NoPhotAccept", "kNoRad", "other"};
211    for(Int_t ibin = 0; ibin < 6; ibin++) trkH->GetXaxis()->SetBinLabel(ibin+1,Form("%i  %s",ibin+1,summary[ibin].Data()));
212    fHistList->Add(trkH);
213
214   
215    TH2F *mod[7], *dq[7];
216    TH1F *q[7];
217
218    //1-7
219    for(Int_t i=0; i< 7 ; i++) {// to preserve the histo sorting
220    mod[i] = new TH2F(Form("mod%i",i),Form("MIP position in chamber %i",i),180,0,180,180,0,180);
221    mod[i]->SetMarkerStyle(8);
222    mod[i]->SetMarkerColor(2);
223    fHistList->Add(mod[i]);
224    }
225    //8-14
226    for(Int_t i=0; i< 7 ; i++) {//to reserve the histo sorting
227    q[i] = new TH1F(Form("q%i",i),Form("MIP charge distribution in chamber %i",i),5000,0,5000);
228    q[i]->SetMarkerColor(2);
229    fHistList->Add(q[i]);
230    }
231    //15-21
232    for(Int_t i=0; i< 7 ; i++) {//to reserve the histo sorting
233    dq[i] = new TH2F(Form("dq%i",i),Form("#Delta(mip-track) vs Q_{mip} in chamber %i",i),1000,0,100,5000,0,5000),
234    dq[i]->SetMarkerStyle(6);
235    fHistList->Add(dq[i]);
236    }
237    //22
238    TH2F *inner = new TH2F("inner","inner track XY",800,-400,400,800,-400,400);
239    inner->SetMarkerStyle(6);
240    inner->SetXTitle("X cm");
241    inner->SetYTitle("Y cm");
242    fHistList->Add(inner); 
243    //23
244    TH2F *outer = new TH2F("outer","outer track XY",800,-400,400,800,-400,400);
245    outer->SetMarkerStyle(6);
246    outer->SetXTitle("X cm");
247    outer->SetYTitle("Y cm");
248    fHistList->Add(outer);
249
250 }
251
252 #endif