]>
Commit | Line | Data |
---|---|---|
f97ea980 | 1 | #include "AliAnalysisTask.h" |
2 | #include "AliAnalysisManager.h" | |
3 | #include "AliAnalysisDataContainer.h" | |
4 | #include "AliITSRecPoint.h" | |
5 | #include "AliESDEvent.h" | |
6 | #include "AliTrackPointArray.h" | |
7 | #include "AliITSgeomTGeo.h" | |
8 | #include "AliESDfriend.h" | |
9 | #include "AliCDBManager.h" | |
10 | #include "AliCDBEntry.h" | |
11 | #include "AliITSCalibrationSDD.h" | |
12 | #include "AliITSresponseSDD.h" | |
13 | #include "AliGeomManager.h" | |
14 | #include <TSystem.h> | |
15 | #include <TTree.h> | |
16 | #include <TH1F.h> | |
17 | #include <TChain.h> | |
18 | #include <TGeoGlobalMagField.h> | |
19 | #include "AliESDInputHandlerRP.h" | |
20 | #include "AliAnalysisTaskSDDRP.h" | |
21 | ||
22 | ClassImp(AliAnalysisTaskSDDRP) | |
23 | //______________________________________________________________________________ | |
24 | AliAnalysisTaskSDDRP::AliAnalysisTaskSDDRP() : AliAnalysisTask("SDD RecPoints",""), | |
25 | fOutput(0), | |
26 | fHistNEvents(0), | |
27 | fRecPMod(0), | |
28 | fTrackPMod(0), | |
29 | fGoodAnMod(0), | |
30 | fRecPLadLay3(0), | |
31 | fRecPLadLay4(0), | |
32 | fTrackPLadLay3(0), | |
33 | fTrackPLadLay4(0), | |
34 | fGoodAnLadLay3(0), | |
35 | fGoodAnLadLay4(0), | |
36 | fESD(0), | |
37 | fESDfriend(0), | |
38 | fResp(0), | |
39 | fRunNumber(0), | |
40 | fMinITSpts(3), | |
41 | fMinPfordEdx(1.5), | |
42 | fOnlyCINT1BTrig(0) | |
43 | { | |
44 | // | |
45 | DefineInput(0, TChain::Class()); | |
46 | DefineOutput(0, TList::Class()); | |
47 | } | |
48 | ||
49 | ||
50 | //___________________________________________________________________________ | |
51 | AliAnalysisTaskSDDRP::~AliAnalysisTaskSDDRP(){ | |
52 | // | |
53 | if (fOutput) { | |
54 | delete fOutput; | |
55 | fOutput = 0; | |
56 | } | |
57 | } | |
58 | ||
59 | //___________________________________________________________________________ | |
60 | void AliAnalysisTaskSDDRP::ConnectInputData(Option_t *) { | |
61 | // | |
62 | TTree* tree = dynamic_cast<TTree*> (GetInputData(0)); | |
63 | if(!tree) { | |
64 | printf("ERROR: Could not read chain from input slot 0\n"); | |
65 | }else{ | |
66 | tree->SetBranchStatus("ESDfriend*", 1); | |
67 | tree->SetBranchAddress("ESDfriend.",&fESDfriend); | |
68 | AliESDInputHandlerRP *hand = dynamic_cast<AliESDInputHandlerRP*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()); | |
69 | fESD=hand->GetEvent(); | |
70 | } | |
71 | return; | |
72 | } | |
73 | //___________________________________________________________________________ | |
74 | ||
75 | void AliAnalysisTaskSDDRP::CreateOutputObjects() { | |
76 | // | |
77 | AliGeomManager::LoadGeometry(fGeomFile.Data()); | |
78 | AliCDBManager* man = AliCDBManager::Instance(); | |
79 | man->SetDefaultStorage("raw://"); | |
80 | man->SetRun(fRunNumber); | |
81 | AliGeomManager::ApplyAlignObjsFromCDB("ITS"); | |
82 | ||
83 | AliCDBEntry* eR=(AliCDBEntry*)man->Get("ITS/Calib/RespSDD"); | |
84 | fResp=(AliITSresponseSDD*)eR->GetObject(); | |
85 | ||
86 | AliCDBEntry* eC=(AliCDBEntry*)man->Get("ITS/Calib/CalibSDD"); | |
87 | TObjArray* calsdd=(TObjArray*)eC->GetObject(); | |
88 | Int_t countGood3[14]; | |
89 | Int_t countGood4[22]; | |
90 | Int_t countGoodMod[260]; | |
91 | for(Int_t ilad=0;ilad<14;ilad++) countGood3[ilad]=0; | |
92 | for(Int_t ilad=0;ilad<22;ilad++) countGood4[ilad]=0; | |
93 | for(Int_t imod=0;imod<260;imod++) countGoodMod[imod]=0; | |
94 | for(Int_t imod=0;imod<260;imod++){ | |
95 | AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)calsdd->At(imod); | |
96 | if(cal->IsBad()) continue; | |
97 | Int_t modId=imod+AliITSgeomTGeo::GetModuleIndex(3,1,1); | |
98 | Int_t lay,lad,det; | |
99 | AliITSgeomTGeo::GetModuleId(modId,lay,lad,det); | |
100 | if(!CheckModule(lay,lad,det)) continue; | |
101 | for(Int_t ian=0; ian<512; ian++){ | |
102 | if(cal->IsBadChannel(ian)) continue; | |
103 | countGoodMod[imod]++; | |
104 | if(lay==3) countGood3[lad-1]++; | |
105 | else if(lay==4) countGood4[lad-1]++; | |
106 | } | |
107 | } | |
108 | ||
109 | ||
110 | fOutput = new TList(); | |
111 | fOutput->SetOwner(); | |
112 | fOutput->SetName("OutputHistos"); | |
113 | ||
114 | fHistNEvents = new TH1F("hNEvents", "Number of processed events",3,-1.5,1.5); | |
115 | fHistNEvents->Sumw2(); | |
116 | fHistNEvents->SetMinimum(0); | |
117 | fOutput->Add(fHistNEvents); | |
118 | ||
119 | // -- Module histos | |
120 | ||
121 | fRecPMod = new TH1F("hRPMod","Rec Points per Module",260,239.5,499.5); | |
122 | fRecPMod->Sumw2(); | |
123 | fRecPMod->SetMinimum(0); | |
124 | fOutput->Add(fRecPMod); | |
125 | ||
126 | fTrackPMod = new TH1F("hTPMod","Track Points per Module",260,239.5,499.5); | |
127 | fTrackPMod->Sumw2(); | |
128 | fTrackPMod->SetMinimum(0); | |
129 | fOutput->Add(fTrackPMod); | |
130 | ||
131 | fGoodAnMod = new TH1F("hGAMod","Good Anodes per Module",260,239.5,499.5); | |
132 | for(Int_t imod=0;imod<260;imod++) fGoodAnMod->SetBinContent(imod+1,countGoodMod[imod]); | |
133 | fGoodAnMod->SetMinimum(0); | |
134 | fOutput->Add(fGoodAnMod); | |
135 | ||
136 | // -- Ladder histos | |
137 | ||
138 | fRecPLadLay3 = new TH1F("hRPLad3","Rec Points per Ladder Layer 3",14,-0.5,13.5); | |
139 | fRecPLadLay3->Sumw2(); | |
140 | fRecPLadLay3->SetMinimum(0); | |
141 | fOutput->Add(fRecPLadLay3); | |
142 | ||
143 | fRecPLadLay4 = new TH1F("hRPLad4","Rec Points per Ladder Layer 4",22,-0.5,21.5); | |
144 | fRecPLadLay4->Sumw2(); | |
145 | fRecPLadLay4->SetMinimum(0); | |
146 | fOutput->Add(fRecPLadLay4); | |
147 | ||
148 | fTrackPLadLay3 = new TH1F("hTPLad3","Track Points per Ladder Layer 3",14,-0.5,13.5); | |
149 | fTrackPLadLay3->Sumw2(); | |
150 | fTrackPLadLay3->SetMinimum(0); | |
151 | fOutput->Add(fTrackPLadLay3); | |
152 | ||
153 | fTrackPLadLay4 = new TH1F("hTPLad4","Track Points per Ladder Layer 4",22,-0.5,21.5); | |
154 | fTrackPLadLay4->Sumw2(); | |
155 | fTrackPLadLay4->SetMinimum(0); | |
156 | fOutput->Add(fTrackPLadLay4); | |
157 | ||
158 | fGoodAnLadLay3 = new TH1F("hGALad3","Good Anodes per Ladder Layer 3",14,-0.5,13.5); | |
159 | for(Int_t ilad=0;ilad<14;ilad++) fGoodAnLadLay3->SetBinContent(ilad+1,countGood3[ilad]); | |
160 | fGoodAnLadLay3->SetMinimum(0); | |
161 | fOutput->Add(fGoodAnLadLay3); | |
162 | ||
163 | fGoodAnLadLay4 = new TH1F("hGALad4","Good Anodes per Ladder Layer 4",22,-0.5,21.5); | |
164 | for(Int_t ilad=0;ilad<22;ilad++) fGoodAnLadLay4->SetBinContent(ilad+1,countGood4[ilad]); | |
165 | fGoodAnLadLay4->SetMinimum(0); | |
166 | fOutput->Add(fGoodAnLadLay4); | |
167 | ||
168 | for(Int_t it=0; it<8; it++){ | |
169 | fSignalTime[it]=new TH1F(Form("hSigTimeInt%d",it),Form("hSigTimeInt%d",it),100,0.,300.); | |
170 | fSignalTime[it]->Sumw2(); | |
171 | fSignalTime[it]->SetMinimum(0); | |
172 | fOutput->Add(fSignalTime[it]); | |
173 | } | |
174 | } | |
175 | //______________________________________________________________________________ | |
176 | void AliAnalysisTaskSDDRP::Exec(Option_t *) | |
177 | { | |
178 | if(!fESD) { | |
179 | printf("AliAnalysisTaskSDDRP::Exec(): bad ESD\n"); | |
180 | return; | |
181 | } | |
182 | if(!fESDfriend) { | |
183 | printf("AliAnalysisTaskSDDRP::Exec(): bad ESDfriend\n"); | |
184 | return; | |
185 | } | |
186 | PostData(0,fOutput); | |
187 | fESD->SetESDfriend(fESDfriend); | |
188 | fHistNEvents->Fill(0); | |
189 | if(fOnlyCINT1BTrig){ | |
190 | if(!fESD->IsTriggerClassFired("CINT1B-ABCE-NOPF-ALL")) return; | |
191 | fHistNEvents->Fill(1); | |
192 | } | |
193 | const AliTrackPointArray *array = 0; | |
194 | Int_t ntracks = fESD->GetNumberOfTracks(); | |
195 | for (Int_t itrack=0; itrack < ntracks; itrack++) { | |
196 | AliESDtrack * track = fESD->GetTrack(itrack); | |
197 | if (!track) continue; | |
198 | if(track->GetNcls(1)>0) continue; | |
199 | if(track->GetNcls(0) < fMinITSpts) continue; | |
200 | Double_t dedx[4]; | |
201 | track->GetITSdEdxSamples(dedx); | |
202 | array = track->GetTrackPointArray(); | |
203 | if(!array) continue; | |
204 | for(Int_t ipt=0; ipt<array->GetNPoints(); ipt++) { | |
205 | AliTrackPoint point; | |
206 | Int_t modId; | |
207 | array->GetPoint(point,ipt); | |
208 | Int_t volId = point.GetVolumeID(); | |
209 | Int_t layerId = AliGeomManager::VolUIDToLayer(volId,modId); | |
210 | if(layerId<3 || layerId>4) continue; | |
211 | modId+=AliITSgeomTGeo::GetModuleIndex(layerId,1,1); | |
212 | Int_t lay,lad,det; | |
213 | AliITSgeomTGeo::GetModuleId(modId,lay,lad,det); | |
214 | if(!CheckModule(lay,lad,det)) continue; | |
215 | fTrackPMod->Fill(modId); | |
216 | Float_t dtime=point.GetDriftTime()-fResp->GetTimeZero(modId); | |
217 | Int_t theBin=int(dtime/6500.*8.); | |
218 | if(layerId==3){ | |
219 | fTrackPLadLay3->Fill(lad-1); | |
220 | if(dedx[0]>0. && track->P()>fMinPfordEdx) fSignalTime[theBin]->Fill(dedx[0]); | |
221 | } | |
222 | if(layerId==4){ | |
223 | fTrackPLadLay4->Fill(lad-1); | |
224 | if(dedx[0]>0.&& track->P()>fMinPfordEdx) fSignalTime[theBin]->Fill(dedx[1]); | |
225 | } | |
226 | } | |
227 | } | |
228 | ||
229 | AliESDInputHandlerRP *hand = dynamic_cast<AliESDInputHandlerRP*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()); | |
230 | TTree* tR = hand->GetTreeR("ITS"); | |
231 | TClonesArray *ITSrec= new TClonesArray("AliITSRecPoint"); | |
232 | TBranch *branch =tR->GetBranch("ITSRecPoints"); | |
233 | branch->SetAddress(&ITSrec); | |
234 | for (Int_t modId=240; modId<500; modId++){ | |
235 | Int_t lay,lad,det; | |
236 | AliITSgeomTGeo::GetModuleId(modId,lay,lad,det); | |
237 | if(!CheckModule(lay,lad,det)) continue; | |
238 | branch->GetEvent(modId); | |
239 | Int_t nrecp = ITSrec->GetEntries(); | |
240 | fRecPMod->Fill(modId,nrecp); | |
241 | if(lay==3) fRecPLadLay3->Fill(lad-1,nrecp); | |
242 | if(lay==4) fRecPLadLay4->Fill(lad-1,nrecp); | |
243 | ||
244 | } | |
245 | ITSrec->Delete(); | |
246 | delete ITSrec; | |
247 | ||
248 | PostData(0,fOutput); | |
249 | ||
250 | } | |
251 | //______________________________________________________________________________ | |
252 | Bool_t AliAnalysisTaskSDDRP::CheckModule(Int_t lay, Int_t lad, Int_t det) const{ | |
253 | // | |
254 | if(lay==4){ | |
255 | if(lad==3 && det==5) return kFALSE; // 1500 V | |
256 | if(lad==3 && det==6) return kFALSE; // 0 V | |
257 | if(lad==3 && det==7) return kFALSE; // 1500 V | |
258 | if(lad==4 && det==1) return kFALSE; // 0 V | |
259 | if(lad==4 && det==2) return kFALSE; // 1500 V | |
260 | if(lad==7 && det==3) return kFALSE; // noisy | |
261 | if(lad==7 && det==5) return kFALSE; // 0 MV + noisy | |
262 | if(lad==9 && det==3) return kFALSE; // 1500 V | |
263 | if(lad==9 && det==4) return kFALSE; // 0 V | |
264 | if(lad==9 && det==5) return kFALSE; // 1500 V | |
265 | if(lad==11 && det==6) return kFALSE; // 1500 V | |
266 | if(lad==11 && det==7) return kFALSE; // 0 V | |
267 | if(lad==11 && det==8) return kFALSE; // 1500 V | |
268 | if(lad==18 && det==5) return kFALSE; // 1500 V | |
269 | if(lad==18 && det==6) return kFALSE; // 0 V | |
270 | if(lad==18 && det==7) return kFALSE; // 1500 V | |
271 | if(lad==22 && det==1) return kFALSE; // 0 V | |
272 | if(lad==22 && det==2) return kFALSE; // 1500 V | |
273 | } | |
274 | if(lay==3){ | |
275 | if(lad==4 && det==4) return kFALSE; // 1500 V | |
276 | if(lad==3) return kFALSE; // swapped in geometry | |
277 | if(lad==5 && det==1) return kFALSE; // noisy | |
278 | if(lad==5 && det==2) return kFALSE; // noisy | |
279 | if(lad==6 && det==1) return kFALSE; // noisy | |
280 | } | |
281 | return kTRUE; | |
282 | } | |
283 | ||
284 | //______________________________________________________________________________ | |
285 | void AliAnalysisTaskSDDRP::Terminate(Option_t */*option*/) | |
286 | { | |
287 | // Terminate analysis | |
288 | fOutput = dynamic_cast<TList*> (GetOutputData(0)); | |
289 | if (!fOutput) { | |
290 | printf("ERROR: fOutput not available\n"); | |
291 | return; | |
292 | } | |
293 | fHistNEvents= dynamic_cast<TH1F*>(fOutput->FindObject("hNEvents")); | |
294 | fRecPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hRPMod")); | |
295 | fTrackPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hTPMod")); | |
296 | fGoodAnMod= dynamic_cast<TH1F*>(fOutput->FindObject("hGAMod")); | |
297 | ||
298 | fRecPLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hRPLad3")); | |
299 | fRecPLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hRPLad4")); | |
300 | fTrackPLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hTPLad3")); | |
301 | fTrackPLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hTPLad4")); | |
302 | fGoodAnLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hGALad3")); | |
303 | fGoodAnLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hGALad4")); | |
304 | ||
305 | for(Int_t it=0; it<8; it++){ | |
306 | fSignalTime[it]= dynamic_cast<TH1F*>(fOutput->FindObject(Form("hSigTimeInt%d",it))); | |
307 | } | |
308 | ||
309 | return; | |
310 | } | |
311 | ||
312 | ||
313 | ||
314 | ||
315 |