1 /**************************************************************************
2 * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // *************************************************************
19 // Checks the quality assurance
20 // by comparing with reference data
22 // -------------------------------------------------------------
23 // W. Ferrarese + P. Cerello Feb 2008
24 // M.Siciliano Aug 2008 QA RecPoints
27 // --- ROOT system ---
29 #include <TProfile2D.h>
37 #include <TDirectory.h>
39 #include <TObjArray.h>
41 // --- Standard library ---
43 // --- AliRoot header files ---
44 #include "AliITSQASDDDataMakerRec.h"
47 #include "AliQAChecker.h"
48 #include "AliRawReader.h"
49 #include "AliITSRawStream.h"
50 #include "AliITSRawStreamSDD.h"
51 #include "AliITSRawStreamSDDCompressed.h"
52 #include "AliITSDetTypeRec.h"
53 #include "AliITSdigit.h"
54 #include "AliITSRecPoint.h"
55 #include "AliITSRecPointContainer.h"
56 #include "AliITSgeomTGeo.h"
57 #include "AliCDBManager.h"
58 #include "AliCDBStorage.h"
59 #include "AliCDBEntry.h"
60 #include "Riostream.h"
61 #include "AliITSdigitSDD.h"
63 #include "AliRunLoader.h"
64 #include "AliITSLoader.h"
65 #include "AliITSDetTypeRec.h"
66 #include "AliITSCalibrationSDD.h"
69 ClassImp(AliITSQASDDDataMakerRec)
71 //____________________________________________________________________________
72 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) :
74 fAliITSQADataMakerRec(aliITSQADataMakerRec),
79 fSDDhRecPointsTask(0),
82 fGenRecPointsOffset(0),
90 //ctor used to discriminate OnLine-Offline analysis
91 if(fLDC < 0 || fLDC > 6) {
92 AliError("Error: LDC number out of range; return\n");
94 fGenRawsOffset = new Int_t[AliRecoParam::kNSpecies];
95 fGenRecPointsOffset = new Int_t[AliRecoParam::kNSpecies];
96 fGenDigitsOffset = new Int_t[AliRecoParam::kNSpecies];
97 for(Int_t i=0; i<AliRecoParam::kNSpecies; i++) {
98 fGenRawsOffset[i] = 0;
99 fGenRecPointsOffset[i] = 0;
100 fGenDigitsOffset[i]=0;
103 InitCalibrationArray();
106 //____________________________________________________________________________
107 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) :
109 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
110 fkOnline(qadm.fkOnline),
112 fSDDhRawsTask(qadm.fSDDhRawsTask),
113 fSDDhDigitsTask(qadm.fSDDhDigitsTask),
114 fSDDhRecPointsTask(qadm.fSDDhRecPointsTask),
115 fGenRawsOffset(qadm.fGenRawsOffset),
116 fGenDigitsOffset(qadm.fGenDigitsOffset),
117 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
118 fTimeBinSize(qadm.fTimeBinSize),
119 fNEvent(qadm.fNEvent),
120 fNEventRP(qadm.fNEventRP),
121 fDDLModuleMap(qadm.fDDLModuleMap),
122 fCalibration(qadm.fCalibration),
123 fHistoCalibration(qadm.fHistoCalibration)
126 fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ;
127 fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
128 //fDDLModuleMap=NULL;
131 //____________________________________________________________________________
132 AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){
134 //if(fDDLModuleMap) delete fDDLModuleMap;
135 if(fHistoCalibration){delete fHistoCalibration; fHistoCalibration=NULL;}
137 //__________________________________________________________________
138 AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac )
141 this->~AliITSQASDDDataMakerRec();
142 new(this) AliITSQASDDDataMakerRec(qac);
146 //____________________________________________________________________________
147 void AliITSQASDDDataMakerRec::StartOfDetectorCycle()
150 if(!fCalibration) {CreateTheCalibration();}
152 //Detector specific actions at start of cycle
153 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SDD Cycle\n");
154 if(fAliITSQADataMakerRec->GetRawsData(0)!=NULL){
155 fAliITSQADataMakerRec->GetRawsData(3+ fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
156 fAliITSQADataMakerRec->GetRawsData(4+ fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
157 fAliITSQADataMakerRec->GetRawsData(5+ fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
160 if(fAliITSQADataMakerRec->GetRecPointsData(0)!=NULL){
161 fAliITSQADataMakerRec->GetRecPointsData(9+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
162 fAliITSQADataMakerRec->GetRecPointsData(10+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
163 fAliITSQADataMakerRec->GetRecPointsData(11+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
167 //____________________________________________________________________________
168 void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray* /*list*/)
170 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n");
171 if(task==AliQAv1::kRAWS){
172 // printf("fNevent %d \n",fNEvent);
174 ((TH1D*)fAliITSQADataMakerRec->GetRawsData(3 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH1D*)fAliITSQADataMakerRec->GetRawsData(0 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH1D*) (fHistoCalibration->At(0))),1.,(Double_t)fNEvent);
176 ((TH2D*)fAliITSQADataMakerRec->GetRawsData(4 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRawsData(1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(1))),1.,(Double_t)fNEvent);
178 ((TH2D*)fAliITSQADataMakerRec->GetRawsData(5 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRawsData(2 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(2))),1.,(Double_t)fNEvent);
182 if(task==AliQAv1::kRECPOINTS){
183 // printf("fNeventRP %d \n",fNEventRP);
185 ((TH1D*)fAliITSQADataMakerRec->GetRecPointsData(9 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH1D*)fAliITSQADataMakerRec->GetRecPointsData(6 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH1D*) (fHistoCalibration->At(0))),1.,(Double_t)fNEventRP);
187 ((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(10+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(7 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(1))),1.,(Double_t)fNEventRP);
189 ((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(11+ fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]))->Divide((TH2D*)fAliITSQADataMakerRec->GetRecPointsData(8 + fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]),((TH2D*)(fHistoCalibration->At(2))),1.,(Double_t)fNEventRP);
195 //____________________________________________________________________________
196 Int_t AliITSQASDDDataMakerRec::InitRaws()
199 // Initialization for RAW data - SDD -
200 const Bool_t expert = kTRUE ;
201 const Bool_t saveCorr = kTRUE ;
202 const Bool_t image = kTRUE ;
209 if(fkOnline){AliInfo("Book Online Histograms for SDD\n");}
210 else {AliInfo("Book Offline Histograms for SDD\n ");}
211 TH1D *h0 = new TH1D("SDDModPattern","HW Modules pattern",fgknSDDmodules,239.5,499.5); //0
212 h0->GetXaxis()->SetTitle("Module Number");
213 h0->GetYaxis()->SetTitle("Counts");
214 rv = fAliITSQADataMakerRec->Add2RawsList(h0,0+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
217 //zPhi distribution using ladder and modules numbers
218 TH2D *hphil3 = new TH2D("SDDphizL3","SDD #varphiz Layer3 ",12,0.5,6.5,14,0.5,14.5);//1
219 hphil3->GetXaxis()->SetTitle("z[Module Number L3 ]");
220 hphil3->GetYaxis()->SetTitle("#varphi[ Ladder Number L3]");
221 rv = fAliITSQADataMakerRec->Add2RawsList(hphil3,1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, saveCorr);
224 TH2D *hphil4 = new TH2D("SDDphizL4","SDD #varphiz Layer4 ",16,0.5,8.5,22,0.5,22.5); //2
225 hphil4->GetXaxis()->SetTitle("z[Module Number L4]");
226 hphil4->GetYaxis()->SetTitle("#varphi[Ladder Number L4]");
227 rv = fAliITSQADataMakerRec->Add2RawsList(hphil4,2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, saveCorr);
230 //normalized histograms
231 TH1D *h0norm = new TH1D("SDDModPatternNORM","NORM HW Modules pattern",fgknSDDmodules,239.5,499.5); //3
232 h0norm->GetXaxis()->SetTitle("Module Number");
233 h0norm->GetYaxis()->SetTitle("Counts");
234 rv = fAliITSQADataMakerRec->Add2RawsList(h0norm,3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
237 //zPhi distribution using ladder and modules numbers
238 TH2D *hphil3norm = new TH2D("SDDphizL3NORM","NORM SDD #varphiz Layer3 ",12,0.5,6.5,14,0.5,14.5);//4
239 hphil3norm->GetXaxis()->SetTitle("z[Module Number L3 ]");
240 hphil3norm->GetYaxis()->SetTitle("#varphi[ Ladder Number L3]");
241 rv = fAliITSQADataMakerRec->Add2RawsList(hphil3norm,4+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
244 TH2D *hphil4norm = new TH2D("SDDphizL4NORM","NORM SDD #varphiz Layer4 ",16,0.5,8.5,22,0.5,22.5); //5
245 hphil4norm->GetXaxis()->SetTitle("z[Module Number L4]");
246 hphil4norm->GetYaxis()->SetTitle("#varphi[Ladder Number L4]");
247 rv = fAliITSQADataMakerRec->Add2RawsList(hphil4norm,5+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
252 TH2D *hddl = new TH2D("SDDDDLPattern","SDD DDL Pattern ",24,-0.5,11.5,24,-0.5,23.5); //6
253 hddl->GetXaxis()->SetTitle("Channel");
254 hddl->GetYaxis()->SetTitle("DDL Number");
255 rv = fAliITSQADataMakerRec->Add2RawsList(hddl,6+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
257 Int_t indexlast1 = 0;
262 indexlast1 = fSDDhRawsTask;
264 for(Int_t i=0; i<3; i++) hname[i]= new char[50];
265 for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
266 for(Int_t iside=0;iside<fgknSide;iside++){
267 AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
268 sprintf(hname[0],"SDDchargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside);
269 sprintf(hname[1],"SDDChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside);
270 // sprintf(hname[2],"SDDhmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);
271 TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
272 fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin");
273 fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode");
274 rv = fAliITSQADataMakerRec->Add2RawsList(fModuleChargeMapFSE,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
280 for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
281 for(Int_t iside=0;iside<fgknSide;iside++){
282 AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
283 sprintf(hname[0],"SDDchargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
284 sprintf(hname[1],"SDDChargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
285 TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
286 fModuleChargeMap->GetXaxis()->SetTitle("Time Bin");
287 fModuleChargeMap->GetYaxis()->SetTitle("Anode Number");
288 rv = fAliITSQADataMakerRec->Add2RawsList(fModuleChargeMap,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
295 TH1F *hsize = new TH1F("SDDEventSize","SDD Event Size ",500,-0.5,199.5);
296 hsize->SetBit(TH1::kCanRebin);
297 hsize->GetXaxis()->SetTitle("Event Size [kB]");
298 hsize->GetYaxis()->SetTitle("Entries");
299 rv = fAliITSQADataMakerRec->Add2RawsList(hsize,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
303 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Raws histograms booked\n",fSDDhRawsTask));
308 //____________________________________________________________________________
309 Int_t AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader)
311 // Fill QA for RAW - SDD -
313 // Check id histograms already created for this Event Specie
315 if(!fDDLModuleMap){CreateTheMap();}
316 if(rawReader->GetType() != 7) return rv; // skips non physical triggers
317 AliDebug(AliQAv1::GetQADebugLevel(),"entering MakeRaws\n");
319 rawReader->Select("ITSSDD");
320 AliITSRawStream *stream=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader);
321 stream->SetDDLModuleMap(fDDLModuleMap);
327 for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
328 for(Int_t iside=0;iside<fgknSide;iside++) {
329 if(fSDDhRawsTask > 7 + index) fAliITSQADataMakerRec->GetRawsData(7 + index +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
330 // 4 because the 2D histos for single events start after the fourth position
340 Int_t coord1, coord2, signal, moduleSDD, activeModule, index1;
343 Int_t prevDDLID = -1;
345 int totalddl=static_cast<int>(fDDLModuleMap->GetNDDLs());
346 Bool_t *ddldata=new Bool_t[totalddl];
347 for(Int_t jddl=0;jddl<totalddl;jddl++){ddldata[jddl]=kFALSE;}
349 while(stream->Next()) {
350 ildcID = rawReader->GetLDCId();
351 iddl = rawReader->GetDDLID();// - fgkDDLIDshift;
352 //printf("----------------------iddl %i\n",iddl);
354 isddmod = fDDLModuleMap->GetModuleNumber(iddl,stream->GetCarlosId());
356 AliDebug(AliQAv1::GetQADebugLevel(),Form("Found module with iddl: %d, stream->GetCarlosId: %d \n",iddl,stream->GetCarlosId()));
359 if(stream->IsCompletedModule()) {
360 AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedModule == KTRUE\n"));
363 if(stream->IsCompletedDDL()) {
365 if ((rawReader->GetDDLID() != prevDDLID)&&(ddldata[iddl])==kFALSE){
366 size += rawReader->GetDataSize();//in bytes
367 prevDDLID = rawReader->GetDDLID();
371 AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedDDL == KTRUE\n"));
375 coord1 = stream->GetCoord1();
376 coord2 = stream->GetCoord2();
377 signal = stream->GetSignal();
379 moduleSDD = isddmod - fgkmodoffset;
381 if(isddmod <fgkmodoffset|| isddmod>fgknSDDmodules+fgkmodoffset-1) {
382 AliDebug(AliQAv1::GetQADebugLevel(),Form( "Module SDD = %d, resetting it to 1 \n",isddmod));
386 AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);
387 Short_t iside = stream->GetChannel();
389 //printf(" \n%i %i %i %i \n ",lay, lad, det,iside );
390 fAliITSQADataMakerRec->GetRawsData( 0 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()] )->Fill(isddmod);
391 if(lay==3) fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);
392 if(lay==4) {fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);}
396 fAliITSQADataMakerRec->GetRawsData(6+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill((stream->GetCarlosId())+0.5*iside -0.5,iddl);
397 // 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));
398 //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));
399 activeModule = moduleSDD;
400 index1 = activeModule * 2 + iside;
403 AliDebug(AliQAv1::GetQADebugLevel(),Form("Wrong index number %d - patched to 0\n",index1));
407 if(fSDDhRawsTask > 7 + index1) {
408 ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(7 + index1 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(coord2, coord1, signal);
409 ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(7 + index1 + 260*2 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(coord2, coord1, signal);
413 if(!(cnt%10000)) AliDebug(AliQAv1::GetQADebugLevel(),Form(" %d raw digits read",cnt));
415 if(fkOnline){((TH1F*)(fAliITSQADataMakerRec->GetRawsData(7 + 260*4 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(size/1024.);//KB
417 AliDebug(AliQAv1::GetQADebugLevel(),Form("Event completed, %d raw digits read",cnt));
422 // AnalyseBNG(); // Analyse Baseline, Noise, Gain
423 // AnalyseINJ(); // Analyse Injectors
431 //____________________________________________________________________________
432 Int_t AliITSQASDDDataMakerRec::InitDigits()
434 // if(!fHistoCalibration)InitCalibrationArray();
435 // Initialization for DIGIT data - SDD -
436 const Bool_t expert = kTRUE ;
437 const Bool_t image = kTRUE ;
439 // fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
440 //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List
441 TH1F* h0=new TH1F("SDD DIGITS Module Pattern","SDD DIGITS Module Pattern",260,239.5,499.5); //hmod
442 h0->GetXaxis()->SetTitle("SDD Module Number");
443 h0->GetYaxis()->SetTitle("# DIGITS");
444 rv = fAliITSQADataMakerRec->Add2DigitsList(h0,fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
446 // printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
447 TH1F* h1=new TH1F("SDD Anode Distribution","DIGITS Anode Distribution",512,-0.5,511.5); //hanocc
448 h1->GetXaxis()->SetTitle("Anode Number");
449 h1->GetYaxis()->SetTitle("# DIGITS");
450 rv = fAliITSQADataMakerRec->Add2DigitsList(h1,1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
452 //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
453 TH1F* h2=new TH1F("SDD Tbin Distribution","DIGITS Tbin Distribution",256,-0.5,255.5); //htbocc
454 h2->GetXaxis()->SetTitle("Tbin Number");
455 h2->GetYaxis()->SetTitle("# DIGITS");
456 rv = fAliITSQADataMakerRec->Add2DigitsList(h2,2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
458 //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
459 TH1F* h3=new TH1F("SDD ADC Counts Distribution","DIGITS ADC Counts Distribution",200,0.,1024.); //hsig
460 h3->GetXaxis()->SetTitle("ADC Value");
461 h3->GetYaxis()->SetTitle("# DIGITS");
462 rv = fAliITSQADataMakerRec->Add2DigitsList(h3,3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
464 //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
465 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Digits histograms booked\n",fSDDhDigitsTask));
469 //____________________________________________________________________________
470 Int_t AliITSQASDDDataMakerRec::MakeDigits(TTree * digits)
473 // Fill QA for DIGIT - SDD -
474 //AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
475 //fITS->SetTreeAddress();
476 //TClonesArray *iITSdigits = fITS->DigitsAddress(1);
481 TBranch *branchD = digits->GetBranch("ITSDigitsSDD");
484 AliError("can't get the branch with the ITS SDD digits !");
487 // Check id histograms already created for this Event Specie
488 // if ( ! fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset) )
489 // rv = InitDigits() ;
491 static TClonesArray statDigits("AliITSdigitSDD");
492 TClonesArray *iITSdigits = &statDigits;
493 branchD->SetAddress(&iITSdigits);
495 for(Int_t i=0; i<260; i++){
497 digits->GetEvent(nmod);
498 Int_t ndigits = iITSdigits->GetEntries();
499 fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nmod,ndigits);
501 for (Int_t idig=0; idig<ndigits; idig++) {
502 AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
503 Int_t iz=dig->GetCoord1(); // cell number z
504 Int_t ix=dig->GetCoord2(); // cell number x
505 Int_t sig=dig->GetSignal();
506 fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(iz);
507 fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(ix);
508 fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(sig);
514 //____________________________________________________________________________
515 Int_t AliITSQASDDDataMakerRec::InitRecPoints()
517 //if(!fHistoCalibration)InitCalibrationArray();
518 //AliInfo("Initialize SDD recpoints histos\n");
519 // Initialization for RECPOINTS - SDD -
520 const Bool_t expert = kTRUE ;
521 const Bool_t image = kTRUE ;
523 // fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
540 //AliInfo(Form("fAliITSQADataMakerRec->GetEventSpecie() %d\n",fAliITSQADataMakerRec->GetEventSpecie()));
541 //AliInfo(Form("fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()] %d\n",fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]));
542 TH1F *h0 = new TH1F("SDDLay3TotCh","Layer 3 total charge",1000/nOnline,-0.5, 499.5); //position number 0
543 //h0->SetBit(TH1::kCanRebin);
544 h0->GetXaxis()->SetTitle("ADC value");
545 h0->GetYaxis()->SetTitle("Entries");
546 rv = fAliITSQADataMakerRec->Add2RecPointsList(h0, 0 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
547 fSDDhRecPointsTask++;
549 TH1F *h1 = new TH1F("SDDLay4TotCh","Layer 4 total charge",1000/nOnline,-0.5, 499.5);//position number 1
550 //h1->SetBit(TH1::kCanRebin);
551 h1->GetXaxis()->SetTitle("ADC value");
552 h1->GetYaxis()->SetTitle("Entries");
553 rv = fAliITSQADataMakerRec->Add2RecPointsList(h1, 1 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
554 fSDDhRecPointsTask++;
557 TH2F *h2 = new TH2F("SDDGlobalCoordDistribYX","YX Global Coord Distrib",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 2
558 h2->GetYaxis()->SetTitle("Y[cm]");
559 h2->GetXaxis()->SetTitle("X[cm]");
560 rv = fAliITSQADataMakerRec->Add2RecPointsList(h2,2+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image
561 fSDDhRecPointsTask++;
563 TH2F *h3 = new TH2F("SDDGlobalCoordDistribRZ","RZ Global Coord Distrib",6400/nOnline3,-32,32,1400/nOnline4,12,26);//position number 3
564 h3->GetYaxis()->SetTitle("R[cm]");
565 h3->GetXaxis()->SetTitle("Z[cm]");
566 rv = fAliITSQADataMakerRec->Add2RecPointsList(h3,3+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image
567 fSDDhRecPointsTask++;
569 TH2F *h4 = new TH2F("SDDGlobalCoordDistribL3PHIZ","#varphi Z Global Coord Distrib L3",4600/nOnline3,-23,23,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 4
570 h4->GetYaxis()->SetTitle("#phi[rad]");
571 h4->GetXaxis()->SetTitle("Z[cm]");
572 rv = fAliITSQADataMakerRec->Add2RecPointsList(h4,4+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
573 fSDDhRecPointsTask++;
575 TH2F *h5 = new TH2F("SDDGlobalCoordDistribL4PHIZ","#varphi Z Global Coord Distrib L4",6200/nOnline6,-31,31,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 5
576 h5->GetYaxis()->SetTitle("#phi[rad]");
577 h5->GetXaxis()->SetTitle("Z[cm]");
578 rv = fAliITSQADataMakerRec->Add2RecPointsList(h5,5+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
579 fSDDhRecPointsTask++;
581 TH1D *h6 = new TH1D("SDDModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); //position number 6
582 h6->GetXaxis()->SetTitle("Module number"); //spd offset = 240
583 h6->GetYaxis()->SetTitle("Entries");
584 rv = fAliITSQADataMakerRec->Add2RecPointsList(h6,6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
585 fSDDhRecPointsTask++;
588 TH2D *h7 = new TH2D("SDDModPatternL3RP","Modules pattern L3 RP",12,0.5,6.5,14,0.5,14.5); //position number 7
589 h7->GetXaxis()->SetTitle("z[#Module L3 ]");
590 h7->GetYaxis()->SetTitle("#varphi[#Ladder L3]");
591 rv = fAliITSQADataMakerRec->Add2RecPointsList(h7,7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
592 fSDDhRecPointsTask++;
594 TH2D *h8 = new TH2D("SDDModPatternL4RP","Modules pattern L4 RP",16,0.5,8.5,22,0.5,22.5); //position number 8
595 h8->GetXaxis()->SetTitle("[#Module L3 ]");
596 h8->GetYaxis()->SetTitle("#varphi[#Ladder L4]");
597 rv = fAliITSQADataMakerRec->Add2RecPointsList(h8,8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
598 fSDDhRecPointsTask++;
600 //------------------------norm--------------------------//
603 TH1D *h9 = new TH1D("SDDModPatternRPNORM","Modules pattern RP NORM",fgknSDDmodules,239.5,499.5); //position number 9
604 h9->GetXaxis()->SetTitle("Module number"); //spd offset = 240
605 h9->GetYaxis()->SetTitle("Entries");
606 rv = fAliITSQADataMakerRec->Add2RecPointsList(h9,9 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
607 fSDDhRecPointsTask++;
610 TH2D *h10 = new TH2D("SDDModPatternL3RPNORM","Modules pattern L3 RP NORM",12,0.5,6.5,14,0.5,14.5); //position number 10
611 h10->GetXaxis()->SetTitle("z[#Module L3 ]");
612 h10->GetYaxis()->SetTitle("#varphi[#Ladder L3]");
613 rv = fAliITSQADataMakerRec->Add2RecPointsList(h10,10 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
614 fSDDhRecPointsTask++;
616 TH2D *h11 = new TH2D("SDDModPatternL4RPNORM","Modules pattern L4 RP NORM",16,0.5,8.5,22,0.5,22.5); //position number 11
617 h11->GetXaxis()->SetTitle("[#Module L3 ]");
618 h10->GetYaxis()->SetTitle("#varphi[#Ladder L4]");
619 rv = fAliITSQADataMakerRec->Add2RecPointsList(h11,11 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
620 fSDDhRecPointsTask++;
622 //--------------------------------------------------------//
624 TH2F *h12 = new TH2F("SDDLocalCoordDistrib","Local Coord Distrib",1000/nOnline,-4,4,1000/nOnline,-4,4);//position number 12
625 h12->GetXaxis()->SetTitle("X local coord, drift, cm");
626 h12->GetYaxis()->SetTitle("Z local coord, anode, cm");
627 rv = fAliITSQADataMakerRec->Add2RecPointsList(h12,12 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
628 fSDDhRecPointsTask++;
630 //AliInfo("Create SDD recpoints histos\n");
632 TH1F *h13 = new TH1F("SDDrdistrib_Layer3" ,"SDD r distribution Layer3" ,100,14.,16.5);//position number 13 (L3)
633 h13->SetBit(TH1::kCanRebin);
634 h13->GetXaxis()->SetTitle("r[cm]");
635 h13->GetXaxis()->CenterTitle();
636 h13->GetYaxis()->SetTitle("Entries");
637 rv = fAliITSQADataMakerRec->Add2RecPointsList(h13,13 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
638 fSDDhRecPointsTask++;
640 TH1F *h14 = new TH1F("SDDrdistrib_Layer4" ,"SDD r distribution Layer4" ,100,23.,25.);// and position number 14 (L4)
641 h14->SetBit(TH1::kCanRebin);
642 h14->GetXaxis()->SetTitle("r[cm]");
643 h14->GetXaxis()->CenterTitle();
644 h14->GetYaxis()->SetTitle("Entries");
645 rv = fAliITSQADataMakerRec->Add2RecPointsList(h14,14 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
646 fSDDhRecPointsTask++;
648 for(Int_t iLay=0; iLay<=1; iLay++){
649 sprintf(hisnam,"SDDphidistrib_Layer%d",iLay+3);
650 TH1F *h15 = new TH1F(hisnam,hisnam,180,-TMath::Pi(),TMath::Pi());//position number 15 (L3) and position number 16 (L4)
651 h15->GetXaxis()->SetTitle("#varphi[rad]");
652 h15->GetXaxis()->CenterTitle();
653 h15->GetYaxis()->SetTitle("Entries");
654 rv = fAliITSQADataMakerRec->Add2RecPointsList(h15,iLay+15+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
655 fSDDhRecPointsTask++;
658 for(Int_t iLay=0; iLay<=1; iLay++){
659 sprintf(hisnam,"SDDdrifttime_Layer%d",iLay+3);
660 TH1F *h17 = new TH1F(hisnam,hisnam,90/nOnline5,-0.5,4499.5);//position number 17 (L3) and position number 18 (L4)
661 h17->SetBit(TH1::kCanRebin);
662 h17->GetXaxis()->SetTitle("drift time[#mus]");
663 h17->GetXaxis()->CenterTitle();
664 h17->GetYaxis()->SetTitle("Entries");
665 rv = fAliITSQADataMakerRec->Add2RecPointsList(h17,iLay+17+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image
666 fSDDhRecPointsTask++;
670 TH2F *h19 = new TH2F("SDDGlobalCoordDistribYXFSE","YX Global Coord Distrib FSE",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 19
671 h19->GetYaxis()->SetTitle("Y[cm]");
672 h19->GetXaxis()->SetTitle("X[cm]");
673 rv = fAliITSQADataMakerRec->Add2RecPointsList(h19,19+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
674 fSDDhRecPointsTask++;
676 TH2F *h20 = new TH2F("SDDGlobalCoordDistribRZFSE","RZ Global Coord Distrib FSE",Int_t(6400/nOnline3),-32,32,1400/nOnline4,12,26);//position number 20
677 h20->GetYaxis()->SetTitle("R[cm]");
678 h20->GetXaxis()->SetTitle("Z[cm]");
679 rv = fAliITSQADataMakerRec->Add2RecPointsList(h20,20+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
680 fSDDhRecPointsTask++;
683 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Recs histograms booked\n",fSDDhRecPointsTask));
688 //____________________________________________________________________________
689 Int_t AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)
691 // Fill QA for RecPoints - SDD -
695 //AliInfo("get the branch with the ITS clusters !\n");
696 AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
697 TClonesArray *recpoints=NULL;
698 if(fkOnline){recpoints = rpcont->FetchClusters(0,clustersTree,fAliITSQADataMakerRec->GetEventNumber());}
699 else{recpoints = rpcont->FetchClusters(0,clustersTree);}
700 if(!rpcont->GetStatusOK() || !rpcont->IsSDDActive()){
701 AliError("can't get SDD clusters !");
706 Float_t cluglo[3]={0.,0.,0.};
708 for(Int_t i=19;i<21;i++){
709 fAliITSQADataMakerRec->GetRecPointsData(i+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
712 Int_t firMod=AliITSgeomTGeo::GetModuleIndex(3,1,1);
713 Int_t lasMod=AliITSgeomTGeo::GetModuleIndex(5,1,1);
714 for(Int_t module=firMod; module<lasMod;module++){
715 //AliInfo(Form("Module %d\n",module));
716 recpoints = rpcont->UncheckedGetClusters(module);
717 npoints += recpoints->GetEntries();
718 //AliInfo(Form("modnumb %d, npoints %d, total points %d\n",module, recpoints->GetEntries(),npoints));
719 //AliITSgeomTGeo::GetModuleId(module, lay, lad, det);
720 //AliInfo(Form("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det));
722 //AliInfo(Form("modnumb %d, entries %d\n",module, recpoints->GetEntries()));
723 for(Int_t j=0;j<recpoints->GetEntries();j++){
724 //AliInfo(Form("modnumb %d, entry %d \n",module, j));
725 AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);
726 Int_t index = recp->GetDetectorIndex();
727 lay=recp->GetLayer();
728 Int_t modnumb=index+AliITSgeomTGeo::GetModuleIndex(lay+1,1,1);
729 //printf("modnumb %i\n",modnumb);
730 AliITSgeomTGeo::GetModuleId(modnumb, lay, lad, det);
731 fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(modnumb);//modpatternrp
732 recp->GetGlobalXYZ(cluglo);
733 Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]);
734 Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);
735 Float_t drifttime=recp->GetDriftTime();
736 fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());//local distribution
737 if(lay==3||lay==4)fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[0],cluglo[1]);//global distribution YX
738 fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],rad);//global distribution rz
740 fAliITSQADataMakerRec->GetRecPointsData(19 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[0],cluglo[1]);//global distribution YX FSE
741 fAliITSQADataMakerRec->GetRecPointsData(20 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],rad);//global distribution rz FSE
743 if(recp->GetLayer() == 2) {
744 fAliITSQADataMakerRec->GetRecPointsData(0 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetQ()) ;//total charge of layer 3
745 //fAliITSQADataMakerRec->GetRecPointsData(7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lad);//lad pattern layer 3
746 Int_t iside=recp->GetDriftSide();
747 //printf("iside =%d\n",iside);
748 fAliITSQADataMakerRec->GetRecPointsData(7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);//mod pattern layer 3
749 fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);//r distribution layer 3
750 fAliITSQADataMakerRec->GetRecPointsData(15 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);// phi distribution layer 3
751 fAliITSQADataMakerRec->GetRecPointsData(4 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],phi);// zphi distribution layer
752 fAliITSQADataMakerRec->GetRecPointsData(17 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(drifttime);// time distribution layer 3
753 } else if(recp->GetLayer() == 3) {
754 fAliITSQADataMakerRec->GetRecPointsData(1 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetQ()) ;//total charge layer 4
755 //fAliITSQADataMakerRec->GetRecPointsData(8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lad);//ladpatternlayer4
756 Int_t iside=recp->GetDriftSide();
757 //printf("iside =%d\n",iside);
758 fAliITSQADataMakerRec->GetRecPointsData(8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);//mod pattern layer 4
759 fAliITSQADataMakerRec->GetRecPointsData(14 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);//r distribution
760 fAliITSQADataMakerRec->GetRecPointsData(16 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);//phi distribution
761 fAliITSQADataMakerRec->GetRecPointsData(5 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],phi);// zphi distribution layer 4
762 fAliITSQADataMakerRec->GetRecPointsData(18 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(drifttime);// time distribution layer 4
769 //_______________________________________________________________
771 Int_t AliITSQASDDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task, Int_t specie)
774 if( task == AliQAv1::kRAWS ){offset=fGenRawsOffset[specie];}
775 else if(task == AliQAv1::kDIGITSR ){offset=fGenDigitsOffset[specie];}
776 else if( task == AliQAv1::kRECPOINTS ){offset=fGenRecPointsOffset[specie];}
780 //_______________________________________________________________
782 void AliITSQASDDDataMakerRec::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset, Int_t specie) {
783 // Returns offset number according to the specified task
784 if( task == AliQAv1::kRAWS ) {fGenRawsOffset[specie]=offset;}
785 else if( task == AliQAv1::kDIGITSR ) {fGenDigitsOffset[specie]=offset;}
786 else if( task == AliQAv1::kRECPOINTS ) {fGenRecPointsOffset[specie]=offset;}
789 //_______________________________________________________________
791 Int_t AliITSQASDDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task)
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");}
802 //_______________________________________________________________
805 void AliITSQASDDDataMakerRec::CreateTheMap()
808 AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
809 Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
811 AliError("Calibration object retrieval failed! SDD will not be processed");
812 fDDLModuleMap = NULL;
815 fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
816 if(!cacheStatus)ddlMapSDD->SetObject(NULL);
817 ddlMapSDD->SetOwner(kTRUE);
818 if(!cacheStatus){ delete ddlMapSDD;}
819 AliInfo("DDL Map Created\n ");
822 //_______________________________________________________________
825 void AliITSQASDDDataMakerRec::CreateTheCalibration()
828 AliCDBEntry *calibSDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
829 Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
832 AliError("Calibration object retrieval failed! SDD will not be processed");
833 fCalibration = NULL;;
835 fCalibration = (TObjArray *)calibSDD->GetObject();
837 if(!cacheStatus)calibSDD->SetObject(NULL);
838 calibSDD->SetOwner(kTRUE);
839 if(!cacheStatus){delete calibSDD;}
841 AliITSCalibrationSDD * cal=NULL;
842 for(Int_t imod=0;imod<fgknSDDmodules;imod++)
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;
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);
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
866 for(Int_t i=0;i<8;i++) //check on bad chips in good modules
869 if(cal->IsChipBad(i)){
870 if(i<4){badmodhisto2side0+=64;}
871 if(i>=4){badmodhisto2side1+=64;}
875 if(cal->IsChipBad(i)){
876 if(i<4){badmodhisto3side0+=64;}
877 if(i>=4){badmodhisto3side1+=64;}
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
885 if(ic<4) badmodhisto2side0++;
886 else if(ic>=4)badmodhisto2side1++;
889 if(ic<4) badmodhisto3side0++;
890 else if(ic>=4)badmodhisto3side1++;
892 }//end if chip module channel
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);
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);
912 }//end else bad module
917 //____________________________________________________________________
919 void AliITSQASDDDataMakerRec::InitCalibrationArray()
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);
926 fHistoCalibration = new TObjArray(3);
927 fHistoCalibration->AddAtAndExpand(pattern1,0);
928 fHistoCalibration->AddAtAndExpand(patternl3,1);
929 fHistoCalibration->AddAtAndExpand(patternl4,2);
931 // printf("Calibration Histograms created!\n");
934 //____________________________________________________________________
936 void AliITSQASDDDataMakerRec::ResetDetector(AliQAv1::TASKINDEX_t task)
939 AliInfo(Form("Reset detector in SDD called for task index %i", task));
940 if(task== AliQAv1::kRAWS ){
945 ((TH1D*)(fHistoCalibration->At(0)))->Reset();
946 ((TH2D*)(fHistoCalibration->At(1)))->Reset();
947 ((TH2D*)(fHistoCalibration->At(2)))->Reset();