]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerRec.cxx
Added QA for digits during reconstruction (Yves)
[u/mrichter/AliRoot.git] / ITS / AliITSQASDDDataMakerRec.cxx
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
24 //  M.Siciliano Aug 2008 QA RecPoints and HLT mode
25 //  INFN Torino
26
27 // --- ROOT system ---
28
29 #include <TProfile2D.h>
30 #include <TH2D.h>
31 #include <TBranch.h>
32 #include <TTree.h>
33 #include <TGaxis.h>
34 #include <TMath.h>
35 #include <TDirectory.h>
36 #include <TSystem.h>
37 // --- Standard library ---
38
39 // --- AliRoot header files ---
40 #include "AliITSQASDDDataMakerRec.h"
41 #include "AliLog.h"
42 #include "AliQAv1.h"
43 #include "AliQAChecker.h"
44 #include "AliRawReader.h"
45 #include "AliITSRawStream.h"
46 #include "AliITSRawStreamSDD.h"
47 #include "AliITSRawStreamSDDCompressed.h"
48 #include "AliITSDetTypeRec.h"
49 #include "AliITSDigit.h"
50 #include "AliITSRecPoint.h"
51 #include "AliITSgeomTGeo.h"
52 #include "AliITSHLTforSDD.h"
53 #include "AliCDBManager.h"
54 #include "AliCDBStorage.h"
55 #include "AliCDBEntry.h"
56 #include "Riostream.h"
57
58 ClassImp(AliITSQASDDDataMakerRec)
59
60 //____________________________________________________________________________ 
61 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) :
62 TObject(),
63 fAliITSQADataMakerRec(aliITSQADataMakerRec),
64 fkOnline(kMode),
65 fLDC(ldc),
66 fSDDhRawsTask(0),
67 fSDDhDigitsTask(0),
68 fSDDhRecPointsTask(0),
69 fGenRawsOffset(0),
70 fGenDigitsOffset(0),
71 fGenRecPointsOffset(0),
72 fTimeBinSize(1),
73 fDDLModuleMap(0),
74 fHLTMode(0),
75 fHLTSDD(0)
76 {
77   //ctor used to discriminate OnLine-Offline analysis
78   if(fLDC < 0 || fLDC > 4) {
79         AliError("Error: LDC number out of range; return\n");
80   }
81   if(!fkOnline){AliInfo("Offline mode: HLT set from AliITSDetTypeRec for SDD\n");}
82   else
83     if(fkOnline){
84       AliInfo("Online mode: HLT set from environment for SDD\n");
85       SetHLTModeFromEnvironment();
86     }
87   //fDDLModuleMap=NULL;
88 }
89
90 //____________________________________________________________________________ 
91 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) :
92 TObject(),
93 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
94 fkOnline(qadm.fkOnline),
95 fLDC(qadm.fLDC),
96 fSDDhRawsTask(qadm.fSDDhRawsTask),
97 fSDDhDigitsTask(qadm.fSDDhDigitsTask),
98 fSDDhRecPointsTask(qadm.fSDDhRecPointsTask),
99 fGenRawsOffset(qadm.fGenRawsOffset),
100 fGenDigitsOffset(qadm.fGenDigitsOffset),
101 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
102 fTimeBinSize(1),
103 fDDLModuleMap(0),
104 fHLTMode(qadm.fHLTMode),
105 fHLTSDD( qadm.fHLTSDD)
106 {
107   //copy ctor 
108   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; 
109   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
110   fDDLModuleMap=NULL;
111 }
112
113 //____________________________________________________________________________ 
114 AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){
115   // destructor
116   //    if(fDDLModuleMap) delete fDDLModuleMap;
117 }
118 //__________________________________________________________________
119 AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac )
120 {
121   // Equal operator.
122   this->~AliITSQASDDDataMakerRec();
123   new(this) AliITSQASDDDataMakerRec(qac);
124   return *this;
125 }
126
127 //____________________________________________________________________________ 
128 void AliITSQASDDDataMakerRec::StartOfDetectorCycle()
129 {
130   //Detector specific actions at start of cycle
131   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SDD Cycle\n");
132 }
133
134 //____________________________________________________________________________ 
135 void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
136 {
137   // launch the QA checking
138   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); 
139 }
140
141 //____________________________________________________________________________ 
142 void AliITSQASDDDataMakerRec::InitRaws()
143
144   // Initialization for RAW data - SDD -
145   const Bool_t expert   = kTRUE ; 
146   const Bool_t saveCorr = kTRUE ; 
147   const Bool_t image    = kTRUE ; 
148   
149   fGenRawsOffset = (fAliITSQADataMakerRec->fRawsQAList[AliRecoParam::kDefault])->GetEntries();
150   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
151   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
152   if(!ddlMapSDD)
153     {
154       AliError("Calibration object retrieval failed! SDD will not be processed");
155       fDDLModuleMap = NULL;
156       return;
157     }
158   fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
159   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
160   ddlMapSDD->SetOwner(kTRUE);
161   if(!cacheStatus)
162     {
163       delete ddlMapSDD;
164     }
165   
166   if(fkOnline==kFALSE){
167     AliInfo("Offline mode: HLTforSDDobject used \n");
168     AliCDBEntry *hltforSDD = AliCDBManager::Instance()->Get("ITS/Calib/HLTforSDD");
169     if(!hltforSDD){
170       AliError("Calibration object retrieval failed! SDD will not be processed");    
171       fHLTSDD=NULL;
172       return;
173     }  
174     fHLTSDD = (AliITSHLTforSDD*)hltforSDD->GetObject();
175     if(!cacheStatus)hltforSDD->SetObject(NULL);
176     hltforSDD->SetOwner(kTRUE);
177     if(!cacheStatus)
178       {
179         delete hltforSDD;
180       }
181   }
182   Int_t lay, lad, det;
183   Int_t indexlast = 0;
184   Int_t index1 = 0;
185
186   if(fkOnline) 
187     {
188       AliInfo("Book Online Histograms for SDD\n");
189     }
190   else 
191     {
192       AliInfo("Book Offline Histograms for SDD\n ");
193     }
194   TH1D *h0 = new TH1D("SDDModPattern","HW Modules pattern",fgknSDDmodules,239.5,499.5); //0
195   h0->GetXaxis()->SetTitle("Module Number");
196   h0->GetYaxis()->SetTitle("Counts");
197   fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h0)),0+fGenRawsOffset, expert, !image, !saveCorr);
198   delete h0;
199   fSDDhRawsTask++;
200   
201   //zPhi distribution using ladder and modules numbers
202   TH2D *hphil3 = new TH2D("SDDphizL3","SDD #varphiz Layer3 ",6,0.5,6.5,14,0.5,14.5);
203   hphil3->GetXaxis()->SetTitle("z[#Module L3 ]");
204   hphil3->GetYaxis()->SetTitle("#varphi[#Ladder L3]");
205   fAliITSQADataMakerRec->Add2RawsList((new TH2D(*hphil3)),1+fGenRawsOffset, !expert, image, saveCorr); 
206   delete hphil3;
207   fSDDhRawsTask++;
208   
209   TH2D *hphil4 = new TH2D("SDDphizL4","SDD #varphiz Layer4 ",8,0.5,8.5,22,0.5,22.5); 
210   hphil4->GetXaxis()->SetTitle("z[#Module L4]");
211   hphil4->GetYaxis()->SetTitle("#varphi[#Ladder L4]");
212   fAliITSQADataMakerRec->Add2RawsList((new TH2D(*hphil4)),2+fGenRawsOffset, !expert, image, saveCorr); 
213   delete hphil4;
214   fSDDhRawsTask++;
215   
216
217   if(fkOnline) 
218     {
219
220       //DDL Pattern 
221       TH2D *hddl = new TH2D("SDDDDLPattern","SDD DDL Pattern ",24,-0.5,23.5,24,-0.5,23.5); 
222       hddl->GetXaxis()->SetTitle("Channel");
223       hddl->GetYaxis()->SetTitle("#DDL");
224       fAliITSQADataMakerRec->Add2RawsList((new TH2D(*hddl)),3+fGenRawsOffset, expert, !image, !saveCorr);
225       delete hddl;
226       fSDDhRawsTask++;
227       Int_t indexlast1 = 0;
228   
229       fTimeBinSize = 4;
230       indexlast = 0;
231       index1 = 0;
232       indexlast1 = fSDDhRawsTask;
233       char *hname[3];
234       for(Int_t i=0; i<3; i++) hname[i]= new char[50];
235       for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
236         for(Int_t iside=0;iside<fgknSide;iside++){
237           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
238           sprintf(hname[0],"SDDchargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside);
239           sprintf(hname[1],"SDDChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside);
240           sprintf(hname[2],"SDDhmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);
241           TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
242           fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin");
243           fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode");
244           fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMapFSE)),indexlast1 + index1 + fGenRawsOffset, expert, !image, !saveCorr);
245           delete fModuleChargeMapFSE;
246           
247           fSDDhRawsTask++;
248           index1++;      
249         }
250       }
251       
252       for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
253         for(Int_t iside=0;iside<fgknSide;iside++){
254           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
255           sprintf(hname[0],"SDDchargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
256           sprintf(hname[1],"SDDChargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
257           TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
258           fModuleChargeMap->GetXaxis()->SetTitle("Time Bin");
259           fModuleChargeMap->GetYaxis()->SetTitle("Anode");
260           fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMap)),indexlast1 + index1 + fGenRawsOffset, expert, !image, !saveCorr);
261           delete fModuleChargeMap;
262           
263           fSDDhRawsTask++;
264           index1++;      
265         }
266       }
267       
268     }  // kONLINE
269   
270   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Raws histograms booked\n",fSDDhRawsTask));
271 }
272
273
274 //____________________________________________________________________________
275 void AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader)
276
277   // Fill QA for RAW - SDD -
278   
279   if(!fDDLModuleMap){
280     AliError("SDD DDL module map not available - skipping SDD QA");
281     return;
282   }
283   if(rawReader->GetType() != 7) return;  // skips non physical triggers
284   AliDebug(AliQAv1::GetQADebugLevel(),"entering MakeRaws\n");                 
285   
286   rawReader->Reset();       
287   AliITSRawStream *stream;
288   
289   if(fkOnline==kTRUE)
290     {
291       if(GetHLTMode()==kTRUE)
292         {
293           //AliInfo("Online  mode: HLT C compressed mode used for SDD\n");
294           stream = new AliITSRawStreamSDDCompressed(rawReader); }
295       else{ 
296         //AliInfo("Online  mode: HLT A mode used for SDD\n");
297         stream = new AliITSRawStreamSDD(rawReader);}     
298     }
299   else 
300     {
301       if(fHLTSDD->IsHLTmodeC()==kTRUE){
302         //AliInfo("Offline  mode: HLT C compressed mode used for SDD\n");
303         stream = new AliITSRawStreamSDDCompressed(rawReader);
304       }else 
305         {
306           //AliInfo("Offline  mode: HLT A mode used for SDD\n");
307           stream = new AliITSRawStreamSDD(rawReader);
308         }
309     }
310   
311   //ckeck on HLT mode
312  
313   //  AliITSRawStreamSDD s(rawReader); 
314   stream->SetDDLModuleMap(fDDLModuleMap);
315   
316   Int_t lay, lad, det; 
317   
318   Int_t index=0;
319   if(fkOnline) {
320     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
321       for(Int_t iside=0;iside<fgknSide;iside++) {
322         if(fSDDhRawsTask > 4 + index) fAliITSQADataMakerRec->GetRawsData(4 + index +fGenRawsOffset)->Reset();   
323         // 4  because the 2D histos for single events start after the fourth position
324         index++;
325       }
326     }
327   }
328   
329   Int_t cnt = 0;
330   Int_t ildcID = -1;
331   Int_t iddl = -1;
332   Int_t isddmod = -1;
333   Int_t coord1, coord2, signal, moduleSDD, activeModule, index1; 
334   
335   while(stream->Next()) {
336     ildcID = rawReader->GetLDCId();
337     iddl = rawReader->GetDDLID() - fgkDDLIDshift;
338     
339     isddmod = fDDLModuleMap->GetModuleNumber(iddl,stream->GetCarlosId());
340     if(isddmod==-1){
341       AliDebug(AliQAv1::GetQADebugLevel(),Form("Found module with iddl: %d, stream->GetCarlosId: %d \n",iddl,stream->GetCarlosId()));
342       continue;
343     }
344     if(stream->IsCompletedModule()) {
345       AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedModule == KTRUE\n"));
346       continue;
347     } 
348     if(stream->IsCompletedDDL()) {
349       AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedDDL == KTRUE\n"));
350       continue;
351     } 
352     
353     coord1 = stream->GetCoord1();
354     coord2 = stream->GetCoord2();
355     signal = stream->GetSignal();
356     
357     moduleSDD = isddmod - fgkmodoffset;
358     
359     if(isddmod <fgkmodoffset|| isddmod>fgknSDDmodules+fgkmodoffset-1) {
360       AliDebug(AliQAv1::GetQADebugLevel(),Form( "Module SDD = %d, resetting it to 1 \n",isddmod));
361       isddmod = 1;
362     }
363     
364     AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);
365
366     
367     fAliITSQADataMakerRec->GetRawsData( 0 + fGenRawsOffset )->Fill(isddmod);   
368     
369     if(lay==3)    fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset)->Fill(det,lad); 
370     if(lay==4) { 
371       fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset)->Fill(det,lad);}  
372     
373     Short_t iside = stream->GetChannel();
374     
375     fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset)->Fill(2*(stream->GetCarlosId())+iside,iddl);
376     
377
378     if(fkOnline) {
379
380       activeModule = moduleSDD;
381       index1 = activeModule * 2 + iside;
382       
383       if(index1<0){
384         AliDebug(AliQAv1::GetQADebugLevel(),Form("Wrong index number %d - patched to 0\n",index1));
385         index1 = 0;
386       }      
387       fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset)->Fill(2*(stream->GetCarlosId())+iside,iddl);
388       if(fSDDhRawsTask > 4 + index1) {                                  
389         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(4 + index1 +fGenRawsOffset)))->Fill(coord2, coord1, signal);     
390         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(4 + index1 + 260*2 +fGenRawsOffset)))->Fill(coord2, coord1, signal); 
391       }
392     }
393     cnt++;
394     if(!(cnt%10000)) AliDebug(AliQAv1::GetQADebugLevel(),Form(" %d raw digits read",cnt));
395   }
396   AliDebug(AliQAv1::GetQADebugLevel(),Form("Event completed, %d raw digits read",cnt)); 
397   delete stream;
398   stream = NULL; 
399 }
400
401 //____________________________________________________________________________ 
402 void AliITSQASDDDataMakerRec::InitDigits()
403
404   // Initialization for DIGIT data - SDD -  
405   const Bool_t expert   = kTRUE ; 
406   const Bool_t image    = kTRUE ;
407   
408   fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
409   //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List
410   TH1F* h0=new TH1F("SDD DIGITS Module Pattern","SDD DIGITS Module Pattern",260,239.5,499.5);       //hmod
411   h0->GetXaxis()->SetTitle("SDD Module Number");
412   h0->GetYaxis()->SetTitle("# DIGITS");
413   fAliITSQADataMakerRec->Add2DigitsList(h0,fGenDigitsOffset, !expert, image);
414   fSDDhDigitsTask ++;
415   TH1F* h1=new TH1F("SDD Anode Distribution","DIGITS Anode Distribution",512,-0.5,511.5);      //hanocc
416   h1->GetXaxis()->SetTitle("Anode Number");
417   h1->GetYaxis()->SetTitle("# DIGITS");
418   fAliITSQADataMakerRec->Add2DigitsList(h1,1+fGenDigitsOffset, !expert, image);
419   fSDDhDigitsTask ++;
420   TH1F* h2=new TH1F("SDD Tbin Distribution","DIGITS Tbin Distribution",256,-0.5,255.5);      //htbocc
421   h2->GetXaxis()->SetTitle("Tbin Number");
422   h2->GetYaxis()->SetTitle("# DIGITS");
423   fAliITSQADataMakerRec->Add2DigitsList(h2,2+fGenDigitsOffset, !expert, image);
424   fSDDhDigitsTask ++;
425   TH1F* h3=new TH1F("SDD ADC Counts Distribution","DIGITS ADC Counts Distribution",200,0.,1024.);          //hsig
426   h3->GetXaxis()->SetTitle("ADC Value");
427   h3->GetYaxis()->SetTitle("# DIGITS");
428   fAliITSQADataMakerRec->Add2DigitsList(h3,3+fGenDigitsOffset, !expert, image);
429   fSDDhDigitsTask ++;
430   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Digits histograms booked\n",fSDDhDigitsTask));
431 }
432
433 //____________________________________________________________________________
434 void AliITSQASDDDataMakerRec::MakeDigits(TTree * digits)
435
436   
437   // Fill QA for DIGIT - SDD -
438 //  AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");
439 //  fITS->SetTreeAddress();
440 //  TClonesArray *iITSdigits  = fITS->DigitsAddress(1);
441   TBranch *branchD = digits->GetBranch("ITS");
442   if (!branchD) { 
443     AliError("can't get the branch with the ITS digits !");
444     return;
445   }
446   static TClonesArray statDigits("AliITSDigit");
447   TClonesArray *iITSdigits = &statDigits;
448   branchD->SetAddress(&iITSdigits);
449   for(Int_t i=0; i<260; i++){
450     Int_t nmod=i+240;
451     digits->GetEvent(nmod);
452     Int_t ndigits = iITSdigits->GetEntries();
453     fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset)->Fill(nmod,ndigits);
454     for (Int_t idig=0; idig<ndigits; idig++) {
455       AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
456       Int_t iz=dig->GetCoord1();  // cell number z
457       Int_t ix=dig->GetCoord2();  // cell number x
458       Int_t sig=dig->GetSignal();
459       fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset)->Fill(iz);
460       fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset)->Fill(ix);
461       fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset)->Fill(sig);
462     }
463   }
464 }
465
466 //____________________________________________________________________________ 
467 void AliITSQASDDDataMakerRec::InitRecPoints()
468 {
469   // Initialization for RECPOINTS - SDD -
470   const Bool_t expert   = kTRUE ; 
471   const Bool_t image    = kTRUE ; 
472   
473   fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
474
475   Int_t nOnline=1;
476   Int_t  nOnline2=1;
477   Int_t  nOnline3=1; 
478   Int_t  nOnline4=1;
479   if(fkOnline)
480     {
481       nOnline=4;
482       nOnline2=28;
483       nOnline3=64;
484       nOnline4=14;
485     }
486
487   
488   TH1F *h0 = new TH1F("SDDLay3TotCh","Layer 3 total charge",1000/nOnline,-0.5, 499.5); //position number 0
489   h0->GetXaxis()->SetTitle("ADC value");
490   h0->GetYaxis()->SetTitle("Entries");
491   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h0)), 0 +fGenRecPointsOffset, !expert, image);
492   delete h0;
493   fSDDhRecPointsTask++;
494  
495   TH1F *h1 = new TH1F("SDDLay4TotCh","Layer 4 total charge",1000/nOnline,-0.5, 499.5);//position number 1
496   h1->GetXaxis()->SetTitle("ADC value");
497   h1->GetYaxis()->SetTitle("Entries");
498   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h1)), 1 +fGenRecPointsOffset, !expert, image);
499   delete h1;
500   fSDDhRecPointsTask++;
501
502   char hisnam[50];
503   TH2F *h2 = new TH2F("SDDGlobalCoordDistribYX","YX Global Coord Distrib",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 2
504   h2->GetYaxis()->SetTitle("Y[cm]");
505   h2->GetXaxis()->SetTitle("X[cm]");
506   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h2)),2+fGenRecPointsOffset, expert, !image);
507   delete h2;
508   fSDDhRecPointsTask++;
509
510   TH2F *h3 = new TH2F("SDDGlobalCoordDistribRZ","RZ Global Coord Distrib",6400/nOnline3,-32,32,1400/nOnline4,12,26);//position number 3
511   h3->GetYaxis()->SetTitle("R[cm]");
512   h3->GetXaxis()->SetTitle("Z[cm]");
513   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h3)),3+fGenRecPointsOffset, expert, !image);
514   delete h3;
515   fSDDhRecPointsTask++;
516   
517   TH2F *h4 = new TH2F("SDDGlobalCoordDistribL3PHIZ","#varphi Z Global Coord Distrib L3",6400/nOnline3,-32,32,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 4
518   h4->GetYaxis()->SetTitle("#phi[rad]");
519   h4->GetXaxis()->SetTitle("Z[cm]");
520   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h4)),4+fGenRecPointsOffset, !expert, image);
521   delete h4;
522   fSDDhRecPointsTask++;
523
524   TH2F *h5 = new TH2F("SDDGlobalCoordDistribL4PHIZ","#varphi Z Global Coord Distrib L4",6400/nOnline3,-32,32,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 5
525   h5->GetYaxis()->SetTitle("#phi[rad]");
526   h5->GetXaxis()->SetTitle("Z[cm]");
527   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h5)),5+fGenRecPointsOffset, !expert, image);
528   delete h5;
529   fSDDhRecPointsTask++;
530   
531   TH1F *h6 = new TH1F("SDDModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); //position number 6
532   h6->GetXaxis()->SetTitle("Module number"); //spd offset = 240
533   h6->GetYaxis()->SetTitle("Entries");
534   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h6)),6 +fGenRecPointsOffset, expert, !image);
535   delete h6;
536   fSDDhRecPointsTask++;
537   TH1F *h7 = new TH1F("SDDLadPatternL3RP","Ladder pattern L3 RP",14,0.5,14.5);  //position number 7
538   h7->GetXaxis()->SetTitle("Ladder #, Layer 3");
539   h7->GetYaxis()->SetTitle("Entries");
540   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h7)),7 +fGenRecPointsOffset, expert, !image);
541   delete h7;
542   fSDDhRecPointsTask++;
543   TH1F *h8 = new TH1F("SDDLadPatternL4RP","Ladder pattern L4 RP",22,0.5,22.5); //position number 8
544   h8->GetXaxis()->SetTitle("Ladder #, Layer 4");
545   h8->GetYaxis()->SetTitle("Entries");
546   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h8)),8 +fGenRecPointsOffset, expert, !image);
547   delete h8;
548   fSDDhRecPointsTask++;
549   TH2F *h9 = new TH2F("SDDLocalCoordDistrib","Local Coord Distrib",1000/nOnline,-4,4,1000/nOnline,-4,4);//position number 9
550   h9->GetXaxis()->SetTitle("X local coord, drift, cm");
551   h9->GetYaxis()->SetTitle("Z local coord, anode, cm");
552   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h9)),9 +fGenRecPointsOffset, expert, !image);
553   delete h9;
554   fSDDhRecPointsTask++;
555
556
557     TH1F *h10 = new TH1F("SDDrdistrib_Layer3" ,"SDD r distribution Layer3" ,100,14.,18.);//position number 10 (L3)
558     h10->GetXaxis()->SetTitle("r[cm]");
559     h10->GetXaxis()->CenterTitle();
560     h10->GetYaxis()->SetTitle("Entries");
561     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h10)),10 +fGenRecPointsOffset, expert, !image);
562     delete h10;
563     fSDDhRecPointsTask++;
564
565     TH1F *h11 = new TH1F("SDDrdistrib_Layer4" ,"SDD r distribution Layer4" ,100,22.,26.);// and position number 11 (L4)
566     h11->GetXaxis()->SetTitle("r[cm]");
567     h11->GetXaxis()->CenterTitle();
568     h11->GetYaxis()->SetTitle("Entries");
569     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h11)),11 +fGenRecPointsOffset, expert, !image);
570     delete h11;
571     fSDDhRecPointsTask++;
572
573   for(Int_t iLay=0; iLay<=1; iLay++){
574     sprintf(hisnam,"SDDphidistrib_Layer%d",iLay+3);
575     TH1F *h12 = new TH1F(hisnam,hisnam,180,-TMath::Pi(),TMath::Pi());//position number 12 (L3) and position number 13 (L4)
576     h12->GetXaxis()->SetTitle("#varphi[rad]");
577     h12->GetXaxis()->CenterTitle();
578     h12->GetYaxis()->SetTitle("Entries");
579     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h12)),iLay+12+fGenRecPointsOffset, expert, !image);
580     delete h12;
581     fSDDhRecPointsTask++;
582   }
583
584   if(fkOnline)
585     {
586       TH2F *h14 = new TH2F("SDDGlobalCoordDistribYXFSE","YX Global Coord Distrib FSE",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 14
587       h14->GetYaxis()->SetTitle("Y[cm]");
588       h14->GetXaxis()->SetTitle("X[cm]");
589       fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h14)),14+fGenRecPointsOffset, expert, !image);
590       delete h14;
591       fSDDhRecPointsTask++;
592       
593       TH2F *h15 = new TH2F("SDDGlobalCoordDistribRZFSE","RZ Global Coord Distrib FSE",Int_t(6400/nOnline3),-32,32,1400/nOnline4,12,26);//position number 15
594       h15->GetYaxis()->SetTitle("R[cm]");
595       h15->GetXaxis()->SetTitle("Z[cm]");
596       fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h15)),15+fGenRecPointsOffset, expert, !image);
597       delete h15;
598       fSDDhRecPointsTask++;
599       
600     }//online
601
602   //printf("%d SDD Recs histograms booked\n",fSDDhRecPointsTask);
603
604
605   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Recs histograms booked\n",fSDDhRecPointsTask));
606
607
608 }
609
610 //____________________________________________________________________________ 
611 void AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)
612 {
613
614
615  // Fill QA for RecPoints - SDD -
616   Int_t lay, lad, det; 
617   TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
618   if (!branchRecP) {
619     AliError("can't get the branch with the ITS clusters !");
620     return;
621   }
622
623   static TClonesArray statRecpoints("AliITSRecPoint") ;
624   TClonesArray *recpoints = &statRecpoints;
625   branchRecP->SetAddress(&recpoints);
626   Int_t npoints = 0;      
627   Float_t cluglo[3]={0.,0.,0.}; 
628   if(fkOnline)
629     {
630       for(Int_t i=14;i<16;i++)
631         {
632           fAliITSQADataMakerRec->GetRecPointsData(i+fGenRecPointsOffset)->Reset();
633         }
634     }
635   for(Int_t module=0; module<clustersTree->GetEntries();module++){
636     branchRecP->GetEvent(module);
637     npoints += recpoints->GetEntries();
638     AliITSgeomTGeo::GetModuleId(module, lay, lad, det);
639     //printf("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det);
640     for(Int_t j=0;j<recpoints->GetEntries();j++){
641       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);    
642       fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset)->Fill(module);//modpatternrp
643       recp->GetGlobalXYZ(cluglo);
644       Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); 
645       Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);
646       fAliITSQADataMakerRec->GetRecPointsData(9 +fGenRecPointsOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());//local distribution
647       fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset)->Fill(cluglo[0],cluglo[1]);//global distribution YX
648       fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset)->Fill(cluglo[2],rad);//global distribution rz
649       if(fkOnline)
650         {
651           fAliITSQADataMakerRec->GetRecPointsData(14 +fGenRecPointsOffset)->Fill(cluglo[0],cluglo[1]);//global distribution YX FSE
652           fAliITSQADataMakerRec->GetRecPointsData(15 +fGenRecPointsOffset)->Fill(cluglo[2],rad);//global distribution rz FSE
653         }
654       if(recp->GetLayer() ==2) {
655         fAliITSQADataMakerRec->GetRecPointsData(0  +fGenRecPointsOffset)->Fill(recp->GetQ()) ;//total charge of layer 3
656         fAliITSQADataMakerRec->GetRecPointsData(7  +fGenRecPointsOffset)->Fill(lad);//lad pattern layer 3
657         fAliITSQADataMakerRec->GetRecPointsData(10 +fGenRecPointsOffset)->Fill(rad);//r distribution layer 3
658         fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset)->Fill(phi);// phi distribution layer 3
659         fAliITSQADataMakerRec->GetRecPointsData(4  +fGenRecPointsOffset)->Fill(cluglo[2],phi);// phi distribution layer 3
660       }
661       else if(recp->GetLayer() ==3) {
662         fAliITSQADataMakerRec->GetRecPointsData(1  +fGenRecPointsOffset)->Fill(recp->GetQ()) ;//total charge layer 4
663         fAliITSQADataMakerRec->GetRecPointsData(8  +fGenRecPointsOffset)->Fill(lad);//ladpatternlayer4
664         fAliITSQADataMakerRec->GetRecPointsData(11 +fGenRecPointsOffset)->Fill(rad);//r distribution
665         fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset)->Fill(phi);//phi distribution
666         fAliITSQADataMakerRec->GetRecPointsData(5  +fGenRecPointsOffset)->Fill(cluglo[2],phi);// phi distribution layer 4
667       }
668     }
669   }
670   statRecpoints.Clear();
671
672 }
673
674 //_______________________________________________________________
675
676 void AliITSQASDDDataMakerRec::SetHLTModeFromEnvironment()
677 {
678
679    Int_t  hltmode= ::atoi(gSystem->Getenv("HLT_MODE"));
680
681    if(hltmode==1)
682      {
683        AliInfo("Online mode: HLT mode A selected from environment for SDD\n");
684        SetHLTMode(kFALSE);
685      }
686    else
687      if(hltmode==2)
688        {
689        AliInfo("Online mode: HLT mode C compressed selected from environment for SDD\n");
690        SetHLTMode(kTRUE);
691        }
692 }
693
694
695 //_______________________________________________________________
696
697 Int_t AliITSQASDDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task)
698 {
699   Int_t offset=0;
700   if( task == AliQAv1::kRAWS )
701     {
702       offset=fGenRawsOffset;  
703     }
704   else
705     if( task == AliQAv1::kRECPOINTS )
706       {
707         offset=fGenRecPointsOffset;   
708       }
709     else AliInfo("No task has been selected. Offset set to zero.\n");
710   return offset;
711 }
712
713 //_______________________________________________________________
714
715 Int_t AliITSQASDDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task)
716 {
717
718   Int_t histotot=0;
719
720   if( task == AliQAv1::kRAWS )
721     {
722       histotot=fSDDhRawsTask ;  
723     }
724   else
725     if( task == AliQAv1::kRECPOINTS )
726       {
727         histotot=fSDDhRecPointsTask;   
728       }
729     else AliInfo("No task has been selected. TaskHisto set to zero.\n");
730   return histotot;
731 }