]>
Commit | Line | Data |
---|---|---|
ebdd0606 | 1 | /************************************************************************** |
2 | * Copyright(c) 2007-2009, 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 | /* $Id$ */ | |
17 | ||
18 | // ************************************************************* | |
19 | // Checks the quality assurance | |
20 | // by comparing with reference data | |
21 | // contained in a DB | |
22 | // ------------------------------------------------------------- | |
23 | // W. Ferrarese + P. Cerello Feb 2008 | |
8b7e858c | 24 | // M.Siciliano Aug 2008 QA RecPoints |
ebdd0606 | 25 | // INFN Torino |
26 | ||
27 | // --- ROOT system --- | |
ad300de9 | 28 | |
ebdd0606 | 29 | #include <TProfile2D.h> |
30 | #include <TH2D.h> | |
58ceb8ca | 31 | #include <TH1F.h> |
ebdd0606 | 32 | #include <TBranch.h> |
33 | #include <TTree.h> | |
ebdd0606 | 34 | #include <TMath.h> |
3f90b3dc | 35 | #include <TObjArray.h> |
36 | ||
ebdd0606 | 37 | // --- Standard library --- |
38 | ||
39 | // --- AliRoot header files --- | |
40 | #include "AliITSQASDDDataMakerRec.h" | |
4e25ac79 | 41 | #include "AliQAv1.h" |
ebdd0606 | 42 | #include "AliRawReader.h" |
e41192d7 | 43 | #include "AliITSRawStream.h" |
ebdd0606 | 44 | #include "AliITSRawStreamSDD.h" |
5192f264 | 45 | #include "AliITSdigit.h" |
ebdd0606 | 46 | #include "AliITSRecPoint.h" |
e62fe478 | 47 | #include "AliITSRecPointContainer.h" |
ebdd0606 | 48 | #include "AliITSgeomTGeo.h" |
ebdd0606 | 49 | #include "AliCDBManager.h" |
ebdd0606 | 50 | #include "AliCDBEntry.h" |
3f90b3dc | 51 | #include "AliITSCalibrationSDD.h" |
fc89f88c | 52 | |
cfb59c70 | 53 | class TGaxis; |
54 | class TF1; | |
55 | class TSystem; | |
56 | class AliLog; | |
57 | class AliQAChecker; | |
58 | class AliITSRawStreamSDDCompressed; | |
59 | class AliCDBStorage; | |
60 | class Riostream; | |
61 | class AliITSdigitSDD; | |
62 | class AliITS; | |
63 | class AliRunLoader; | |
64 | class AliITSLoader; | |
65 | class AliITSDetTypeRec; | |
66 | ||
67 | ||
ebdd0606 | 68 | |
69 | ClassImp(AliITSQASDDDataMakerRec) | |
70 | ||
71 | //____________________________________________________________________________ | |
72 | AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) : | |
73 | TObject(), | |
74 | fAliITSQADataMakerRec(aliITSQADataMakerRec), | |
75 | fkOnline(kMode), | |
76 | fLDC(ldc), | |
ad300de9 | 77 | fSDDhRawsTask(0), |
44ed7a66 | 78 | fSDDhDigitsTask(0), |
ad300de9 | 79 | fSDDhRecPointsTask(0), |
ad300de9 | 80 | fGenRawsOffset(0), |
44ed7a66 | 81 | fGenDigitsOffset(0), |
ad300de9 | 82 | fGenRecPointsOffset(0), |
ebdd0606 | 83 | fTimeBinSize(1), |
3f90b3dc | 84 | fNEvent(0), |
85 | fNEventRP(0), | |
86 | fDDLModuleMap(0), | |
87 | fCalibration(0), | |
88 | fHistoCalibration(0) | |
ebdd0606 | 89 | { |
90 | //ctor used to discriminate OnLine-Offline analysis | |
8b7e858c | 91 | if(fLDC < 0 || fLDC > 6) { |
ebdd0606 | 92 | AliError("Error: LDC number out of range; return\n"); |
93 | } | |
8b7e858c | 94 | fGenRawsOffset = new Int_t[AliRecoParam::kNSpecies]; |
95 | fGenRecPointsOffset = new Int_t[AliRecoParam::kNSpecies]; | |
4a903927 | 96 | fGenDigitsOffset = new Int_t[AliRecoParam::kNSpecies]; |
8b7e858c | 97 | for(Int_t i=0; i<AliRecoParam::kNSpecies; i++) { |
98 | fGenRawsOffset[i] = 0; | |
99 | fGenRecPointsOffset[i] = 0; | |
4a903927 | 100 | fGenDigitsOffset[i]=0; |
8b7e858c | 101 | } |
3d0c40f9 | 102 | |
103 | InitCalibrationArray(); | |
ebdd0606 | 104 | } |
105 | ||
106 | //____________________________________________________________________________ | |
107 | AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) : | |
108 | TObject(), | |
109 | fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec), | |
110 | fkOnline(qadm.fkOnline), | |
111 | fLDC(qadm.fLDC), | |
ad300de9 | 112 | fSDDhRawsTask(qadm.fSDDhRawsTask), |
44ed7a66 | 113 | fSDDhDigitsTask(qadm.fSDDhDigitsTask), |
ad300de9 | 114 | fSDDhRecPointsTask(qadm.fSDDhRecPointsTask), |
ad300de9 | 115 | fGenRawsOffset(qadm.fGenRawsOffset), |
44ed7a66 | 116 | fGenDigitsOffset(qadm.fGenDigitsOffset), |
ad300de9 | 117 | fGenRecPointsOffset(qadm.fGenRecPointsOffset), |
3f90b3dc | 118 | fTimeBinSize(qadm.fTimeBinSize), |
119 | fNEvent(qadm.fNEvent), | |
120 | fNEventRP(qadm.fNEventRP), | |
121 | fDDLModuleMap(qadm.fDDLModuleMap), | |
122 | fCalibration(qadm.fCalibration), | |
123 | fHistoCalibration(qadm.fHistoCalibration) | |
58ceb8ca | 124 | { |
ebdd0606 | 125 | //copy ctor |
126 | fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; | |
127 | fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle()); | |
3d0c40f9 | 128 | //fDDLModuleMap=NULL; |
ebdd0606 | 129 | } |
130 | ||
131 | //____________________________________________________________________________ | |
132 | AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){ | |
133 | // destructor | |
58ceb8ca | 134 | //if(fDDLModuleMap) delete fDDLModuleMap; |
3f90b3dc | 135 | if(fHistoCalibration){delete fHistoCalibration; fHistoCalibration=NULL;} |
ebdd0606 | 136 | } |
137 | //__________________________________________________________________ | |
138 | AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac ) | |
139 | { | |
7a0e5776 | 140 | // Equal operator. |
ebdd0606 | 141 | this->~AliITSQASDDDataMakerRec(); |
142 | new(this) AliITSQASDDDataMakerRec(qac); | |
143 | return *this; | |
144 | } | |
145 | ||
146 | //____________________________________________________________________________ | |
147 | void AliITSQASDDDataMakerRec::StartOfDetectorCycle() | |
148 | { | |
cfb59c70 | 149 | |
150 | //Start of a QA cycle | |
151 | ||
1e7991d8 | 152 | AliDebug(AliQAv1::GetQADebugLevel(),Form("Start of SDD Cycle with event specie %s for task %s\n",AliRecoParam::GetEventSpecieName(fAliITSQADataMakerRec->GetEventSpecie()),AliQAv1::GetTaskName(fAliITSQADataMakerRec->GetTaskIndexSelected()).Data())); |
3d0c40f9 | 153 | if(!fCalibration) {CreateTheCalibration();} |
154 | ||
1e7991d8 | 155 | if(fAliITSQADataMakerRec->GetEventSpecie()==0) return;//not the kDefault EventSpecie |
156 | //Detector specific actions at start of cycle | |
157 | if(fAliITSQADataMakerRec->GetTaskIndexSelected()==AliQAv1::kRAWS){ | |
158 | AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SDD Cycle\n"); | |
159 | if(fAliITSQADataMakerRec->ListExists(AliQAv1::kRAWS)==kFALSE)return; | |
3f90b3dc | 160 | |
1e7991d8 | 161 | AliDebug(AliQAv1::GetQADebugLevel(),Form("Reset of Raw Data normalized histograms with eventspecie %s ",AliRecoParam::GetEventSpecieName(fAliITSQADataMakerRec->GetEventSpecie()))); |
162 | ||
163 | fAliITSQADataMakerRec->GetRawsData(3+ fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); | |
164 | fAliITSQADataMakerRec->GetRawsData(4+ fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); | |
165 | fAliITSQADataMakerRec->GetRawsData(5+ fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); | |
166 | ||
167 | } | |
168 | if(fAliITSQADataMakerRec->GetTaskIndexSelected()==AliQAv1::kRECPOINTS){ | |
169 | if(fAliITSQADataMakerRec->ListExists(AliQAv1::kRECPOINTS)==kFALSE)return; | |
170 | ||
171 | AliDebug(AliQAv1::GetQADebugLevel(),Form("Reset of RecPoints normalized histograms with eventspecie %s ",AliRecoParam::GetEventSpecieName(fAliITSQADataMakerRec->GetEventSpecie()))); | |
172 | ||
173 | fAliITSQADataMakerRec->GetRecPointsData(9+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); | |
174 | fAliITSQADataMakerRec->GetRecPointsData(10+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); | |
175 | fAliITSQADataMakerRec->GetRecPointsData(11+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); | |
176 | } | |
ebdd0606 | 177 | } |
178 | ||
179 | //____________________________________________________________________________ | |
3f90b3dc | 180 | void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray* /*list*/) |
ebdd0606 | 181 | { |
cfb59c70 | 182 | //end of a QA cycle |
8b7e858c | 183 | AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); |
3f90b3dc | 184 | if(task==AliQAv1::kRAWS){ |
fad6e2fc | 185 | // printf("fNevent %d \n",fNEvent); |
3f90b3dc | 186 | if(fNEvent!=0){ |
187 | ((TH1D*)fAliITSQADataMakerRec->GetRawsData(3 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH1D*)fAliITSQADataMakerRec->GetRawsData(0 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH1D*) (fHistoCalibration->At(0))),1.,(Double_t)fNEvent); | |
188 | ||
189 | ((TH2D*)fAliITSQADataMakerRec->GetRawsData(4 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRawsData(1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(1))),1.,(Double_t)fNEvent); | |
190 | ||
191 | ((TH2D*)fAliITSQADataMakerRec->GetRawsData(5 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRawsData(2 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(2))),1.,(Double_t)fNEvent); | |
192 | } | |
193 | }//end raws | |
194 | ||
3d0c40f9 | 195 | if(task==AliQAv1::kRECPOINTS){ |
fad6e2fc | 196 | // printf("fNeventRP %d \n",fNEventRP); |
3f90b3dc | 197 | if(fNEventRP!=0){ |
198 | ((TH1D*)fAliITSQADataMakerRec->GetRecPointsData(9 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH1D*)fAliITSQADataMakerRec->GetRecPointsData(6 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH1D*) (fHistoCalibration->At(0))),1.,(Double_t)fNEventRP); | |
199 | ||
200 | ((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(10+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(7 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(1))),1.,(Double_t)fNEventRP); | |
201 | ||
202 | ((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(11+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(8 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(2))),1.,(Double_t)fNEventRP); | |
203 | } | |
204 | }//end recpoints | |
205 | ||
ebdd0606 | 206 | } |
207 | ||
208 | //____________________________________________________________________________ | |
eca4fa66 | 209 | Int_t AliITSQASDDDataMakerRec::InitRaws() |
ebdd0606 | 210 | { |
3d0c40f9 | 211 | |
ebdd0606 | 212 | // Initialization for RAW data - SDD - |
7d297381 | 213 | const Bool_t expert = kTRUE ; |
214 | const Bool_t saveCorr = kTRUE ; | |
215 | const Bool_t image = kTRUE ; | |
58ceb8ca | 216 | |
eca4fa66 | 217 | Int_t rv = 0 ; |
ebdd0606 | 218 | Int_t lay, lad, det; |
ebdd0606 | 219 | Int_t indexlast = 0; |
220 | Int_t index1 = 0; | |
221 | ||
3f90b3dc | 222 | if(fkOnline){AliInfo("Book Online Histograms for SDD\n");} |
223 | else {AliInfo("Book Offline Histograms for SDD\n ");} | |
7a0e5776 | 224 | TH1D *h0 = new TH1D("SDDModPattern","HW Modules pattern",fgknSDDmodules,239.5,499.5); //0 |
ebdd0606 | 225 | h0->GetXaxis()->SetTitle("Module Number"); |
226 | h0->GetYaxis()->SetTitle("Counts"); | |
ce8c37f2 | 227 | rv = fAliITSQADataMakerRec->Add2RawsList(h0,0+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); |
ad300de9 | 228 | fSDDhRawsTask++; |
7a0e5776 | 229 | |
230 | //zPhi distribution using ladder and modules numbers | |
3f90b3dc | 231 | TH2D *hphil3 = new TH2D("SDDphizL3","SDD #varphiz Layer3 ",12,0.5,6.5,14,0.5,14.5);//1 |
5af1d436 | 232 | hphil3->GetXaxis()->SetTitle("z[Module Number L3 ]"); |
233 | hphil3->GetYaxis()->SetTitle("#varphi[ Ladder Number L3]"); | |
ce8c37f2 | 234 | rv = fAliITSQADataMakerRec->Add2RawsList(hphil3,1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, saveCorr); |
7a0e5776 | 235 | fSDDhRawsTask++; |
236 | ||
3f90b3dc | 237 | TH2D *hphil4 = new TH2D("SDDphizL4","SDD #varphiz Layer4 ",16,0.5,8.5,22,0.5,22.5); //2 |
5af1d436 | 238 | hphil4->GetXaxis()->SetTitle("z[Module Number L4]"); |
239 | hphil4->GetYaxis()->SetTitle("#varphi[Ladder Number L4]"); | |
ce8c37f2 | 240 | rv = fAliITSQADataMakerRec->Add2RawsList(hphil4,2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, saveCorr); |
7a0e5776 | 241 | fSDDhRawsTask++; |
242 | ||
3f90b3dc | 243 | //normalized histograms |
244 | TH1D *h0norm = new TH1D("SDDModPatternNORM","NORM HW Modules pattern",fgknSDDmodules,239.5,499.5); //3 | |
245 | h0norm->GetXaxis()->SetTitle("Module Number"); | |
246 | h0norm->GetYaxis()->SetTitle("Counts"); | |
247 | rv = fAliITSQADataMakerRec->Add2RawsList(h0norm,3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); | |
248 | fSDDhRawsTask++; | |
249 | ||
250 | //zPhi distribution using ladder and modules numbers | |
251 | TH2D *hphil3norm = new TH2D("SDDphizL3NORM","NORM SDD #varphiz Layer3 ",12,0.5,6.5,14,0.5,14.5);//4 | |
252 | hphil3norm->GetXaxis()->SetTitle("z[Module Number L3 ]"); | |
253 | hphil3norm->GetYaxis()->SetTitle("#varphi[ Ladder Number L3]"); | |
254 | rv = fAliITSQADataMakerRec->Add2RawsList(hphil3norm,4+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); | |
255 | fSDDhRawsTask++; | |
256 | ||
257 | TH2D *hphil4norm = new TH2D("SDDphizL4NORM","NORM SDD #varphiz Layer4 ",16,0.5,8.5,22,0.5,22.5); //5 | |
258 | hphil4norm->GetXaxis()->SetTitle("z[Module Number L4]"); | |
259 | hphil4norm->GetYaxis()->SetTitle("#varphi[Ladder Number L4]"); | |
260 | rv = fAliITSQADataMakerRec->Add2RawsList(hphil4norm,5+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); | |
261 | fSDDhRawsTask++; | |
ebdd0606 | 262 | |
3f90b3dc | 263 | if(fkOnline){ |
7a0e5776 | 264 | //DDL Pattern |
3f90b3dc | 265 | TH2D *hddl = new TH2D("SDDDDLPattern","SDD DDL Pattern ",24,-0.5,11.5,24,-0.5,23.5); //6 |
ad300de9 | 266 | hddl->GetXaxis()->SetTitle("Channel"); |
5af1d436 | 267 | hddl->GetYaxis()->SetTitle("DDL Number"); |
3f90b3dc | 268 | rv = fAliITSQADataMakerRec->Add2RawsList(hddl,6+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); |
ad300de9 | 269 | fSDDhRawsTask++; |
7a0e5776 | 270 | Int_t indexlast1 = 0; |
271 | ||
272 | fTimeBinSize = 4; | |
273 | indexlast = 0; | |
274 | index1 = 0; | |
275 | indexlast1 = fSDDhRawsTask; | |
276 | char *hname[3]; | |
277 | for(Int_t i=0; i<3; i++) hname[i]= new char[50]; | |
278 | for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){ | |
279 | for(Int_t iside=0;iside<fgknSide;iside++){ | |
280 | AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det); | |
281 | sprintf(hname[0],"SDDchargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside); | |
282 | sprintf(hname[1],"SDDChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside); | |
fc89f88c | 283 | // sprintf(hname[2],"SDDhmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside); |
7a0e5776 | 284 | TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5); |
285 | fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin"); | |
286 | fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode"); | |
80b9610c | 287 | rv = fAliITSQADataMakerRec->Add2RawsList(fModuleChargeMapFSE,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); |
7a0e5776 | 288 | fSDDhRawsTask++; |
289 | index1++; | |
ebdd0606 | 290 | } |
7a0e5776 | 291 | } |
292 | ||
293 | for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){ | |
294 | for(Int_t iside=0;iside<fgknSide;iside++){ | |
295 | AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det); | |
296 | sprintf(hname[0],"SDDchargeMap_L%d_%d_%d_%d",lay,lad,det,iside); | |
297 | sprintf(hname[1],"SDDChargeMap_L%d_%d_%d_%d",lay,lad,det,iside); | |
298 | TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5); | |
299 | fModuleChargeMap->GetXaxis()->SetTitle("Time Bin"); | |
5af1d436 | 300 | fModuleChargeMap->GetYaxis()->SetTitle("Anode Number"); |
80b9610c | 301 | rv = fAliITSQADataMakerRec->Add2RawsList(fModuleChargeMap,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); |
7a0e5776 | 302 | fSDDhRawsTask++; |
303 | index1++; | |
ebdd0606 | 304 | } |
7a0e5776 | 305 | } |
306 | ||
58ceb8ca | 307 | //Event Size |
3f90b3dc | 308 | TH1F *hsize = new TH1F("SDDEventSize","SDD Event Size ",500,-0.5,199.5); |
309 | hsize->SetBit(TH1::kCanRebin); | |
58ceb8ca | 310 | hsize->GetXaxis()->SetTitle("Event Size [kB]"); |
311 | hsize->GetYaxis()->SetTitle("Entries"); | |
ce8c37f2 | 312 | rv = fAliITSQADataMakerRec->Add2RawsList(hsize,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr); |
58ceb8ca | 313 | fSDDhRawsTask++; |
8b7e858c | 314 | |
7a0e5776 | 315 | } // kONLINE |
5379c4a3 | 316 | AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Raws histograms booked\n",fSDDhRawsTask)); |
eca4fa66 | 317 | return rv ; |
ebdd0606 | 318 | } |
319 | ||
320 | ||
321 | //____________________________________________________________________________ | |
eca4fa66 | 322 | Int_t AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader) |
ebdd0606 | 323 | { |
324 | // Fill QA for RAW - SDD - | |
3f90b3dc | 325 | Int_t rv = 0; |
7555afef | 326 | // Check id histograms already created for this Event Specie |
3f90b3dc | 327 | fNEvent++; |
328 | if(!fDDLModuleMap){CreateTheMap();} | |
eca4fa66 | 329 | if(rawReader->GetType() != 7) return rv; // skips non physical triggers |
5379c4a3 | 330 | AliDebug(AliQAv1::GetQADebugLevel(),"entering MakeRaws\n"); |
e41192d7 | 331 | rawReader->Reset(); |
58ceb8ca | 332 | rawReader->Select("ITSSDD"); |
8b7e858c | 333 | AliITSRawStream *stream=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader); |
334 | stream->SetDDLModuleMap(fDDLModuleMap); | |
e41192d7 | 335 | |
ebdd0606 | 336 | Int_t lay, lad, det; |
7a0e5776 | 337 | |
ebdd0606 | 338 | Int_t index=0; |
339 | if(fkOnline) { | |
340 | for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){ | |
7a0e5776 | 341 | for(Int_t iside=0;iside<fgknSide;iside++) { |
3f90b3dc | 342 | if(fSDDhRawsTask > 7 + index) fAliITSQADataMakerRec->GetRawsData(7 + index +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); |
7a0e5776 | 343 | // 4 because the 2D histos for single events start after the fourth position |
344 | index++; | |
345 | } | |
ebdd0606 | 346 | } |
347 | } | |
e41192d7 | 348 | |
ebdd0606 | 349 | Int_t cnt = 0; |
350 | Int_t ildcID = -1; | |
351 | Int_t iddl = -1; | |
352 | Int_t isddmod = -1; | |
7a0e5776 | 353 | Int_t coord1, coord2, signal, moduleSDD, activeModule, index1; |
58ceb8ca | 354 | //if(fkOnline) |
355 | //{ | |
356 | Int_t prevDDLID = -1; | |
357 | UInt_t size = 0; | |
358 | int totalddl=static_cast<int>(fDDLModuleMap->GetNDDLs()); | |
359 | Bool_t *ddldata=new Bool_t[totalddl]; | |
360 | for(Int_t jddl=0;jddl<totalddl;jddl++){ddldata[jddl]=kFALSE;} | |
361 | //} | |
e41192d7 | 362 | while(stream->Next()) { |
ebdd0606 | 363 | ildcID = rawReader->GetLDCId(); |
58ceb8ca | 364 | iddl = rawReader->GetDDLID();// - fgkDDLIDshift; |
3f90b3dc | 365 | //printf("----------------------iddl %i\n",iddl); |
58ceb8ca | 366 | |
e41192d7 | 367 | isddmod = fDDLModuleMap->GetModuleNumber(iddl,stream->GetCarlosId()); |
ebdd0606 | 368 | if(isddmod==-1){ |
5379c4a3 | 369 | AliDebug(AliQAv1::GetQADebugLevel(),Form("Found module with iddl: %d, stream->GetCarlosId: %d \n",iddl,stream->GetCarlosId())); |
ebdd0606 | 370 | continue; |
371 | } | |
e41192d7 | 372 | if(stream->IsCompletedModule()) { |
5379c4a3 | 373 | AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedModule == KTRUE\n")); |
ebdd0606 | 374 | continue; |
375 | } | |
de075dae | 376 | if(stream->IsCompletedDDL()) { |
58ceb8ca | 377 | if(fkOnline){ |
378 | if ((rawReader->GetDDLID() != prevDDLID)&&(ddldata[iddl])==kFALSE){ | |
379 | size += rawReader->GetDataSize();//in bytes | |
380 | prevDDLID = rawReader->GetDDLID(); | |
381 | ddldata[iddl]=kTRUE; | |
382 | } | |
383 | } | |
5379c4a3 | 384 | AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedDDL == KTRUE\n")); |
de075dae | 385 | continue; |
386 | } | |
ebdd0606 | 387 | |
e41192d7 | 388 | coord1 = stream->GetCoord1(); |
389 | coord2 = stream->GetCoord2(); | |
390 | signal = stream->GetSignal(); | |
ebdd0606 | 391 | |
7a0e5776 | 392 | moduleSDD = isddmod - fgkmodoffset; |
393 | ||
b62d7123 | 394 | if(isddmod <fgkmodoffset|| isddmod>fgknSDDmodules+fgkmodoffset-1) { |
5379c4a3 | 395 | AliDebug(AliQAv1::GetQADebugLevel(),Form( "Module SDD = %d, resetting it to 1 \n",isddmod)); |
b62d7123 | 396 | isddmod = 1; |
ebdd0606 | 397 | } |
7a0e5776 | 398 | |
ebdd0606 | 399 | AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det); |
e41192d7 | 400 | Short_t iside = stream->GetChannel(); |
ad300de9 | 401 | |
58ceb8ca | 402 | //printf(" \n%i %i %i %i \n ",lay, lad, det,iside ); |
403 | fAliITSQADataMakerRec->GetRawsData( 0 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()] )->Fill(isddmod); | |
404 | if(lay==3) fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad); | |
3f90b3dc | 405 | if(lay==4) {fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);} |
58ceb8ca | 406 | |
ebdd0606 | 407 | if(fkOnline) { |
fc89f88c | 408 | |
3f90b3dc | 409 | fAliITSQADataMakerRec->GetRawsData(6+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill((stream->GetCarlosId())+0.5*iside -0.5,iddl); |
410 | // printf("content ddlmap %d, %d = %f \n",(stream->GetCarlosId()+0.5*iside -0.5),iddl,fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetBinContent(1+(det-1)*2;lad)); | |
411 | //printf("content ddlmap %d, %d = %f \n",(stream->GetCarlosId())+0.5*iside -0.5,iddl,fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetBinContent(1+(stream->GetCarlosId()-1)*2,iddl)); | |
7a0e5776 | 412 | activeModule = moduleSDD; |
413 | index1 = activeModule * 2 + iside; | |
414 | ||
415 | if(index1<0){ | |
5379c4a3 | 416 | AliDebug(AliQAv1::GetQADebugLevel(),Form("Wrong index number %d - patched to 0\n",index1)); |
7a0e5776 | 417 | index1 = 0; |
418 | } | |
78f8430c | 419 | |
3f90b3dc | 420 | if(fSDDhRawsTask > 7 + index1) { |
421 | ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(7 + index1 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(coord2, coord1, signal); | |
422 | ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(7 + index1 + 260*2 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(coord2, coord1, signal); | |
ebdd0606 | 423 | } |
58ceb8ca | 424 | }//online |
7a0e5776 | 425 | cnt++; |
5379c4a3 | 426 | if(!(cnt%10000)) AliDebug(AliQAv1::GetQADebugLevel(),Form(" %d raw digits read",cnt)); |
58ceb8ca | 427 | }//end next() |
3f90b3dc | 428 | if(fkOnline){((TH1F*)(fAliITSQADataMakerRec->GetRawsData(7 + 260*4 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(size/1024.);//KB |
429 | } | |
5379c4a3 | 430 | AliDebug(AliQAv1::GetQADebugLevel(),Form("Event completed, %d raw digits read",cnt)); |
e41192d7 | 431 | delete stream; |
432 | stream = NULL; | |
8b7e858c | 433 | |
434 | // if(fkOnline) { | |
435 | // AnalyseBNG(); // Analyse Baseline, Noise, Gain | |
436 | // AnalyseINJ(); // Analyse Injectors | |
437 | // } | |
438 | ||
58ceb8ca | 439 | delete []ddldata; |
440 | ddldata=NULL; | |
eca4fa66 | 441 | return rv ; |
7a0e5776 | 442 | } |
ebdd0606 | 443 | |
44ed7a66 | 444 | //____________________________________________________________________________ |
eca4fa66 | 445 | Int_t AliITSQASDDDataMakerRec::InitDigits() |
44ed7a66 | 446 | { |
3d0c40f9 | 447 | // if(!fHistoCalibration)InitCalibrationArray(); |
44ed7a66 | 448 | // Initialization for DIGIT data - SDD - |
449 | const Bool_t expert = kTRUE ; | |
450 | const Bool_t image = kTRUE ; | |
eca4fa66 | 451 | Int_t rv = 0 ; |
452 | // fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries(); | |
44ed7a66 | 453 | //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List |
454 | TH1F* h0=new TH1F("SDD DIGITS Module Pattern","SDD DIGITS Module Pattern",260,239.5,499.5); //hmod | |
455 | h0->GetXaxis()->SetTitle("SDD Module Number"); | |
456 | h0->GetYaxis()->SetTitle("# DIGITS"); | |
4a903927 | 457 | rv = fAliITSQADataMakerRec->Add2DigitsList(h0,fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image); |
44ed7a66 | 458 | fSDDhDigitsTask ++; |
4a903927 | 459 | // printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask ); |
44ed7a66 | 460 | TH1F* h1=new TH1F("SDD Anode Distribution","DIGITS Anode Distribution",512,-0.5,511.5); //hanocc |
461 | h1->GetXaxis()->SetTitle("Anode Number"); | |
462 | h1->GetYaxis()->SetTitle("# DIGITS"); | |
4a903927 | 463 | rv = fAliITSQADataMakerRec->Add2DigitsList(h1,1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image); |
44ed7a66 | 464 | fSDDhDigitsTask ++; |
4a903927 | 465 | //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask ); |
44ed7a66 | 466 | TH1F* h2=new TH1F("SDD Tbin Distribution","DIGITS Tbin Distribution",256,-0.5,255.5); //htbocc |
467 | h2->GetXaxis()->SetTitle("Tbin Number"); | |
468 | h2->GetYaxis()->SetTitle("# DIGITS"); | |
4a903927 | 469 | rv = fAliITSQADataMakerRec->Add2DigitsList(h2,2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image); |
44ed7a66 | 470 | fSDDhDigitsTask ++; |
4a903927 | 471 | //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask ); |
44ed7a66 | 472 | TH1F* h3=new TH1F("SDD ADC Counts Distribution","DIGITS ADC Counts Distribution",200,0.,1024.); //hsig |
473 | h3->GetXaxis()->SetTitle("ADC Value"); | |
474 | h3->GetYaxis()->SetTitle("# DIGITS"); | |
4a903927 | 475 | rv = fAliITSQADataMakerRec->Add2DigitsList(h3,3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image); |
44ed7a66 | 476 | fSDDhDigitsTask ++; |
4a903927 | 477 | //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask ); |
44ed7a66 | 478 | AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Digits histograms booked\n",fSDDhDigitsTask)); |
eca4fa66 | 479 | return rv ; |
44ed7a66 | 480 | } |
481 | ||
482 | //____________________________________________________________________________ | |
eca4fa66 | 483 | Int_t AliITSQASDDDataMakerRec::MakeDigits(TTree * digits) |
44ed7a66 | 484 | { |
7555afef | 485 | |
44ed7a66 | 486 | // Fill QA for DIGIT - SDD - |
fc89f88c | 487 | //AliITS *fITS = (AliITS*)gAlice->GetModule("ITS"); |
488 | //fITS->SetTreeAddress(); | |
489 | //TClonesArray *iITSdigits = fITS->DigitsAddress(1); | |
490 | ||
491 | ||
eca4fa66 | 492 | Int_t rv = 0 ; |
493 | ||
fc89f88c | 494 | TBranch *branchD = digits->GetBranch("ITSDigitsSDD"); |
495 | ||
496 | if (!branchD) { | |
497 | AliError("can't get the branch with the ITS SDD digits !"); | |
eca4fa66 | 498 | return rv ; |
44ed7a66 | 499 | } |
eca4fa66 | 500 | // Check id histograms already created for this Event Specie |
8b7e858c | 501 | // if ( ! fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset) ) |
502 | // rv = InitDigits() ; | |
eca4fa66 | 503 | |
fc89f88c | 504 | static TClonesArray statDigits("AliITSdigitSDD"); |
44ed7a66 | 505 | TClonesArray *iITSdigits = &statDigits; |
506 | branchD->SetAddress(&iITSdigits); | |
fc89f88c | 507 | |
44ed7a66 | 508 | for(Int_t i=0; i<260; i++){ |
509 | Int_t nmod=i+240; | |
510 | digits->GetEvent(nmod); | |
511 | Int_t ndigits = iITSdigits->GetEntries(); | |
4a903927 | 512 | fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nmod,ndigits); |
7555afef | 513 | |
44ed7a66 | 514 | for (Int_t idig=0; idig<ndigits; idig++) { |
515 | AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig); | |
516 | Int_t iz=dig->GetCoord1(); // cell number z | |
517 | Int_t ix=dig->GetCoord2(); // cell number x | |
518 | Int_t sig=dig->GetSignal(); | |
4a903927 | 519 | fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(iz); |
520 | fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(ix); | |
521 | fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(sig); | |
44ed7a66 | 522 | } |
523 | } | |
eca4fa66 | 524 | return rv ; |
44ed7a66 | 525 | } |
526 | ||
ebdd0606 | 527 | //____________________________________________________________________________ |
eca4fa66 | 528 | Int_t AliITSQASDDDataMakerRec::InitRecPoints() |
ebdd0606 | 529 | { |
3d0c40f9 | 530 | //if(!fHistoCalibration)InitCalibrationArray(); |
8b7e858c | 531 | //AliInfo("Initialize SDD recpoints histos\n"); |
ebdd0606 | 532 | // Initialization for RECPOINTS - SDD - |
7d297381 | 533 | const Bool_t expert = kTRUE ; |
534 | const Bool_t image = kTRUE ; | |
eca4fa66 | 535 | Int_t rv = 0 ; |
536 | // fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries(); | |
ad300de9 | 537 | |
8b7e858c | 538 | //AliInfo(Form("fAliITSQADataMakerRec->GetEventSpecie() %d\n",fAliITSQADataMakerRec->GetEventSpecie())); |
539 | //AliInfo(Form("fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()] %d\n",fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])); | |
be2d913a | 540 | TH1F *h0 = new TH1F("SDDLay3TotCh","Layer 3 total charge",250,-0.5, 499.5); //position number 0 |
3f90b3dc | 541 | //h0->SetBit(TH1::kCanRebin); |
cfb59c70 | 542 | h0->GetXaxis()->SetTitle("KeV"); |
ebdd0606 | 543 | h0->GetYaxis()->SetTitle("Entries"); |
80b9610c | 544 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h0, 0 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image |
ad300de9 | 545 | fSDDhRecPointsTask++; |
ebdd0606 | 546 | |
be2d913a | 547 | TH1F *h1 = new TH1F("SDDLay4TotCh","Layer 4 total charge",250,-0.5, 499.5);//position number 1 |
3f90b3dc | 548 | //h1->SetBit(TH1::kCanRebin); |
cfb59c70 | 549 | h1->GetXaxis()->SetTitle("Kev"); |
ebdd0606 | 550 | h1->GetYaxis()->SetTitle("Entries"); |
80b9610c | 551 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h1, 1 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image |
ad300de9 | 552 | fSDDhRecPointsTask++; |
ebdd0606 | 553 | |
554 | char hisnam[50]; | |
be2d913a | 555 | TH2F *h2 = new TH2F("SDDGlobalCoordDistribYX","YX Global Coord Distrib",56,-28,28,56,-28,28);//position number 2 |
ebdd0606 | 556 | h2->GetYaxis()->SetTitle("Y[cm]"); |
557 | h2->GetXaxis()->SetTitle("X[cm]"); | |
80b9610c | 558 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h2,2+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image |
ad300de9 | 559 | fSDDhRecPointsTask++; |
ebdd0606 | 560 | |
be2d913a | 561 | TH2F *h3 = new TH2F("SDDGlobalCoordDistribRZ","RZ Global Coord Distrib",128,-32,32,56,12,26);//position number 3 |
ebdd0606 | 562 | h3->GetYaxis()->SetTitle("R[cm]"); |
563 | h3->GetXaxis()->SetTitle("Z[cm]"); | |
80b9610c | 564 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h3,3+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image |
ad300de9 | 565 | fSDDhRecPointsTask++; |
ebdd0606 | 566 | |
be2d913a | 567 | TH2F *h4 = new TH2F("SDDGlobalCoordDistribL3PHIZ","#varphi Z Global Coord Distrib L3",46,-23,23,90,-TMath::Pi(),TMath::Pi());//position number 4 |
ebdd0606 | 568 | h4->GetYaxis()->SetTitle("#phi[rad]"); |
569 | h4->GetXaxis()->SetTitle("Z[cm]"); | |
80b9610c | 570 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h4,4+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image |
ad300de9 | 571 | fSDDhRecPointsTask++; |
ebdd0606 | 572 | |
be2d913a | 573 | TH2F *h5 = new TH2F("SDDGlobalCoordDistribL4PHIZ","#varphi Z Global Coord Distrib L4",62,-31,31,90,-TMath::Pi(),TMath::Pi());//position number 5 |
ebdd0606 | 574 | h5->GetYaxis()->SetTitle("#phi[rad]"); |
575 | h5->GetXaxis()->SetTitle("Z[cm]"); | |
80b9610c | 576 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h5,5+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image |
ad300de9 | 577 | fSDDhRecPointsTask++; |
ebdd0606 | 578 | |
3f90b3dc | 579 | TH1D *h6 = new TH1D("SDDModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); //position number 6 |
ebdd0606 | 580 | h6->GetXaxis()->SetTitle("Module number"); //spd offset = 240 |
581 | h6->GetYaxis()->SetTitle("Entries"); | |
80b9610c | 582 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h6,6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image |
ad300de9 | 583 | fSDDhRecPointsTask++; |
58ceb8ca | 584 | |
3f90b3dc | 585 | |
586 | TH2D *h7 = new TH2D("SDDModPatternL3RP","Modules pattern L3 RP",12,0.5,6.5,14,0.5,14.5); //position number 7 | |
5af1d436 | 587 | h7->GetXaxis()->SetTitle("z[#Module L3 ]"); |
588 | h7->GetYaxis()->SetTitle("#varphi[#Ladder L3]"); | |
3f90b3dc | 589 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h7,7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image |
ad300de9 | 590 | fSDDhRecPointsTask++; |
58ceb8ca | 591 | |
3f90b3dc | 592 | TH2D *h8 = new TH2D("SDDModPatternL4RP","Modules pattern L4 RP",16,0.5,8.5,22,0.5,22.5); //position number 8 |
5af1d436 | 593 | h8->GetXaxis()->SetTitle("[#Module L3 ]"); |
594 | h8->GetYaxis()->SetTitle("#varphi[#Ladder L4]"); | |
3f90b3dc | 595 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h8,8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image |
596 | fSDDhRecPointsTask++; | |
597 | ||
598 | //------------------------norm--------------------------// | |
599 | ||
600 | ||
601 | TH1D *h9 = new TH1D("SDDModPatternRPNORM","Modules pattern RP NORM",fgknSDDmodules,239.5,499.5); //position number 9 | |
602 | h9->GetXaxis()->SetTitle("Module number"); //spd offset = 240 | |
603 | h9->GetYaxis()->SetTitle("Entries"); | |
604 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h9,9 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
ad300de9 | 605 | fSDDhRecPointsTask++; |
58ceb8ca | 606 | |
3f90b3dc | 607 | |
608 | TH2D *h10 = new TH2D("SDDModPatternL3RPNORM","Modules pattern L3 RP NORM",12,0.5,6.5,14,0.5,14.5); //position number 10 | |
609 | h10->GetXaxis()->SetTitle("z[#Module L3 ]"); | |
610 | h10->GetYaxis()->SetTitle("#varphi[#Ladder L3]"); | |
611 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h10,10 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
612 | fSDDhRecPointsTask++; | |
613 | ||
614 | TH2D *h11 = new TH2D("SDDModPatternL4RPNORM","Modules pattern L4 RP NORM",16,0.5,8.5,22,0.5,22.5); //position number 11 | |
615 | h11->GetXaxis()->SetTitle("[#Module L3 ]"); | |
616 | h10->GetYaxis()->SetTitle("#varphi[#Ladder L4]"); | |
617 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h11,11 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
618 | fSDDhRecPointsTask++; | |
619 | ||
620 | //--------------------------------------------------------// | |
621 | ||
be2d913a | 622 | TH2F *h12 = new TH2F("SDDLocalCoordDistrib","Local Coord Distrib",160,-4,4,160,-4,4);//position number 12 |
3f90b3dc | 623 | h12->GetXaxis()->SetTitle("X local coord, drift, cm"); |
624 | h12->GetYaxis()->SetTitle("Z local coord, anode, cm"); | |
625 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h12,12 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
ad300de9 | 626 | fSDDhRecPointsTask++; |
58ceb8ca | 627 | |
628 | //AliInfo("Create SDD recpoints histos\n"); | |
629 | ||
3f90b3dc | 630 | TH1F *h13 = new TH1F("SDDrdistrib_Layer3" ,"SDD r distribution Layer3" ,100,14.,16.5);//position number 13 (L3) |
631 | h13->SetBit(TH1::kCanRebin); | |
632 | h13->GetXaxis()->SetTitle("r[cm]"); | |
633 | h13->GetXaxis()->CenterTitle(); | |
634 | h13->GetYaxis()->SetTitle("Entries"); | |
635 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h13,13 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
58ceb8ca | 636 | fSDDhRecPointsTask++; |
637 | ||
3f90b3dc | 638 | TH1F *h14 = new TH1F("SDDrdistrib_Layer4" ,"SDD r distribution Layer4" ,100,23.,25.);// and position number 14 (L4) |
639 | h14->SetBit(TH1::kCanRebin); | |
640 | h14->GetXaxis()->SetTitle("r[cm]"); | |
641 | h14->GetXaxis()->CenterTitle(); | |
642 | h14->GetYaxis()->SetTitle("Entries"); | |
643 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h14,14 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
58ceb8ca | 644 | fSDDhRecPointsTask++; |
645 | ||
ebdd0606 | 646 | for(Int_t iLay=0; iLay<=1; iLay++){ |
647 | sprintf(hisnam,"SDDphidistrib_Layer%d",iLay+3); | |
3f90b3dc | 648 | TH1F *h15 = new TH1F(hisnam,hisnam,180,-TMath::Pi(),TMath::Pi());//position number 15 (L3) and position number 16 (L4) |
649 | h15->GetXaxis()->SetTitle("#varphi[rad]"); | |
650 | h15->GetXaxis()->CenterTitle(); | |
651 | h15->GetYaxis()->SetTitle("Entries"); | |
652 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h15,iLay+15+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
ad300de9 | 653 | fSDDhRecPointsTask++; |
ebdd0606 | 654 | } |
58ceb8ca | 655 | |
656 | for(Int_t iLay=0; iLay<=1; iLay++){ | |
657 | sprintf(hisnam,"SDDdrifttime_Layer%d",iLay+3); | |
be2d913a | 658 | TH1F *h17 = new TH1F(hisnam,hisnam,45,-0.5,4499.5);//position number 17 (L3) and position number 18 (L4) |
3f90b3dc | 659 | h17->SetBit(TH1::kCanRebin); |
1e7991d8 | 660 | h17->GetXaxis()->SetTitle("drift time[ns]"); |
3f90b3dc | 661 | h17->GetXaxis()->CenterTitle(); |
662 | h17->GetYaxis()->SetTitle("Entries"); | |
663 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h17,iLay+17+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image | |
58ceb8ca | 664 | fSDDhRecPointsTask++; |
665 | } | |
666 | ||
3f90b3dc | 667 | if(fkOnline){ |
be2d913a | 668 | TH2F *h19 = new TH2F("SDDGlobalCoordDistribYXFSE","YX Global Coord Distrib FSE",112,-28,28,112,-28,28);//position number 19 |
3f90b3dc | 669 | h19->GetYaxis()->SetTitle("Y[cm]"); |
670 | h19->GetXaxis()->SetTitle("X[cm]"); | |
671 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h19,19+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
ad300de9 | 672 | fSDDhRecPointsTask++; |
ebdd0606 | 673 | |
be2d913a | 674 | TH2F *h20 = new TH2F("SDDGlobalCoordDistribRZFSE","RZ Global Coord Distrib FSE",128,-32,32,56,12,26);//position number 20 |
3f90b3dc | 675 | h20->GetYaxis()->SetTitle("R[cm]"); |
676 | h20->GetXaxis()->SetTitle("Z[cm]"); | |
677 | rv = fAliITSQADataMakerRec->Add2RecPointsList(h20,20+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image | |
678 | fSDDhRecPointsTask++; | |
ebdd0606 | 679 | }//online |
58ceb8ca | 680 | |
5379c4a3 | 681 | AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Recs histograms booked\n",fSDDhRecPointsTask)); |
58ceb8ca | 682 | |
eca4fa66 | 683 | return rv ; |
ebdd0606 | 684 | } |
685 | ||
686 | //____________________________________________________________________________ | |
eca4fa66 | 687 | Int_t AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree) |
7a0e5776 | 688 | { |
ebdd0606 | 689 | // Fill QA for RecPoints - SDD - |
3f90b3dc | 690 | Int_t rv = 0 ; |
691 | fNEventRP++; | |
ebdd0606 | 692 | Int_t lay, lad, det; |
8b7e858c | 693 | //AliInfo("get the branch with the ITS clusters !\n"); |
e62fe478 | 694 | AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance(); |
aa444532 | 695 | TClonesArray *recpoints=NULL; |
3f90b3dc | 696 | if(fkOnline){recpoints = rpcont->FetchClusters(0,clustersTree,fAliITSQADataMakerRec->GetEventNumber());} |
697 | else{recpoints = rpcont->FetchClusters(0,clustersTree);} | |
e62fe478 | 698 | if(!rpcont->GetStatusOK() || !rpcont->IsSDDActive()){ |
699 | AliError("can't get SDD clusters !"); | |
700 | return rv; | |
ebdd0606 | 701 | } |
702 | ||
ebdd0606 | 703 | Int_t npoints = 0; |
704 | Float_t cluglo[3]={0.,0.,0.}; | |
3f90b3dc | 705 | if(fkOnline){ |
706 | for(Int_t i=19;i<21;i++){ | |
8b7e858c | 707 | fAliITSQADataMakerRec->GetRecPointsData(i+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset(); |
ebdd0606 | 708 | } |
709 | } | |
e62fe478 | 710 | Int_t firMod=AliITSgeomTGeo::GetModuleIndex(3,1,1); |
711 | Int_t lasMod=AliITSgeomTGeo::GetModuleIndex(5,1,1); | |
712 | for(Int_t module=firMod; module<lasMod;module++){ | |
58ceb8ca | 713 | //AliInfo(Form("Module %d\n",module)); |
e62fe478 | 714 | recpoints = rpcont->UncheckedGetClusters(module); |
58ceb8ca | 715 | npoints += recpoints->GetEntries(); |
716 | //AliInfo(Form("modnumb %d, npoints %d, total points %d\n",module, recpoints->GetEntries(),npoints)); | |
80b9610c | 717 | //AliITSgeomTGeo::GetModuleId(module, lay, lad, det); |
58ceb8ca | 718 | //AliInfo(Form("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det)); |
719 | ||
720 | //AliInfo(Form("modnumb %d, entries %d\n",module, recpoints->GetEntries())); | |
721 | for(Int_t j=0;j<recpoints->GetEntries();j++){ | |
722 | //AliInfo(Form("modnumb %d, entry %d \n",module, j)); | |
80b9610c | 723 | AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j); |
724 | Int_t index = recp->GetDetectorIndex(); | |
725 | lay=recp->GetLayer(); | |
726 | Int_t modnumb=index+AliITSgeomTGeo::GetModuleIndex(lay+1,1,1); | |
727 | //printf("modnumb %i\n",modnumb); | |
728 | AliITSgeomTGeo::GetModuleId(modnumb, lay, lad, det); | |
729 | fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(modnumb);//modpatternrp | |
58ceb8ca | 730 | recp->GetGlobalXYZ(cluglo); |
731 | Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); | |
732 | Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]); | |
733 | Float_t drifttime=recp->GetDriftTime(); | |
3f90b3dc | 734 | fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());//local distribution |
80b9610c | 735 | if(lay==3||lay==4)fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[0],cluglo[1]);//global distribution YX |
58ceb8ca | 736 | fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],rad);//global distribution rz |
737 | if(fkOnline) { | |
3f90b3dc | 738 | fAliITSQADataMakerRec->GetRecPointsData(19 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[0],cluglo[1]);//global distribution YX FSE |
739 | fAliITSQADataMakerRec->GetRecPointsData(20 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],rad);//global distribution rz FSE | |
58ceb8ca | 740 | } |
741 | if(recp->GetLayer() == 2) { | |
742 | fAliITSQADataMakerRec->GetRecPointsData(0 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetQ()) ;//total charge of layer 3 | |
5af1d436 | 743 | //fAliITSQADataMakerRec->GetRecPointsData(7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lad);//lad pattern layer 3 |
3f90b3dc | 744 | Int_t iside=recp->GetDriftSide(); |
745 | //printf("iside =%d\n",iside); | |
746 | fAliITSQADataMakerRec->GetRecPointsData(7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);//mod pattern layer 3 | |
747 | fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);//r distribution layer 3 | |
748 | fAliITSQADataMakerRec->GetRecPointsData(15 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);// phi distribution layer 3 | |
58ceb8ca | 749 | fAliITSQADataMakerRec->GetRecPointsData(4 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],phi);// zphi distribution layer |
3f90b3dc | 750 | fAliITSQADataMakerRec->GetRecPointsData(17 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(drifttime);// time distribution layer 3 |
58ceb8ca | 751 | } else if(recp->GetLayer() == 3) { |
752 | fAliITSQADataMakerRec->GetRecPointsData(1 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetQ()) ;//total charge layer 4 | |
5af1d436 | 753 | //fAliITSQADataMakerRec->GetRecPointsData(8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lad);//ladpatternlayer4 |
3f90b3dc | 754 | Int_t iside=recp->GetDriftSide(); |
755 | //printf("iside =%d\n",iside); | |
756 | fAliITSQADataMakerRec->GetRecPointsData(8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);//mod pattern layer 4 | |
757 | fAliITSQADataMakerRec->GetRecPointsData(14 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);//r distribution | |
758 | fAliITSQADataMakerRec->GetRecPointsData(16 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);//phi distribution | |
58ceb8ca | 759 | fAliITSQADataMakerRec->GetRecPointsData(5 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],phi);// zphi distribution layer 4 |
3f90b3dc | 760 | fAliITSQADataMakerRec->GetRecPointsData(18 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(drifttime);// time distribution layer 4 |
58ceb8ca | 761 | } |
762 | } | |
763 | } | |
58ceb8ca | 764 | return rv ; |
ebdd0606 | 765 | } |
766 | ||
ad300de9 | 767 | //_______________________________________________________________ |
768 | ||
cfb59c70 | 769 | Int_t AliITSQASDDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task, Int_t specie)const |
7a0e5776 | 770 | { |
cfb59c70 | 771 | // Returns offset number according to the specified task |
ad300de9 | 772 | Int_t offset=0; |
3f90b3dc | 773 | if( task == AliQAv1::kRAWS ){offset=fGenRawsOffset[specie];} |
774 | else if(task == AliQAv1::kDIGITSR ){offset=fGenDigitsOffset[specie];} | |
775 | else if( task == AliQAv1::kRECPOINTS ){offset=fGenRecPointsOffset[specie];} | |
ad300de9 | 776 | return offset; |
777 | } | |
778 | ||
779 | //_______________________________________________________________ | |
780 | ||
8b7e858c | 781 | void AliITSQASDDDataMakerRec::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset, Int_t specie) { |
cfb59c70 | 782 | // Set offset number according to the specified task |
3f90b3dc | 783 | if( task == AliQAv1::kRAWS ) {fGenRawsOffset[specie]=offset;} |
784 | else if( task == AliQAv1::kDIGITSR ) {fGenDigitsOffset[specie]=offset;} | |
785 | else if( task == AliQAv1::kRECPOINTS ) {fGenRecPointsOffset[specie]=offset;} | |
eca4fa66 | 786 | } |
787 | ||
788 | //_______________________________________________________________ | |
789 | ||
4e25ac79 | 790 | Int_t AliITSQASDDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task) |
7a0e5776 | 791 | { |
cfb59c70 | 792 | //return the number of histo booked for a given Task |
7a0e5776 | 793 | Int_t histotot=0; |
3f90b3dc | 794 | if( task == AliQAv1::kRAWS ){ histotot=fSDDhRawsTask ;} |
795 | else if(task == AliQAv1::kDIGITSR) { histotot=fSDDhDigitsTask;} | |
796 | else if( task == AliQAv1::kRECPOINTS ){ histotot=fSDDhRecPointsTask;} | |
797 | else {AliInfo("No task has been selected. TaskHisto set to zero.\n");} | |
ad300de9 | 798 | return histotot; |
e41192d7 | 799 | } |
8b7e858c | 800 | |
58ceb8ca | 801 | |
802 | //_______________________________________________________________ | |
803 | ||
804 | ||
805 | void AliITSQASDDDataMakerRec::CreateTheMap() | |
806 | { | |
cfb59c70 | 807 | //Create the SDD DDL Module Map |
58ceb8ca | 808 | AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD"); |
809 | Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag(); | |
3f90b3dc | 810 | if(!ddlMapSDD){ |
58ceb8ca | 811 | AliError("Calibration object retrieval failed! SDD will not be processed"); |
812 | fDDLModuleMap = NULL; | |
813 | //return rv; | |
814 | } | |
815 | fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject(); | |
816 | if(!cacheStatus)ddlMapSDD->SetObject(NULL); | |
817 | ddlMapSDD->SetOwner(kTRUE); | |
3f90b3dc | 818 | if(!cacheStatus){ delete ddlMapSDD;} |
58ceb8ca | 819 | AliInfo("DDL Map Created\n "); |
820 | } | |
5af1d436 | 821 | |
3f90b3dc | 822 | //_______________________________________________________________ |
823 | ||
824 | ||
825 | void AliITSQASDDDataMakerRec::CreateTheCalibration() | |
826 | { | |
cfb59c70 | 827 | //Take from the OCDB the calibration information for the SDD |
3f90b3dc | 828 | AliCDBEntry *calibSDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD"); |
829 | Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag(); | |
830 | if(!calibSDD) | |
831 | { | |
832 | AliError("Calibration object retrieval failed! SDD will not be processed"); | |
833 | fCalibration = NULL;; | |
834 | } | |
835 | fCalibration = (TObjArray *)calibSDD->GetObject(); | |
836 | ||
837 | if(!cacheStatus)calibSDD->SetObject(NULL); | |
838 | calibSDD->SetOwner(kTRUE); | |
839 | if(!cacheStatus){delete calibSDD;} | |
840 | ||
841 | AliITSCalibrationSDD * cal=NULL; | |
842 | for(Int_t imod=0;imod<fgknSDDmodules;imod++) | |
843 | { | |
844 | //cal=NULL; | |
845 | Int_t fillmodhisto1=fgkTotalNumberSDDAnodes; | |
846 | Int_t fillmodhisto2side0=fgkNumberOfSDDAnodesperSide; | |
847 | Int_t fillmodhisto2side1=fgkNumberOfSDDAnodesperSide; | |
848 | Int_t fillmodhisto3side0=fgkNumberOfSDDAnodesperSide; | |
849 | Int_t fillmodhisto3side1=fgkNumberOfSDDAnodesperSide; | |
850 | ||
851 | Int_t badmodhisto1=0; | |
852 | Int_t badmodhisto2side0=0; | |
853 | Int_t badmodhisto2side1=0; | |
854 | Int_t badmodhisto3side0=0; | |
855 | Int_t badmodhisto3side1=0; | |
856 | //printf("imod %i\t ==== \t",imod); | |
857 | Int_t module=imod + 240; | |
858 | //printf("module %i\t ==== \t",module); | |
859 | cal=(AliITSCalibrationSDD*)fCalibration->At(imod); | |
860 | Int_t lay,lad,det; | |
861 | AliITSgeomTGeo::GetModuleId(module,lay,lad,det); | |
862 | Int_t index=1+(det-1)*2; | |
863 | if(cal==0){continue;} | |
864 | if (cal->IsBad()){continue;}//bad module check | |
865 | else{ | |
866 | for(Int_t i=0;i<8;i++) //check on bad chips in good modules | |
867 | { | |
868 | if(lay==3){ | |
869 | if(cal->IsChipBad(i)){ | |
870 | if(i<4){badmodhisto2side0+=64;} | |
871 | if(i>=4){badmodhisto2side1+=64;} | |
872 | }//end if chip | |
873 | }//end if layer3 | |
874 | else if(lay==4){ | |
875 | if(cal->IsChipBad(i)){ | |
876 | if(i<4){badmodhisto3side0+=64;} | |
877 | if(i>=4){badmodhisto3side1+=64;} | |
878 | }//end if chip | |
879 | }//ens if layer4 | |
880 | }//end for chip | |
881 | for(Int_t iAn=0; iAn<512; iAn++){//anodes loop | |
882 | Int_t ic=cal->GetChip(iAn);//chip with this anode number | |
883 | if(!cal->IsChipBad(ic) && !cal->IsBad() && cal->IsBadChannel(iAn)){// good chip good module bad channel | |
884 | if(lay==3){ | |
885 | if(ic<4) badmodhisto2side0++; | |
886 | else if(ic>=4)badmodhisto2side1++; | |
887 | }//end if layer 3 | |
888 | else if(lay==4){ | |
889 | if(ic<4) badmodhisto3side0++; | |
890 | else if(ic>=4)badmodhisto3side1++; | |
891 | }//end if layer 4 | |
892 | }//end if chip module channel | |
893 | }//end for anodes | |
894 | if(lay==3){ | |
895 | badmodhisto1=badmodhisto2side0+badmodhisto2side1; | |
896 | fillmodhisto1-=badmodhisto1; | |
897 | fillmodhisto2side0-=badmodhisto2side0; | |
898 | fillmodhisto2side1-=badmodhisto2side1; | |
899 | ((TH1D*)(fHistoCalibration->At(0)))->SetBinContent(imod+1,fillmodhisto1); | |
900 | ((TH2D*)(fHistoCalibration->At(1)))->SetBinContent(index,lad,fillmodhisto2side0); | |
901 | ((TH2D*)(fHistoCalibration->At(1)))->SetBinContent(index+1,lad,fillmodhisto2side1); | |
902 | }//end layer 3 | |
903 | else if(lay==4){ | |
904 | badmodhisto1=badmodhisto3side0+badmodhisto3side1; | |
905 | fillmodhisto1-=badmodhisto1; | |
906 | fillmodhisto3side0-=badmodhisto3side0; | |
907 | fillmodhisto3side1-=badmodhisto3side1; | |
908 | ((TH1D*)(fHistoCalibration->At(0)))->SetBinContent(imod+1,fillmodhisto1); | |
909 | ((TH2D*)(fHistoCalibration->At(2)))->SetBinContent(index,lad,fillmodhisto3side0); | |
910 | ((TH2D*)(fHistoCalibration->At(2)))->SetBinContent(index+1,lad,fillmodhisto2side1); | |
911 | }//end layer 4 | |
912 | }//end else bad module | |
913 | }//end module for | |
914 | ||
915 | } | |
916 | ||
917 | //____________________________________________________________________ | |
918 | ||
919 | void AliITSQASDDDataMakerRec::InitCalibrationArray() | |
920 | { | |
cfb59c70 | 921 | //create the histograms with the calibration informations. The histograms are stored in a TObjArray |
3f90b3dc | 922 | TH1D *pattern1 = new TH1D("CALSDDModPattern","Calibration HW Modules pattern",fgknSDDmodules,239.5,499.5); |
923 | TH2D *patternl3 = new TH2D("CALSDDphizL3","Calibration SDD #varphiz Layer3 ",12,0.5,6.5,14,0.5,14.5); | |
924 | TH2D *patternl4 = new TH2D("CALSDDphizL4"," Calibration SDD #varphiz Layer4 ",16,0.5,8.5,22,0.5,22.5); | |
925 | ||
926 | fHistoCalibration = new TObjArray(3); | |
927 | fHistoCalibration->AddAtAndExpand(pattern1,0); | |
928 | fHistoCalibration->AddAtAndExpand(patternl3,1); | |
929 | fHistoCalibration->AddAtAndExpand(patternl4,2); | |
930 | ||
fad6e2fc | 931 | // printf("Calibration Histograms created!\n"); |
3f90b3dc | 932 | } |
5af1d436 | 933 | |
934 | //____________________________________________________________________ | |
935 | ||
936 | void AliITSQASDDDataMakerRec::ResetDetector(AliQAv1::TASKINDEX_t task) | |
937 | { | |
cfb59c70 | 938 | //reset the SDD calibration histograms |
5af1d436 | 939 | AliInfo(Form("Reset detector in SDD called for task index %i", task)); |
940 | if(task== AliQAv1::kRAWS ){ | |
5af1d436 | 941 | fDDLModuleMap=NULL; |
3d0c40f9 | 942 | } |
3f90b3dc | 943 | fCalibration=NULL; |
3d0c40f9 | 944 | |
3f90b3dc | 945 | ((TH1D*)(fHistoCalibration->At(0)))->Reset(); |
946 | ((TH2D*)(fHistoCalibration->At(1)))->Reset(); | |
947 | ((TH2D*)(fHistoCalibration->At(2)))->Reset(); | |
3d0c40f9 | 948 | |
5af1d436 | 949 | } |
be2d913a | 950 |