]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDOfflineCalibTask.cxx
Separate trees for ITSsa and global tracks
[u/mrichter/AliRoot.git] / PMD / AliPMDOfflineCalibTask.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  * AliPMD offline calibration task
18  *
19  * Satyajit Jena, IIT Bombay
20  * sjena@cern.ch
21  * Fri Feb 12 13:30:19 IST 2010
22  *
23  ************************************/
24
25 #include "TChain.h"
26 #include "TList.h"
27 #include "TFile.h"
28 #include "TTree.h"
29 #include "TH1F.h"
30 #include "TCanvas.h"
31 #include "AliAnalysisTask.h"
32 #include "AliAnalysisManager.h"
33 #include "AliVEvent.h"
34 #include "AliESD.h"
35 #include "AliESDEvent.h"
36 #include "AliAODEvent.h"
37 #include "TObjArray.h"
38
39 #include "AliPMDOfflineCalibTask.h"
40
41 ClassImp(AliPMDOfflineCalibTask)
42
43 //________________________________________________________________________
44 AliPMDOfflineCalibTask::AliPMDOfflineCalibTask(const char *name) 
45 : AliAnalysisTaskSE(name),
46   fListOfHistos(0),
47   fPmdCalibAdcP(0),
48   fPmdCalibAdcC(0),
49   fPmdCalibEntP(0),
50   fPmdCalibEntC(0),
51   fNEvents(0),
52   fSelectedTrigger(new TObjArray()),
53   fRejected(kTRUE)
54 {
55   // Constructor
56   
57   DefineOutput(1, TList::Class());
58 }
59
60 //________________________________________________________________________
61 void AliPMDOfflineCalibTask::UserCreateOutputObjects()
62 {
63   fListOfHistos = new TList();
64
65   fNEvents = new TH1I("hNEvents","Events Statistic", 3, 0, 5);
66
67   fPmdCalibAdcP = new TH1F("fPmdCalibAdcP", "PMD Calibration ADC fill for PRE", 234795, 0, 234795);
68   fPmdCalibAdcP->GetYaxis()->SetTitle("ADC");
69   fPmdCalibAdcP->GetXaxis()->SetTitle("Cells in SMN-ROW-COL");
70
71   fPmdCalibAdcC = new TH1F("fPmdCalibAdcC", "PMD Calibration ADC fill for CPV", 234795, 0, 234795);
72   fPmdCalibAdcC->GetYaxis()->SetTitle("ADC");
73   fPmdCalibAdcC->GetXaxis()->SetTitle("Cells in SMN-ROW-COL");
74
75   fPmdCalibEntP = new TH1F("fPmdCalibEntP", "PMD Calibration Entries fill for PRE", 234795, 0, 234795);
76   fPmdCalibEntP->GetYaxis()->SetTitle("Frequescy");
77   fPmdCalibEntP->GetXaxis()->SetTitle("Cells in SMN-ROW-COL");
78
79   fPmdCalibEntC = new TH1F("fPmdCalibEntC", "PMD Calibration Entries fill for CPV", 234795, 0, 234795);
80   fPmdCalibEntC->GetYaxis()->SetTitle("Frequescy ");
81   fPmdCalibEntC->GetXaxis()->SetTitle("Cells in SMN-ROW-COL");
82  
83   fListOfHistos->Add(fPmdCalibAdcP);
84   fListOfHistos->Add(fPmdCalibAdcC);
85   fListOfHistos->Add(fPmdCalibEntP);
86   fListOfHistos->Add(fPmdCalibEntC);
87   fListOfHistos->Add(fNEvents);
88 }
89
90 //________________________________________________________________________
91 void AliPMDOfflineCalibTask::UserExec(Option_t *) 
92 {
93
94   AliVEvent *event = InputEvent();
95   if (!event) {
96      Printf("ERROR: Could not retrieve event");
97      return;
98   }
99
100   fNEvents->Fill(1);
101
102   AliESDEvent* pmdesd = dynamic_cast<AliESDEvent*>(event);
103   if(!pmdesd) return; 
104   fNEvents->Fill(2);
105
106   Bool_t pass = kTRUE;
107   Int_t numberOfTriggerSelected = fSelectedTrigger->GetEntriesFast();
108   
109   if(fRejected) 
110     {
111       pass = kTRUE;
112       for(Int_t k = 0; k < numberOfTriggerSelected; k++)
113         {
114           const TObjString *const obString = (TObjString*)fSelectedTrigger->At(k);
115         
116           const TString tString = obString->GetString();
117           if(pmdesd->IsTriggerClassFired((const char*)tString)) 
118             {
119               pass = kFALSE;
120             }
121         }
122     }
123   else 
124     {
125       pass = kFALSE;
126       for(Int_t k = 0; k < numberOfTriggerSelected; k++)
127         {
128           const TObjString *const obString=(TObjString*)fSelectedTrigger->At(k);
129           
130           const TString tString = obString->GetString();
131           if(pmdesd->IsTriggerClassFired((const char*)tString)) 
132             {
133               pass = kTRUE;
134             }
135         }
136     }
137   
138   if(!pass) 
139     {
140       PostData(1, fListOfHistos);
141       return;
142     }
143
144   fNEvents->Fill(3);
145   
146   Int_t npmdcl = pmdesd->GetNumberOfPmdTracks();
147   if(npmdcl < 1) fNEvents->Fill(4);
148
149   while (npmdcl--)
150     {
151       AliESDPmdTrack *pmdtr = pmdesd->GetPmdTrack(npmdcl);
152       
153       Int_t   det     = pmdtr->GetDetector();
154       Int_t   smn     = pmdtr->GetSmn();
155       Float_t adc     = pmdtr->GetClusterADC();
156       Float_t sTag    = pmdtr->GetClusterSigmaX();
157       Float_t sRowCol = pmdtr->GetClusterSigmaY();
158       
159       
160       Float_t rc      = smn*10000 + sRowCol;
161       if(adc > 1200.) continue;
162       
163       if(sTag > 999. && sTag < 1000.) 
164         {
165           if(det == 0) 
166             {
167               fPmdCalibAdcP->Fill(rc,adc);
168               fPmdCalibEntP->Fill(rc);
169             }
170           else if(det == 1) 
171             {
172               fPmdCalibAdcC->Fill(rc,adc);
173               fPmdCalibEntC->Fill(rc);
174             }
175         }
176     }
177   
178
179   PostData(1, fListOfHistos);
180 }      
181
182 //________________________________________________________________________
183 void AliPMDOfflineCalibTask::Terminate(Option_t *) 
184 {
185   fListOfHistos = dynamic_cast<TList*>(GetOutputData(1));
186     
187     fPmdCalibAdcP = dynamic_cast<TH1F*>(fListOfHistos->At(0));
188     if(!fPmdCalibAdcP) {
189       printf("ERROR: No ADC File Generated for PMD-PRE ");
190       return;
191     }
192     
193     fPmdCalibAdcC = dynamic_cast<TH1F*>(fListOfHistos->At(1));  
194     if(!fPmdCalibAdcC) {
195       printf("ERROR: No ADC File Generated for PMD-CPV ");
196       return;
197     }
198     
199     fPmdCalibEntP = dynamic_cast<TH1F*>(fListOfHistos->At(2));
200     if(!fPmdCalibEntP) {
201       printf("ERROR: No Entry File Generated for PMD-PRE ");
202       printf(" No fhXyPRE ");
203       return;
204     }
205
206     fPmdCalibEntC = dynamic_cast<TH1F*>(fListOfHistos->At(3));
207     if(!fPmdCalibEntC) {
208       printf("ERROR: No Entry File Generated for PMD-CPV ");
209       return;
210     }
211    
212     fNEvents = dynamic_cast<TH1I*>(fListOfHistos->At(4));
213     if(!fNEvents) {
214       printf("ERROR: No Entry File Generated for Event Counter ");
215       return;
216     }
217
218  
219     Info("AliPMDOfflineCalibTask","PMD offline Calibration Successfully finished");
220
221
222 }