]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PMD/anal/AliPmdCustomCalibTask.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PMD / anal / AliPmdCustomCalibTask.cxx
CommitLineData
233407cb 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#include "Riostream.h"
17#include "TChain.h"
18#include "TTree.h"
19#include "TH1F.h"
20#include "TH2F.h"
21#include "TCanvas.h"
22#include "TList.h"
23
24#include "AliAnalysisTaskSE.h"
25#include "AliAnalysisManager.h"
26#include "AliStack.h"
27#include "AliESDtrackCuts.h"
28#include "AliESDEvent.h"
29#include "AliVEvent.h"
30#include "AliESDInputHandler.h"
31#include "AliAODEvent.h"
32#include "AliMCEvent.h"
33#include "AliESDPmdTrack.h"
34
35#include "AliPmdCustomCalibTask.h"
36
37/*------------------------------------------------------------------------
38 . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
39 .
40 . PMD Offiline Calibration Task to run Over
41 . ESDs - can be used in user task
42 . - runs for single module in both plan
43 . - Centrality is not needed for PbPb
44 . - (Let me know if you need any updates)
45 .
46 . Satyajit Jena, IIT Bombay
47 . sjena@cern.ch
48 . 3/8/2011
49 .
50 .
51 . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . .
52 ------------------------------------------------------------------------*/
53
54ClassImp(AliPmdCustomCalibTask)
55
56//________________________________________________________________________
57AliPmdCustomCalibTask::AliPmdCustomCalibTask() : AliAnalysisTaskSE(),
58 fOutput(0),
59 fSmn(10),
60 fHistAdcPre(0),
61 fHistAdcCpv(0),
62 fHistClusterXYPre(0),
63 fHistClusterXYCpv(0) {
64 for(Int_t i =0;i<48;i++) {
65 for(Int_t j=0;j<96;j++) {
66 fHistAdcPreRC[i][j] = NULL;
67 fHistAdcCpvRC[i][j] = NULL;
68 }
69 }
70}
71
72//________________________________________________________________________
73AliPmdCustomCalibTask::AliPmdCustomCalibTask(const char *name)
74 :AliAnalysisTaskSE(name),
75 fOutput(0),
76 fSmn(10),
77 fHistAdcPre(0),
78 fHistAdcCpv(0),
79 fHistClusterXYPre(0),
80 fHistClusterXYCpv(0) {
81 for(Int_t i =0;i<48;i++) {
82 for(Int_t j=0;j<96;j++) {
83 fHistAdcPreRC[i][j] = NULL;
84 fHistAdcCpvRC[i][j] = NULL;
85 }
86 }
87 DefineOutput(1, TList::Class()); }
88
89//________________________________________________________________________
90AliPmdCustomCalibTask::~AliPmdCustomCalibTask() {
91 if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) delete fOutput;
92}
93
94//________________________________________________________________________
95void AliPmdCustomCalibTask::UserCreateOutputObjects() {
96 fOutput = new TList();
97 fOutput->SetOwner(kTRUE);
98
99 // Module wise ADC
100 fHistAdcPre = new TH1F("fHistAdcPre", "ADC Distribution for PRE module", 120, 0, 1200);
101 fHistAdcCpv = new TH1F("fHistAdcCpv", "ADC Distribution for CPV module", 120, 0, 1200);
102
103 // XY QA
104 fHistClusterXYPre = new TH2F("fHistClusterXYPre", "ClusterX-Y Plot for PRE module", 400,-200,200,400,-200,200);
105 fHistClusterXYCpv = new TH2F("fHistClusterXYCpv", "ClusterX-Y Plot for CPV module", 400,-200,200,400,-200,200);
106
107 fOutput->Add(fHistAdcPre);
108 fOutput->Add(fHistAdcCpv);
109 fOutput->Add(fHistClusterXYPre);
110 fOutput->Add(fHistClusterXYCpv);
111
112 Char_t *name = new Char_t[100];
113 Char_t *title = new Char_t[200];
114
115 // Cell-wise ADC
116 for(Int_t i = 0; i < 48; i++) {
117 for(Int_t j = 0; j < 96; j++) {
118 sprintf(name,"fHistPreAdcRaw%02dCol%02d",i,j);
119 sprintf(title,"ADC Distribution for PRE | SM %d - Raw%02dCol%02d",fSmn,i,j);
120 fHistAdcPreRC[i][j] = new TH1F(name,title,120,0,1200);
121 fOutput->Add(fHistAdcPreRC[i][j]);
122
123 sprintf(name,"fHistCpvAdcRaw%02dCol%02d",i,j);
124 sprintf(title,"ADC Distribution for CPV | SM %d - Raw%02dCol%02d",fSmn,i,j);
125 fHistAdcCpvRC[i][j] = new TH1F(name,title,120,0,1200);
126 fOutput->Add(fHistAdcCpvRC[i][j]);
127 }
128 }
129
130 PostData(1, fOutput);
131}
132
133//________________________________________________________________________
134void AliPmdCustomCalibTask::UserExec(Option_t *) {
135 AliVEvent *event = InputEvent();
136
137 if (!event) { Printf("ERROR: Could not retrieve event"); return; }
138 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(event);
139 if (!esd) {
140 AliError("Cannot get the ESD event");
141 return;
142 }
143
144
145
146 Int_t npmdcl = esd->GetNumberOfPmdTracks();
147
148 while (npmdcl--) {
149 AliESDPmdTrack *pmdtr = esd->GetPmdTrack(npmdcl);
150
151 Int_t smn = pmdtr->GetSmn();
152
153 if (smn != fSmn) continue; // if not the SM go back
154
155 Float_t adc = pmdtr->GetClusterADC();
156 if(adc>1200) continue;
157
158 Float_t sTag = pmdtr->GetClusterSigmaX();
159 if(sTag < 999. || sTag > 1000.) continue; // Isolated Cell Section
160
161 Int_t det = pmdtr->GetDetector();
162
163 Float_t clsX = pmdtr->GetClusterX();
164 Float_t clsY = pmdtr->GetClusterY();
165
166 Float_t sRowCol = pmdtr->GetClusterSigmaY();
167 Int_t row = (Int_t)sRowCol/100;
168 Int_t col = (Int_t)sRowCol - row*100;
169
170 if(det == 0) {
171 fHistAdcPre->Fill(adc);
172 fHistAdcPreRC[row][col]->Fill(adc);
173 fHistClusterXYPre->Fill(clsX,clsY);
174 } else if(det == 1) {
175 fHistAdcCpv->Fill(adc);
176 fHistAdcCpvRC[row][col]->Fill(adc);
177 fHistClusterXYCpv->Fill(clsX,clsY);
178 }
179 }
180
181 PostData(1, fOutput);
182}
183
184//________________________________________________________________________
185void AliPmdCustomCalibTask::Terminate(Option_t *) {
186 fOutput = dynamic_cast<TList*> (GetOutputData(1));
187 if(!fOutput) { Printf("ERROR: could not retrieve TList fOutput"); return; }
188}
189
190