]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASSDDataMakerRec.cxx
147b0b1d59d206c5e730afe6b681092a174d14ba
[u/mrichter/AliRoot.git] / ITS / AliITSQASSDDataMakerRec.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 /* $Id$    */
16 //  *************************************************************
17 //  Checks the quality assurance 
18 //  by comparing with reference data
19 //  contained in a DB
20 //  -------------------------------------------------------------
21 //  W. Ferrarese + P. Cerello Feb 2008
22 //  INFN Torino
23
24 // --- ROOT system ---
25 #include <TH2D.h>
26 #include <TTree.h>
27 #include <TMath.h>
28 #include <TString.h>
29 #include <TSystem.h>
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliITSQASSDDataMakerRec.h"
35 #include "AliQADataMakerRec.h"
36 #include "AliLog.h"
37 #include "AliQAv1.h"
38 #include "AliQAChecker.h"
39 #include "AliRawReader.h"
40 #include "AliRawReaderRoot.h"
41 #include "AliITSRawStreamSSD.h"
42 #include "AliITSgeomTGeo.h"
43 #include "AliRawEventHeaderBase.h"
44 #include "AliITSRecPoint.h"
45 #include "AliITSdigitSSD.h"
46 #include "AliITSBadChannelsSSDv2.h"
47
48 #include "AliCDBManager.h"
49 #include "AliCDBEntry.h"
50
51 ClassImp(AliITSQASSDDataMakerRec)
52
53 AliITSQASSDDataMakerRec::AliITSQASSDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Int_t ldc) :
54 TObject(),
55 fAliITSQADataMakerRec(aliITSQADataMakerRec),
56 fSSDEvent(0),
57 fSSDEventPerCycle(0),
58 fkOnline(kMode),
59 fLDC(ldc),
60 fSSDRawsOffset(0), fSSDRawsDAOffset(0),
61 fSSDRawsCommonLevelOffset(0),
62 fSSDhRawsTask(0),
63 fSSDhDigitsTask(0),
64 fSSDhRecPointsTask(0),
65 fGenRawsOffset(0),
66 fGenDigitsOffset(0),
67 fGenRecPointsOffset(0),
68 fCDBManager(0) {
69   // Default constructor   
70   //initilize the raw signal vs strip number histograms
71   if(fkOnline) {
72     fCDBManager = AliCDBManager::Instance();
73     //fCDBManager->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
74     fCDBManager->SetDefaultStorage(gSystem->Getenv("AMORE_CDB_URI"));
75     Int_t runNumber = atoi(gSystem->Getenv("DATE_RUN_NUMBER"));
76     if(!runNumber) 
77       AliWarning("DATE_RUN_NUMBER not defined!!!\n");
78     
79     fCDBManager->SetRun(runNumber);
80     AliCDBEntry *geomGRP = fCDBManager->Get("GRP/Geometry/Data");
81     if(!geomGRP) AliWarning("GRP geometry not found!!!\n");
82     
83     Int_t gLayer = 0,gLadder = 0, gModule = 0;
84     Int_t gHistCounter = 0;
85     TString gTitle; 
86
87     for(Int_t iModule = 500; iModule < fgkSSDMODULES + 500; iModule++) {
88       AliITSgeomTGeo::GetModuleId(iModule,gLayer,gLadder,gModule);
89       gTitle = "SSD_RawSignal_Layer"; gTitle += gLayer;
90       gTitle += "_Ladder"; gTitle += gLadder;
91       gTitle += "_Module"; gTitle += gModule;
92       fHistSSDRawSignalModule[gHistCounter] = new TH1D(gTitle.Data(),gTitle.Data(),
93                                                        2*fgkNumberOfPSideStrips,0,2*fgkNumberOfPSideStrips);
94       gHistCounter += 1;
95       
96       for(Int_t iStrip = 0; iStrip < 2*fgkNumberOfPSideStrips; iStrip++)
97         fOccupancyMatrix[iModule-500][iStrip] = 0; 
98     }//module loop
99   }//online flag
100   else {
101     for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) 
102       fHistSSDRawSignalModule[iModule]=NULL;
103     fCDBManager = NULL;
104   }
105 }
106
107 //____________________________________________________________________________ 
108 AliITSQASSDDataMakerRec::AliITSQASSDDataMakerRec(const AliITSQASSDDataMakerRec& qadm) :
109 TObject(),
110 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
111 fSSDEvent(qadm.fSSDEvent),
112 fSSDEventPerCycle(qadm.fSSDEventPerCycle),
113 fkOnline(qadm.fkOnline),
114 fLDC(qadm.fLDC),
115 fSSDRawsOffset(qadm.fSSDRawsOffset), fSSDRawsDAOffset(qadm.fSSDRawsDAOffset),
116 fSSDRawsCommonLevelOffset(qadm.fSSDRawsCommonLevelOffset),
117 fSSDhRawsTask(qadm.fSSDhRawsTask),
118 fSSDhDigitsTask(qadm.fSSDhDigitsTask),
119 fSSDhRecPointsTask(qadm.fSSDhRecPointsTask),
120 fGenRawsOffset(qadm.fGenRawsOffset),
121 fGenDigitsOffset(qadm.fGenDigitsOffset),
122 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
123 fCDBManager(qadm.fCDBManager) {
124   //copy ctor 
125   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; 
126   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
127 }
128
129 //__________________________________________________________________
130 AliITSQASSDDataMakerRec& AliITSQASSDDataMakerRec::operator = (const AliITSQASSDDataMakerRec& qac )
131 {
132   // Equal operator.
133   this->~AliITSQASSDDataMakerRec();
134   new(this) AliITSQASSDDataMakerRec(qac);
135   return *this;
136 }
137
138 //__________________________________________________________________
139 AliITSQASSDDataMakerRec::~AliITSQASSDDataMakerRec() {
140   // destructor
141   for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++)
142     if(fHistSSDRawSignalModule[iModule]) delete fHistSSDRawSignalModule[iModule];
143   if(fCDBManager) delete fCDBManager;
144 }
145
146 //____________________________________________________________________________ 
147 void AliITSQASSDDataMakerRec::StartOfDetectorCycle()
148 {
149   if (  fAliITSQADataMakerRec->GetRawsData(0) == NULL ) // Raws not defined
150     return ; 
151
152   //Detector specific actions at start of cycle
153   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SSD Cycle\n");    
154
155   //Data size per DDL
156   ((TH1D *)(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+4)))->Reset();
157   //Data size per LDC
158   ((TH1D *)(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+22)))->Reset();
159
160   //online part
161   if(fkOnline) {
162     for(Int_t iModule = 500; iModule < fgkSSDMODULES + 500; iModule++) {
163       for(Int_t iStrip = 0; iStrip < 2*fgkNumberOfPSideStrips; iStrip++)
164         fOccupancyMatrix[iModule-500][iStrip] = 0; 
165     }//module loop
166
167     Int_t gHistPositionOccupancyPerLadder = 0;
168     Int_t gLayer = 0, gLadder = 0, gModule = 0;
169     for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
170       AliITSgeomTGeo::GetModuleId(iModule+500,gLayer,gLadder,gModule);
171       
172       gHistPositionOccupancyPerLadder = (gLayer == 5) ? 2*(gLadder - 1) : 2*(gLadder - 1 + fgkSSDLADDERSLAYER5);
173       
174       //P-SIDE OCCUPANCY
175       fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder)->Reset();
176       //N-SIDE OCCUPANCY
177       fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder+1)->Reset();
178
179       ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6))->Reset();
180       ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+1))->Reset();
181     }//module loop
182   }//online flag
183 }
184
185 //____________________________________________________________________________ 
186 void AliITSQASSDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
187 {
188   // launch the QA checking
189   // launch the QA checking
190   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); 
191   AliDebug(AliQAv1::GetQADebugLevel(), Form("Offset: %d\n",fGenRawsOffset));
192
193   if (fAliITSQADataMakerRec->GetRawsData(0) != NULL ) {
194     //Data size per DDL
195     for(Int_t i = 0; i < fgkNumOfDDLs; i++) {
196       Double_t gSizePerDDL = ((TH1D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+5+i))->GetMean();
197       //cout<<"DDL: "<<i+2<<" - Size: "<<gSizePerDDL<<" - Mean: "<<
198       //(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+5+i))->GetMean()<<endl;
199       ((TH1D *)(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+4)))->SetBinContent(i+2,gSizePerDDL);
200     }
201     
202     //Data size per LDC
203     for(Int_t i = 0; i < fgkNumOfLDCs; i++) {
204       Double_t gSizePerLDC = ((TH1D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+23+i))->GetMean();
205       ((TH1D *)(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+22)))->SetBinContent(i+6,gSizePerLDC);
206     }
207   }//raw data end of cycle
208
209   //online part
210   if(fkOnline) {
211     //Output of the DA
212     MonitorOCDBObjects();
213
214     Int_t gHistPositionOccupancyPerModule = 0;
215     Int_t gLayer = 0, gLadder = 0, gModule = 0;
216     //occupancy per module
217     for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
218       AliITSgeomTGeo::GetModuleId(iModule+500,gLayer,gLadder,gModule);
219
220       gHistPositionOccupancyPerModule = (gLayer == 5) ? ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + gModule - 1) : ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + gModule + fgkSSDMODULESLAYER5 - 1);
221       for(Int_t iBins = 1; iBins < fHistSSDRawSignalModule[iModule]->GetXaxis()->GetNbins(); iBins++)
222         fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule)->SetBinContent(iBins,fOccupancyMatrix[iModule][iBins-1]);
223       
224       if(fSSDEventPerCycle != 0)
225         ((TH1F *)(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule)))->Scale(100./fSSDEventPerCycle);
226     }//module loop
227
228     //occupancy per ladder
229     Int_t gHistPositionOccupancyPerLadder = 0;
230     Int_t lLadderLocationY = 0;
231     Double_t occupancy = 0.0, occupancyThreshold = 0.0, occupancyAverage = 0.0;
232     for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
233       AliITSgeomTGeo::GetModuleId(iModule+500,gLayer,gLadder,gModule);
234       
235       gHistPositionOccupancyPerModule = (gLayer == 5) ? ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + gModule - 1) : ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + gModule + fgkSSDMODULESLAYER5 - 1);
236       gHistPositionOccupancyPerLadder = (gLayer == 5) ? 2*(gLadder - 1) : 2*(gLadder - 1 + fgkSSDLADDERSLAYER5);
237       
238       //P-SIDE OCCUPANCY
239       occupancy = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),0,0,0);
240       occupancyThreshold = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),0,1,3);
241       occupancyAverage = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),0,2,0);
242
243       fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder)->Fill(gModule,occupancy);
244       lLadderLocationY = 3*gLadder; // sideP=1 sideN=0 
245       if(gLayer == 5) {
246         //occupancy per module - no threshold
247         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6))->SetBinContent(gModule,lLadderLocationY,occupancy);
248         //occupancy per module - threshold @ 3%
249         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+2))->SetBinContent(gModule,lLadderLocationY,occupancyThreshold);
250         //average occupancy per module
251         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+4))->SetBinContent(gModule,lLadderLocationY,occupancyAverage);
252       }
253       else if(gLayer == 6) {
254         //occupancy per module - no threshold
255         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+1))->SetBinContent(gModule,lLadderLocationY,occupancy);
256         //occupancy per module - threshold @ 3%
257         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+3))->SetBinContent(gModule,lLadderLocationY,occupancyThreshold);
258         //average occupancy per module
259         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+5))->SetBinContent(gModule,lLadderLocationY,occupancyAverage);
260       }
261
262       //N-SIDE OCCUPANCY
263       occupancy = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),1,0,0);   
264       occupancyThreshold = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),1,1,3);   
265       occupancyAverage = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),1,2,0);   
266
267       fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder+1)->Fill(gModule,occupancy);
268       if(gLayer == 5) {
269         //occupancy per module - no threshold
270         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6))->SetBinContent(gModule,lLadderLocationY-1,occupancy);
271         //occupancy per module - threshold @ 3%
272         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+2))->SetBinContent(gModule,lLadderLocationY-1,occupancyThreshold);
273         //average occupancy per module
274         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+4))->SetBinContent(gModule,lLadderLocationY-1,occupancyAverage);
275       }
276       else if(gLayer == 6) {
277         //occupancy per module - no threshold
278         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+1))->SetBinContent(gModule,lLadderLocationY-1,occupancy);
279         //occupancy per module - threshold @ 3%
280         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+3))->SetBinContent(gModule,lLadderLocationY-1,occupancyThreshold);
281         //average occupancy per module
282         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+5))->SetBinContent(gModule,lLadderLocationY-1,occupancyAverage);
283       }
284     }//module loop
285   }//online flag for SSD
286
287   fSSDEventPerCycle = 0;
288   
289  // AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);
290 }
291
292 //____________________________________________________________________________ 
293 Int_t AliITSQASSDDataMakerRec::InitRaws() {  
294   // Initialization for RAW data - SSD -
295   const Bool_t expert   = kTRUE ; 
296   const Bool_t saveCorr = kTRUE ; 
297   const Bool_t image    = kTRUE ; 
298   Int_t rv = 0 ; 
299 //  fGenRawsOffset = (fAliITSQADataMakerRec->fRawsQAList[AliRecoParam::kDefault])->GetEntries();
300
301   if(fkOnline) {
302     AliDebug(AliQAv1::GetQADebugLevel(), "Book Online Histograms for SSD\n");
303   }
304   else {
305     AliDebug(AliQAv1::GetQADebugLevel(), "Book Offline Histograms for SSD\n ");
306   }
307   AliDebug(AliQAv1::GetQADebugLevel(), Form("Number of histograms (SPD+SDD): %d\n",fGenRawsOffset));
308   TString gTitle;
309   TString gName;
310   //book online-offline QA histos
311   TH1D *fHistSSDEventType = new TH1D("fHistSSDEventType",
312                                      "SSD Event Type;Event type;Events",
313                                      31,-1,30);
314   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDEventType, 
315                                       fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
316   fSSDRawsOffset += 1;
317   TH1D *fHistSSDDataSize = new TH1D("fHistSSDDataSize",
318                                     "SSD Data Size;(SSD data size) [KB];Events",
319                                     1000,0,500);
320   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSize, 
321                                       fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
322   fSSDRawsOffset += 1;
323   TH1D *fHistSSDDataSizePercentage = new TH1D("fHistSSDDataSizePercentage",
324                                               "SSD Data Size Percentage;SSD data size [%];Events",
325                                               1000,0,100);
326   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePercentage, 
327                                       fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
328   fSSDRawsOffset += 1;
329   TH1D *fHistSSDDDLId = new TH1D("fHistSSDDDLId",
330                                  "SSD DDL Id;DDL id;Events",20,510.5,530.5);
331   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDDLId, 
332                                       fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
333   fSSDRawsOffset += 1;
334   TH1D *fHistSSDDataSizePerDDL = new TH1D("fHistSSDDataSizePerDDL",
335                                           "SSD Data Size Per DDL;DDL id;<SSD data size> [KB]",
336                                           20,510.5,530.5);
337   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePerDDL, 
338                                       fGenRawsOffset+fSSDRawsOffset, !expert, image, !saveCorr);
339   fSSDRawsOffset += 1;
340   TH1D *fHistSSDDataSizeDDL[fgkNumOfDDLs];
341   for(Int_t i = 1; i < fgkNumOfDDLs+1; i++) {
342     gName = Form("fHistSSDDataSizeDDL%d", i+511) ;
343     gTitle = Form("SSD Data Size DDL %d", i+511) ;
344     fHistSSDDataSizeDDL[i-1] = new TH1D(gName.Data(),
345                                         Form("%s;(SSD data size) [KB];Events", gTitle.Data()),
346                                         1000,0,50);
347     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizeDDL[i-1], 
348                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
349     fSSDRawsOffset += 1;
350   }
351   
352   TH1D *fHistSSDLDCId = new TH1D("fHistSSDLDCId","SSD LDC Id;LDC id;Events",10,0.5,10.5);
353   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDLDCId, 
354                                       fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
355   fSSDRawsOffset += 1;
356   TH1D *fHistSSDDataSizePerLDC = new TH1D("fHistSSDDataSizePerLDC",
357                                           "SSD Data Size Per LDC;LDC id;<SSD data size> [KB]",
358                                           20,0.5,20.5);
359   rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePerLDC, 
360                                       fGenRawsOffset+fSSDRawsOffset, !expert, image, !saveCorr);
361   fSSDRawsOffset += 1;
362   TH1D *fHistSSDDataSizeLDC[fgkNumOfLDCs];
363   for(Int_t i = 1; i < fgkNumOfLDCs+1; i++) {
364     gName = "fHistSSDDataSizeLDC"; 
365     if(i == 1) gName += "082";
366     if(i == 2) gName += "086";
367     if(i == 3) gName += "085";
368     
369     gTitle = "SSD Data Size LDC "; 
370     if(i == 1) gTitle += "082";
371     if(i == 2) gTitle += "086";
372     if(i == 3) gTitle += "085";
373     fHistSSDDataSizeLDC[i-1] = new TH1D(gName.Data(),
374                                         Form("%s;SSD data size [KB];Events", gTitle.Data()),
375                                         1000,0,100);
376     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizeLDC[i-1], 
377                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
378     fSSDRawsOffset += 1;
379   }
380   fSSDRawsCommonLevelOffset = fSSDRawsOffset;
381
382   if(fkOnline) {
383     Int_t gLayer = 0, gLadder = 0, gModule = 0;
384     //occupancy per SSD module
385     TH1D *fHistSSDOccupancyModule[fgkSSDMODULES]; 
386     for(Int_t i = 500; i < fgkSSDMODULES + 500; i++) {
387       AliITSgeomTGeo::GetModuleId(i,gLayer,gLadder,gModule);
388       gName = "fHistSSD_Occupancy_Layer";
389       gTitle = "SSD Occupancy Layer";
390       if(gLayer == 5) {
391         gName += gLayer; gName += "_Ladder"; 
392         gName += 499+gLadder;
393         gTitle += gLayer; gTitle += "_Ladder"; 
394         gTitle += 499+gLadder;
395       }
396       if(gLayer == 6) {
397         gName += gLayer; gName += "_Ladder"; 
398         gName += 599+gLadder;
399         gTitle += gLayer; gTitle += "_Ladder"; 
400         gTitle += 599+gLadder;
401       }
402       gName += "_Module"; gName += gModule; 
403       gTitle += "_Module"; gTitle += gModule; 
404
405       fHistSSDOccupancyModule[i-500] = new TH1D(gName.Data(),Form("%s;N_{strip};Occupancy [%]", gTitle.Data()),
406                                                 2*fgkNumberOfPSideStrips,0,2*fgkNumberOfPSideStrips);
407       fHistSSDOccupancyModule[i-500]->GetXaxis()->SetTitleColor(1);
408       rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyModule[i-500], 
409                                           fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
410       fSSDRawsOffset += 1;
411     }
412
413     //Occupancy per SSD ladder
414     Int_t occupancyCounter = 0;
415     TH1D *fHistSSDOccupancyLadder[2*(fgkSSDLADDERSLAYER5 + fgkSSDLADDERSLAYER6)];
416     for(Int_t iLayer = 5; iLayer < 7; iLayer++) {
417       for(Int_t iLadder = 1; iLadder < AliITSgeomTGeo::GetNLadders(iLayer) + 1; iLadder++) {
418         //P-side occupancy plots
419         gName  = "fHistSSD_Occupancy_Layer"; 
420         gTitle = "SSD Occupancy Layer"; 
421         if(iLayer == 5) {
422           gName += iLayer; gName += "_Ladder"; 
423           gName += 499+iLadder;
424           gTitle += iLayer; gTitle += "_Ladder"; 
425           gTitle += 499+iLadder;
426         }
427         if(iLayer == 6) {
428           gName += iLayer; gName += "_Ladder"; 
429           gName += 599+iLadder;
430           gTitle += iLayer; gTitle += "_Ladder"; 
431           gTitle += 599+iLadder;
432         }
433         gName += "_PSide";
434         gTitle += "_PSide";
435         fHistSSDOccupancyLadder[occupancyCounter] = new TH1D(gName.Data(),
436                                                              Form("%s;Module number;Occupancy [%]", gTitle.Data()),
437                                                              AliITSgeomTGeo::GetNDetectors(iLayer),
438                                                              0.5,AliITSgeomTGeo::GetNDetectors(iLayer)+0.5);
439         fHistSSDOccupancyLadder[occupancyCounter]->GetXaxis()->SetTitleColor(1);
440         rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLadder[occupancyCounter], 
441                                             fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
442         occupancyCounter += 1; fSSDRawsOffset += 1;
443         //N-side occupancy plots
444         gName = "fHistSSD_Occupancy_Layer"; 
445         gTitle = "SSD Occupancy Layer"; 
446         if(iLayer == 5) {
447           gName += iLayer; gName += "_Ladder"; 
448           gName += 499+iLadder;
449           gTitle += iLayer; gTitle += "_Ladder"; 
450           gTitle += 499+iLadder;
451         }
452         if(iLayer == 6) {
453           gName += iLayer; gName += "_Ladder"; 
454           gName += 599+iLadder;
455           gTitle += iLayer; gTitle += "_Ladder"; 
456           gTitle += 599+iLadder;
457         }
458         gName += "_NSide";
459         gTitle += "_NSide";
460         fHistSSDOccupancyLadder[occupancyCounter] = new TH1D(gName.Data(),
461                                                              Form("%s;Module number;Occupancy [%]", gTitle.Data()),
462                                                              AliITSgeomTGeo::GetNDetectors(iLayer),
463                                                              0.5,AliITSgeomTGeo::GetNDetectors(iLayer)+0.5);
464         fHistSSDOccupancyLadder[occupancyCounter]->GetXaxis()->SetTitleColor(1);
465         rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLadder[occupancyCounter], 
466                                             fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
467         occupancyCounter += 1; fSSDRawsOffset += 1;
468       }//ladder loop
469     }//layer loop
470
471     //top level occupancy plots
472     //occupancy per module - no threshold
473     TH2D *fHistSSDOccupancyLayer5 = new TH2D("fHistSSDOccupancyLayer5",
474                                              "SSD Occupancy (Layer 5) - No threshold;N_{modules};N_{Ladders}",
475                                              fgkSSDMODULESPERLADDERLAYER5,
476                                              0,fgkSSDMODULESPERLADDERLAYER5,
477                                              3*fgkSSDLADDERSLAYER5,
478                                              5000,500+fgkSSDLADDERSLAYER5);  
479     fHistSSDOccupancyLayer5->GetZaxis()->SetRangeUser(0.0,100.0);
480     Char_t fLabel[3];
481     for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER5 + 1; iBin++){
482       sprintf(fLabel,"%d",iBin);
483       fHistSSDOccupancyLayer5->GetXaxis()->SetBinLabel(iBin,fLabel);
484     }
485     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLayer5, 
486                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
487     fSSDRawsOffset += 1;
488     TH2D *fHistSSDOccupancyLayer6 = new TH2D("fHistSSDOccupancyLayer6",
489                                              "Occupancy per module (Layer 6) - No threshold;N_{modules};N_{Ladders}",
490                                              fgkSSDMODULESPERLADDERLAYER6,
491                                              0,fgkSSDMODULESPERLADDERLAYER6,
492                                              3*fgkSSDLADDERSLAYER6,
493                                              600,600+fgkSSDLADDERSLAYER6);
494     fHistSSDOccupancyLayer6->GetZaxis()->SetRangeUser(0.0,100.0);
495     for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER6 + 1; iBin++){
496       sprintf(fLabel,"%d",iBin);
497       fHistSSDOccupancyLayer6->GetXaxis()->SetBinLabel(iBin,fLabel);
498     }
499     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLayer6, 
500                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
501     fSSDRawsOffset += 1;
502
503     //occupancy per module - threshold @ 3%
504     TH2D *fHistSSDOccupancyThresholdLayer5 = new TH2D("fHistSSDOccupancyThresholdLayer5",
505                                                       "Occupancy per module (Layer 5) - Threshold 3%;N_{modules};N_{Ladders};Entries",
506                                                       fgkSSDMODULESPERLADDERLAYER5,
507                                                       0,fgkSSDMODULESPERLADDERLAYER5,
508                                                       3*fgkSSDLADDERSLAYER5,
509                                                       500,500+fgkSSDLADDERSLAYER5);  
510     fHistSSDOccupancyThresholdLayer5->GetZaxis()->SetRangeUser(3.0,10.0);
511     for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER5 + 1; iBin++){
512       sprintf(fLabel,"%d",iBin);
513       fHistSSDOccupancyThresholdLayer5->GetXaxis()->SetBinLabel(iBin,fLabel);
514     }
515     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyThresholdLayer5, 
516                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
517     fSSDRawsOffset += 1;
518     TH2D *fHistSSDOccupancyThresholdLayer6 = new TH2D("fHistSSDOccupancyThresholdLayer6",
519                                                       "Occupancy per module (Layer 6) - Threshold 3%;N_{modules};N_{Ladders}",
520                                                       fgkSSDMODULESPERLADDERLAYER6,
521                                                       0,fgkSSDMODULESPERLADDERLAYER6,
522                                                       3*fgkSSDLADDERSLAYER6,
523                                                       600,600+fgkSSDLADDERSLAYER6);
524     fHistSSDOccupancyThresholdLayer6->GetZaxis()->SetRangeUser(3.0,10.0);
525     for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER6 + 1; iBin++){
526       sprintf(fLabel,"%d",iBin);
527       fHistSSDOccupancyThresholdLayer6->GetXaxis()->SetBinLabel(iBin,fLabel);
528     }
529     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyThresholdLayer6, 
530                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
531     fSSDRawsOffset += 1;
532
533     //Average occupancy per module
534     TH2D *fHistSSDAverageOccupancyLayer5 = new TH2D("fHistSSDAverageOccupancyLayer5",
535                                                     "Average occupancy per module (Layer 5);N_{modules};N_{Ladders}",
536                                                     fgkSSDMODULESPERLADDERLAYER5,
537                                                     0,fgkSSDMODULESPERLADDERLAYER5,
538                                                     3*fgkSSDLADDERSLAYER5,
539                                                     500,500+fgkSSDLADDERSLAYER5);  
540     fHistSSDAverageOccupancyLayer5->GetZaxis()->SetRangeUser(0.0,5.0);
541     for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER5 + 1; iBin++){
542       sprintf(fLabel,"%d",iBin);
543       fHistSSDAverageOccupancyLayer5->GetXaxis()->SetBinLabel(iBin,fLabel);
544     }
545     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDAverageOccupancyLayer5, 
546                                         fGenRawsOffset+fSSDRawsOffset);
547     fSSDRawsOffset += 1;
548     TH2D *fHistSSDAverageOccupancyLayer6 = new TH2D("fHistSSDAverageOccupancyLayer6",
549                                                     "Average occupancy per module (Layer 6);N_{modules};N_{Ladders}",
550                                                     fgkSSDMODULESPERLADDERLAYER6,
551                                                     0,fgkSSDMODULESPERLADDERLAYER6,
552                                                     3*fgkSSDLADDERSLAYER6,
553                                                     600,600+fgkSSDLADDERSLAYER6);
554     fHistSSDAverageOccupancyLayer6->GetZaxis()->SetRangeUser(0.0,5.0);
555     for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER6 + 1; iBin++){
556       sprintf(fLabel,"%d",iBin);
557       fHistSSDAverageOccupancyLayer6->GetXaxis()->SetBinLabel(iBin,fLabel);
558     }
559     rv = fAliITSQADataMakerRec->Add2RawsList(fHistSSDAverageOccupancyLayer6, 
560                                         fGenRawsOffset+fSSDRawsOffset, !expert, image, !saveCorr);
561     fSSDRawsOffset += 1;
562
563     //Output of the DA
564     TH2D *fHistPSideBadChannelMapLayer5 = new TH2D("fHistPSideBadChannelMapLayer5",
565                                                    "Layer 5;N_{module};N_{ladder}",
566                                                    22,1,23,
567                                                    34,500,534);
568     fHistPSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
569     fHistPSideBadChannelMapLayer5->SetStats(kFALSE);
570     fHistPSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
571     fHistPSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
572     fHistPSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
573     fHistPSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
574     fHistPSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
575     fHistPSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
576     fHistPSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
577     rv = fAliITSQADataMakerRec->Add2RawsList(fHistPSideBadChannelMapLayer5, 
578                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
579     fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
580
581     TH2D *fHistNSideBadChannelMapLayer5 = new TH2D("fHistNSideBadChannelMapLayer5",
582                                                    "Layer 5;N_{module};N_{ladder}",
583                                                    22,1,23,
584                                                    34,500,534);
585     fHistNSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
586     fHistNSideBadChannelMapLayer5->SetStats(kFALSE);
587     fHistNSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
588     fHistNSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
589     fHistNSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
590     fHistNSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
591     fHistNSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
592     fHistNSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
593     fHistNSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
594     rv = fAliITSQADataMakerRec->Add2RawsList(fHistNSideBadChannelMapLayer5, 
595                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
596     fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
597     
598     TH2D *fHistPSideBadChannelMapLayer6 = new TH2D("fHistPSideBadChannelMapLayer6",
599                                                    "Layer 6;N_{module};N_{ladder}",
600                                                    25,1,26,
601                                                    38,600,638);
602     fHistPSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
603     fHistPSideBadChannelMapLayer6->SetStats(kFALSE);
604     fHistPSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
605     fHistPSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
606     fHistPSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
607     fHistPSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
608     fHistPSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
609     fHistPSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
610     fHistPSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
611     rv = fAliITSQADataMakerRec->Add2RawsList(fHistPSideBadChannelMapLayer6, 
612                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
613     fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
614     
615     TH2D *fHistNSideBadChannelMapLayer6 = new TH2D("fHistNSideBadChannelMapLayer6",
616                                                    "Layer 6;N_{module};N_{ladder}",
617                                                    25,1,26,
618                                                    38,600,638);
619     fHistNSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
620     fHistNSideBadChannelMapLayer6->SetStats(kFALSE);
621     fHistNSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
622     fHistNSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
623     fHistNSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
624     fHistNSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
625     fHistNSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
626     fHistNSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
627     fHistNSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
628     rv = fAliITSQADataMakerRec->Add2RawsList(fHistNSideBadChannelMapLayer6, 
629                                         fGenRawsOffset+fSSDRawsOffset, expert, !image, !saveCorr);
630     fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
631   }//online flag
632   
633   fSSDhRawsTask = fSSDRawsOffset;
634   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD Raws histograms booked\n",fSSDhRawsTask));
635   AliDebug(AliQAv1::GetQADebugLevel(), Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenRawsOffset+fSSDhRawsTask));  
636   AliDebug(AliQAv1::GetQADebugLevel(),Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenRawsOffset+fSSDRawsOffset));
637   
638   /*
639    fSSDhTask = fSSDRawsOffset;
640    AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD Raws histograms booked\n",fSSDhTask));
641    AliDebug(AliQAv1::GetQADebugLevel(), Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenRawsOffset+fSSDhTask));  
642    AliDebug(AliQAv1::GetQADebugLevel(),Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenRawsOffset+fSSDRawsOffset));
643    */
644   return rv ; 
645 }
646
647 //____________________________________________________________________________
648 Int_t AliITSQASSDDataMakerRec::MakeRaws(AliRawReader* rawReader) { 
649   // Fill QA for RAW - SSD -
650   Int_t rv = 0 ; 
651   // Check id histograms already created for this Event Specie
652   if ( ! fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset) )
653     rv = InitRaws() ;
654
655   Int_t gStripNumber;
656   Int_t gHistPosition;
657   Int_t gLayer = 0,gLadder = 0, gModule = 0;
658   
659   Double_t gSizePerDDL[fgkNumOfDDLs] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
660   Double_t gSizePerLDC[fgkNumOfLDCs] = {0.,0.,0.};
661   Double_t sumSSDDataSize = 0.0;
662   Double_t eventSize = -1.0;
663
664   if(fkOnline) {
665     //reset the signal vs strip number histograms
666     for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++)
667       fHistSSDRawSignalModule[iModule]->Reset();
668   }//online flag
669
670   rawReader->Select("ITSSSD",-1,-1);  
671   rawReader->Reset(); //rawReader->NextEvent();   
672   (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset))->Fill(rawReader->GetType());
673   if(rawReader->GetType() == 7) {
674     fSSDEvent += 1;
675     fSSDEventPerCycle += 1;
676   }
677
678   AliITSRawStreamSSD gSSDStream(rawReader);    
679   AliRawReaderRoot *rootReader = (AliRawReaderRoot *)rawReader;
680   if(rootReader) {
681     const AliRawEventHeaderBase *header = rootReader->GetEventHeader();
682     if(header)
683       eventSize = header->GetEventSize();
684   }
685   while (gSSDStream.Next()) {
686     if(gSSDStream.GetModuleID() < 0) continue;
687     /*cout<<"DDL: "<<rawReader->GetDDLID()<<
688       " - LDC: "<<rawReader->GetLDCId()<<
689       " - Size: "<<rawReader->GetDataSize()<<
690       " - Equipment size: "<<rawReader->GetEquipmentSize()<<endl;*/
691     gSizePerDDL[rawReader->GetDDLID()] = rawReader->GetDataSize();
692     gSizePerLDC[rawReader->GetLDCId()-6] = rawReader->GetDataSize();
693     AliITSgeomTGeo::GetModuleId(gSSDStream.GetModuleID(),gLayer,gLadder,gModule);
694     if(gSSDStream.GetStrip() < 0) continue;
695     gStripNumber = (gSSDStream.GetSideFlag() == 0) ? gSSDStream.GetStrip() : -gSSDStream.GetStrip() + 2*fgkNumberOfPSideStrips;
696     gHistPosition = (gLayer == 5) ? ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + gModule - 1) : ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + gModule + fgkSSDMODULESLAYER5 - 1);
697     //AliDebug(AliQAv1::GetQADebugLevel(), Form("ModulePosition: %d - Layer: %d - Ladder: %d - Module: %d\n",gHistPosition,gLayer,gLadder,gModule));
698     if(fkOnline)
699       fHistSSDRawSignalModule[gHistPosition]->Fill(gStripNumber,gSSDStream.GetSignal());
700     //fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+gHistPosition+fSSDRawsCommonLevelOffset)->Fill(gStripNumber,gSSDStream.GetSignal());
701   }//streamer loop   
702   
703   //event size calculation and filling info
704   for(Int_t i = 0; i < fgkNumOfDDLs; i++) {
705     sumSSDDataSize += gSizePerDDL[i];
706     if(gSizePerDDL[i] > 0) {
707       (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+3))->Fill(i+512);
708       (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+5+i))->Fill(gSizePerDDL[i]/1e+03);
709     }
710     //(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+4))->Fill(i+512,gSizePerDDL[i]/1e+06);
711   }
712   for(Int_t i = 0; i < fgkNumOfLDCs; i++) {
713     if(gSizePerLDC[i] > 0) {
714       (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+21))->Fill(i+6);
715       //LDC 082
716       if(i == 0)
717         gSizePerLDC[i] = gSizePerDDL[8] + gSizePerDDL[9] + gSizePerDDL[10] +
718           gSizePerDDL[11] + gSizePerDDL[12] + gSizePerDDL[13];
719       //LDC 086
720       if(i == 1)
721         gSizePerLDC[i] = gSizePerDDL[3] + gSizePerDDL[4] + gSizePerDDL[5] +
722           gSizePerDDL[6] + gSizePerDDL[7];
723       //LDC 085
724       if(i == 2)
725         gSizePerLDC[i] = gSizePerDDL[0] + gSizePerDDL[1] + gSizePerDDL[2] +
726           gSizePerDDL[14] + gSizePerDDL[15];
727       (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+23+i))->Fill(gSizePerLDC[i]/1e+03);
728       //cout<<"Event: "<<fSSDEventPerCycle<<" - LDC: "<<i+6<<
729       //" - Data size: "<<gSizePerLDC[i]<<endl;
730     }
731     //(fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+22))->Fill(i+6,gSizePerLDC[i]/1e+06);
732   }
733   if(sumSSDDataSize) 
734     (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+1))->Fill(sumSSDDataSize/1e+03);
735   if(eventSize)
736     (fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+2))->Fill(100.*sumSSDDataSize/eventSize);
737
738   //Occupancy calculation
739   if(fkOnline) {
740     for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
741       GetOccupancyStrip(fHistSSDRawSignalModule[iModule],fOccupancyMatrix[iModule]);
742       //if(iModule == 156) {    
743         //cout<<"========================================"<<endl;
744         //AliITSgeomTGeo::GetModuleId(656,gLayer,gLadder,gModule);
745         /*for(Int_t iBin = 1; iBin < fHistSSDRawSignalModule[iModule]->GetXaxis()->GetNbins(); iBin++) {
746           if((iBin >= 750)&&(iBin <= 780))
747           cout<<"Event: "<<fSSDEventPerCycle<<" - Ladder: "<<gLadder+499<<
748           " - Module: "<<gModule<<" - Strip: "<<iBin<<
749           " - Signal: "<<fHistSSDRawSignalModule[iModule]->GetBinContent(iBin)<<endl;
750           }*///strip loop --> to be removed
751       //}//module cut --> to be removed
752     }//module loop
753   }//online flag for SSD
754   return rv ; 
755 }
756
757 //____________________________________________________________________________ //
758 void AliITSQASSDDataMakerRec::GetOccupancyStrip(TH1 *lHisto, Int_t *occupancyMatrix) { 
759   //Increments the entries in the occupancy matrix based 
760   //on whether the signal for each strip is larger than the cutoff
761   //Currently the cutoff is at 0 which means that if ZS
762   //works, whatever comes from the FEROM is considered as "signal"
763   //Double_t cutoff = 0.0;
764   TString histname = lHisto->GetName();
765   //cout<<histname.Data()<<endl;
766   //if(histname.Contains("Layer5_Ladder8_Module3")) {
767   for(Int_t iBin = 1; iBin < lHisto->GetXaxis()->GetNbins(); iBin++) {
768     Double_t y = lHisto->GetBinContent(iBin);
769     if(y) {
770       occupancyMatrix[iBin-1] += 1;
771     }
772     //if((iBin >= 750)&&(iBin <= 780))
773     //cout<<"Strip: "<<iBin<<
774     //" - Signal: "<<y<<
775     //" - Occupancy: "<<occupancyMatrix[iBin-1]<<endl;
776     
777   }
778   //}
779 }
780
781 //____________________________________________________________________________ 
782 Double_t AliITSQASSDDataMakerRec::GetOccupancyModule(TH1 *lHisto, 
783                                                      Int_t stripside,
784                                                      Int_t mode,
785                                                      Double_t threshold) {
786   //Mode 0: calculates the occupancy of a module 
787   //        without a threshold in the strip occupancy
788   //Mode 1: calculates the occupancy of a module 
789   //        with the set threshold in the strip occupancy
790   //Mode 2: calculates the average occupancy of a module 
791
792   // bo: TDC >0 or # of sigmas wrt noise ?
793   //stripside == 0 --> P-side
794   //stripside == 1 --> N-side
795   Int_t lNumFiredBins = 0;
796   Double_t sumOccupancy = 0.0;
797   TString histname = lHisto->GetName();
798   for(Int_t iBin = 1 + stripside*fgkNumberOfPSideStrips; iBin < fgkNumberOfPSideStrips*(1 + stripside); iBin++){
799     /*if(histname.Contains("Layer5_Ladder507_Module3")) {
800       cout<<lHisto->GetName()<<
801       " - Strip: "<<iBin<<
802       " - Bin content: "<<lHisto->GetBinContent(iBin)<<endl;
803       }*/
804     if (lHisto->GetBinContent(iBin) > threshold) {
805       lNumFiredBins++; 
806       sumOccupancy = lHisto->GetBinContent(iBin);
807     }
808   }
809   
810   Double_t lOccupancy = 0.0;
811   if((mode == 0)||(mode == 1))
812     lOccupancy = (100.*lNumFiredBins)/fgkNumberOfPSideStrips; // percentage  
813   else if(mode == 2)
814     lOccupancy = 100.*sumOccupancy/fgkNumberOfPSideStrips;
815
816   /*if(histname.Contains("Layer5_Ladder507_Module3"))
817     cout<<"Fired strips: "<<lNumFiredBins<<
818     " - Occupancy: "<<lOccupancy<<endl;*/
819   //AliDebug(AliQAv1::GetQADebugLevel(), Form("Fired strips: %d - Total strips: %d - Occupancy :%lf\n",lNumFiredBins,lHisto->GetNbinsX(),lOccupancy));
820   
821   return lOccupancy;
822 }
823
824 //____________________________________________________________________________
825 void AliITSQASSDDataMakerRec::MonitorOCDBObjects() { 
826   //Monitor in AMORE the output of the DA
827   //Currently only the bad channel list is monitored
828   //Todo: Noise - Pedestal
829   ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset))->Reset();
830   ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+1))->Reset();
831   ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+2))->Reset();
832   ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+3))->Reset();
833
834   AliCDBEntry *entryBadChannelsSSD = fCDBManager->Get("ITS/Calib/BadChannelsSSD");
835   if(!entryBadChannelsSSD) 
836     AliError("OCDB entry for the bad channel list is not valid!"); 
837   AliITSBadChannelsSSDv2 *badchannelsSSD = (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
838   if(!badchannelsSSD)
839     AliError("Bad channel list object is not a valid AliITSBadChannelsSSD object!");
840
841   //_____________________________________________________________________________//                       
842   Int_t nBadPSideChannels = 0, nBadNSideChannels = 0;
843   Int_t layer = 0, ladder = 0, module = 0;
844   Int_t nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
845   Int_t nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
846   //_____________________________________________________________________________//                      
847
848   for(Int_t i = 0; i < fgkSSDMODULES; i++) {
849     AliITSgeomTGeo::GetModuleId(i+500,layer,ladder,module);
850     nBadPSideChannels = 0, nBadNSideChannels = 0;
851     nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
852     nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
853
854     Int_t badChannel = 0;
855     for(Int_t j = 0; j < fgkNumberOfPSideStrips; j++) {
856       badChannel = (Int_t)(badchannelsSSD->GetBadChannelP(i,j));
857       if(badChannel != 0) {
858         if(layer == 5)
859           nPSideChannelsLayer5 += 1;
860         if(layer == 6)
861           nPSideChannelsLayer6 += 1;
862         nBadPSideChannels += 1;
863       }//badchannel flag != 0
864       badChannel = (Int_t)(badchannelsSSD->GetBadChannelN(i,j));
865       if(badChannel != 0) {
866         if(layer == 5)
867           nNSideChannelsLayer5 += 1;
868         if(layer == 6)
869           nNSideChannelsLayer6 += 1;
870         nBadNSideChannels += 1;
871       }//badchannel flag != 0
872     }//loop over strips
873     if(layer == 5) {
874       if(nPSideChannelsLayer5 > 0)
875         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset))->Fill(module,499+ladder,
876                                                                             100.*nPSideChannelsLayer5/fgkNumberOfPSideStrips);
877       else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset))->Fill(module,499+ladder,0.0001);
878       if(nNSideChannelsLayer5 > 0)
879         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+1))->Fill(module,499+ladder,
880                                                                             100.*nNSideChannelsLayer5/fgkNumberOfPSideStrips);
881       else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+1))->Fill(module,499+ladder,0.0001);
882     }//layer 5                                                                                                                      
883     if(layer == 6) {
884       if(nPSideChannelsLayer6 > 0)
885         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+2))->Fill(module,599+ladder,
886                                                                             100.*nPSideChannelsLayer6/fgkNumberOfPSideStrips);
887       else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+2))->Fill(module,599+ladder,0.0001);
888       if(nNSideChannelsLayer6 > 0)
889         ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+3))->Fill(module,599+ladder,
890                                                                             100.*nNSideChannelsLayer6/fgkNumberOfPSideStrips);
891       else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset+fSSDRawsOffset-fSSDRawsDAOffset+3))->Fill(module,599+ladder,0.0001);
892     }//layer 6                                                                                                                      
893   }//module loop
894 }
895
896 //____________________________________________________________________________ 
897 Int_t AliITSQASSDDataMakerRec::InitDigits() { 
898   // Initialization for DIGIT data - SSD -
899   const Bool_t expert   = kTRUE ; 
900   const Bool_t image    = kTRUE ; 
901   Int_t rv = 0 ; 
902 //  fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
903   
904   // custom code here
905   TH1F *fHistSSDModule = new TH1F("fHistSSDDigitsModule",
906                                   ";SSD Module Number;N_{DIGITS}",
907                                   1698,499.5,2197.5);  
908  rv =  fAliITSQADataMakerRec->Add2DigitsList(fHistSSDModule,
909                                         fGenDigitsOffset + 0, !expert, image);
910   fSSDhDigitsTask += 1;
911   TH2F *fHistSSDModuleStrip = new TH2F("fHistSSDDigitsModuleStrip",
912                                        ";N_{Strip};N_{Module}",
913                                        1540,0,1540,1698,499.5,2197.5);  
914   rv = fAliITSQADataMakerRec->Add2DigitsList(fHistSSDModuleStrip,
915                                         fGenDigitsOffset + 1, !expert, image);
916   fSSDhDigitsTask += 1;
917   
918   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD Digits histograms booked\n",fSSDhDigitsTask));
919   return rv ; 
920 }
921
922 //____________________________________________________________________________
923 Int_t AliITSQASSDDataMakerRec::MakeDigits(TTree *digits) { 
924   // Fill QA for DIGIT - SSD -
925 //  AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");
926 //  fITS->SetTreeAddress();
927 //  TClonesArray *iSSDdigits  = fITS->DigitsAddress(2);
928   Int_t rv = 0 ; 
929   TBranch *branchD = digits->GetBranch("ITSDigitsSSD");
930   if (!branchD) { 
931     AliError("can't get the branch with the SSD ITS digits !");
932     return rv;
933   }
934   
935   // Check id histograms already created for this Event Specie
936   if ( ! fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset) )
937     rv = InitDigits() ;
938
939   static TClonesArray statDigits("AliITSdigitSSD");
940   TClonesArray *iSSDdigits = &statDigits;
941   branchD->SetAddress(&iSSDdigits);  
942   for(Int_t iModule = 500; iModule < 2198; iModule++) {
943     iSSDdigits->Clear();
944     digits->GetEvent(iModule);    
945     Int_t ndigits = iSSDdigits->GetEntries();
946     fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset + 0)->Fill(iModule,ndigits);
947     if(ndigits != 0)
948       AliDebug(AliQAv1::GetQADebugLevel(),Form("Module: %d - Digits: %d",iModule,ndigits));
949     
950     for (Int_t iDigit = 0; iDigit < ndigits; iDigit++) {
951       AliITSdigit *dig = (AliITSdigit*)iSSDdigits->UncheckedAt(iDigit);
952       Int_t fStripNumber = (dig->GetCoord1() == 0) ? dig->GetCoord2() : dig->GetCoord2() + fgkNumberOfPSideStrips;
953       ((TH2F *)fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset + 1))->Fill(fStripNumber,iModule,dig->GetSignal());
954     }//digit loop
955   }//module loop
956   return rv ; 
957 }
958
959 //____________________________________________________________________________ 
960 Int_t AliITSQASSDDataMakerRec::InitRecPoints()
961 {
962   // Initialization for RECPOINTS - SSD -
963   const Bool_t expert   = kTRUE ; 
964   const Bool_t image    = kTRUE ; 
965   Int_t rv = 0 ; 
966 //  fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
967   //AliDebug(AliQAv1::GetQADebugLevel(), Form("**-------*-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints offset %d \t %d \n",fGenOffset,fGenRecPointsOffset));
968   Int_t nModuleOffset = 500;
969   Int_t nITSTotalModules = AliITSgeomTGeo::GetNModules();
970
971   TH1F *fHistSSDModuleIdLayer5 = new TH1F("fHistSSDModuleIdLayer5",
972                                           "Module Id - Layer 5;Module Id;Entries",
973                                           fgkSSDMODULESLAYER5,
974                                           nModuleOffset - 0.5,
975                                           nITSTotalModules-fgkSSDMODULESLAYER6+0.5);
976   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDModuleIdLayer5, 
977                                            fGenRecPointsOffset + 0, expert, !image);
978   fSSDhRecPointsTask += 1;
979   TH1F *fHistSSDModuleIdLayer6 = new TH1F("fHistSSDModuleIdLayer6",
980                                           "Module Id - Layer 6;Module Id;Entries",
981                                           fgkSSDMODULESLAYER6,
982                                           nModuleOffset+fgkSSDMODULESLAYER5 - 0.5,
983                                           nITSTotalModules + 0.5);
984   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDModuleIdLayer6, 
985                                            fGenRecPointsOffset + 1, expert, !image);
986   fSSDhRecPointsTask += 1;
987   TH1F *fHistSSDClusterPerEventLayer5 = new TH1F("fHistSSDClusterPerEventLayer5",
988                                                  "N_{clusters} - Layer 5;N_{clusters};Entries;",
989                                                  100,0.1,5000);
990   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDClusterPerEventLayer5,
991                                            fGenRecPointsOffset + 2, expert, !image);
992   fSSDhRecPointsTask += 1;
993   TH1F *fHistSSDClusterPerEventLayer6 = new TH1F("fHistSSDClusterPerEventLayer6",
994                                                  "N_{clusters} - Layer 6;N_{clusters};Entries;",
995                                                  100,0.1,5000);
996   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDClusterPerEventLayer6,
997                                            fGenRecPointsOffset + 3, expert, !image);
998   fSSDhRecPointsTask += 1;
999   TH1F *fHistSSDLocalXLayer5 = new TH1F("fHistSSDLocalXLayer5",
1000                                         "Local x coord.- Layer 5;x [cm];Entries;",
1001                                         100,-4.,4.);
1002   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDLocalXLayer5,
1003                                            fGenRecPointsOffset + 4, expert, !image);
1004   fSSDhRecPointsTask += 1;
1005   TH1F *fHistSSDLocalXLayer6 = new TH1F("fHistSSDLocalXLayer6",
1006                                         "Local x coord.- Layer 6;x [cm];Entries;",
1007                                         100,-4.,4.);
1008   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDLocalXLayer6, 
1009                                            fGenRecPointsOffset + 5, expert, !image);
1010   fSSDhRecPointsTask += 1;
1011   TH1F *fHistSSDLocalZLayer5 = new TH1F("fHistSSDLocalZLayer5",
1012                                         "Local z coord.- Layer 5;z [cm];Entries;",
1013                                         100,-4.,4.);
1014   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDLocalZLayer5, 
1015                                            fGenRecPointsOffset + 6, expert, !image);
1016   fSSDhRecPointsTask += 1;
1017   TH1F *fHistSSDLocalZLayer6 = new TH1F("fHistSSDLocalZLayer6",
1018                                         "Local z coord.- Layer 6;z [cm];Entries;",
1019                                         100,-4.,4.);
1020   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDLocalZLayer6, 
1021                                            fGenRecPointsOffset + 7, expert, !image);
1022   fSSDhRecPointsTask += 1;
1023   TH1F *fHistSSDGlobalXLayer5 = new TH1F("fHistSSDGlobalXLayer5",
1024                                          "Global x - Layer 5;x [cm];Entries;",
1025                                          100,-40.,40.);
1026   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDGlobalXLayer5, 
1027                                            fGenRecPointsOffset + 8, expert, !image);
1028   fSSDhRecPointsTask += 1;
1029   TH1F *fHistSSDGlobalXLayer6 = new TH1F("fHistSSDGlobalXLayer6",
1030                                          "Global x - Layer 6;x [cm];Entries;",
1031                                          100,-45.,45.);
1032   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDGlobalXLayer6, 
1033                                            fGenRecPointsOffset + 9, expert, !image);
1034   fSSDhRecPointsTask += 1;
1035   TH1F *fHistSSDGlobalYLayer5 = new TH1F("fHistSSDGlobalYLayer5",
1036                                          "Global y - Layer 5;y [cm];Entries;",
1037                                          100,-40.,40);
1038   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDGlobalYLayer5, 
1039                                            fGenRecPointsOffset + 10, expert, !image);
1040   fSSDhRecPointsTask += 1;
1041   TH1F *fHistSSDGlobalYLayer6 = new TH1F("fHistSSDGlobalYLayer6",
1042                                          "Global y - Layer 6;y [cm];Entries;",
1043                                          100,-45.,45.);
1044   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDGlobalYLayer6, 
1045                                            fGenRecPointsOffset + 11, expert, !image);
1046   fSSDhRecPointsTask += 1;
1047   TH1F *fHistSSDGlobalZLayer5 = new TH1F("fHistSSDGlobalZLayer5",
1048                                          "Global z - Layer 5;z [cm];Entries;",
1049                                          100,-45.,45);
1050   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDGlobalZLayer5, 
1051                                            fGenRecPointsOffset + 12, expert, !image);
1052   fSSDhRecPointsTask += 1;
1053   TH1F *fHistSSDGlobalZLayer6 = new TH1F("fHistSSDGlobalZLayer6",
1054                                          "Global z - Layer 6;z [cm];Entries;",
1055                                          100,-55.,55.);
1056   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDGlobalZLayer6, 
1057                                            fGenRecPointsOffset + 13, expert, !image);
1058   fSSDhRecPointsTask += 1;
1059   TH1F *fHistSSDPhiLayer5 = new TH1F("fHistSSDPhiLayer5",
1060                                      "#phi - Layer 5;#phi [rad];Entries;",
1061                                      100,-TMath::Pi(),TMath::Pi());
1062   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDPhiLayer5, 
1063                                            fGenRecPointsOffset + 14, expert, !image);
1064   fSSDhRecPointsTask += 1;
1065   TH1F *fHistSSDPhiLayer6 = new TH1F("fHistSSDPhiLayer6",
1066                                      "#phi - Layer 6;#phi [rad];Entries;",
1067                                      100,-TMath::Pi(),TMath::Pi());
1068   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDPhiLayer6, 
1069                                            fGenRecPointsOffset + 15, expert, !image);
1070   fSSDhRecPointsTask += 1;
1071   TH1F *fHistSSDThetaLayer5 = new TH1F("fHistSSDThetaLayer5",
1072                                        "#theta - Layer 5;#theta [rad];Entries;",
1073                                        100,-TMath::Pi(),TMath::Pi());
1074   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDThetaLayer5, 
1075                                            fGenRecPointsOffset + 16, expert, !image);
1076   fSSDhRecPointsTask += 1;
1077   TH1F *fHistSSDThetaLayer6 = new TH1F("fHistSSDThetaLayer6",
1078                                        "#theta - Layer 6;#theta [rad];Entries;",
1079                                        100,-TMath::Pi(),TMath::Pi());
1080   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDThetaLayer6, 
1081                                            fGenRecPointsOffset + 17, expert, !image);
1082   fSSDhRecPointsTask += 1;
1083   TH1F *fHistSSDRadiusLayer5 = new TH1F("fHistSSDRadiusLayer5",
1084                                         "r - Layer 5;r [cm];Entries;",
1085                                         100,35.,50.);
1086   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDRadiusLayer5, 
1087                                            fGenRecPointsOffset + 18, expert, !image);
1088   fSSDhRecPointsTask += 1;
1089   TH1F *fHistSSDRadiusLayer6 = new TH1F("fHistSSDRadiusLayer6",
1090                                         "r - Layer 6;r [cm];Entries;",
1091                                         100,35.,50.);
1092   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDRadiusLayer6, 
1093                                            fGenRecPointsOffset + 19, expert, !image);
1094   fSSDhRecPointsTask += 1;
1095   TH1F *fHistSSDClusterTypeLayer5 = new TH1F("fHistSSDClusterTypeLayer5",
1096                                              "CL type - Layer 5;Cluster type;Entries;",
1097                                              150,0,150);
1098   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDClusterTypeLayer5, 
1099                                            fGenRecPointsOffset + 20, expert, !image);
1100   fSSDhRecPointsTask += 1;
1101   TH1F *fHistSSDClusterTypeLayer6 = new TH1F("fHistSSDClusterTypeLayer6",
1102                                              "CL type - Layer 6;Cluster type;Entries;",
1103                                              150,0,150);
1104   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDClusterTypeLayer6, 
1105                                            fGenRecPointsOffset + 21, expert, !image);
1106   fSSDhRecPointsTask += 1;
1107   TH1F *fHistSSDChargeRatioLayer5 = new TH1F("fHistSSDChargeRatioLayer5",
1108                                              "Charge ratio - Layer 5;q_{ratio};Entries;",
1109                                              100,-2.0,2.0);
1110   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeRatioLayer5, 
1111                                            fGenRecPointsOffset + 22, expert, !image);
1112   fSSDhRecPointsTask += 1;
1113   TH1F *fHistSSDChargeRatioLayer6 = new TH1F("fHistSSDChargeRatioLayer6",
1114                                              "Charge ratio - Layer 6;q_{ratio};Entries;",
1115                                              100,-2.0,2.0);
1116   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeRatioLayer6, 
1117                                            fGenRecPointsOffset + 23, expert, !image);
1118   fSSDhRecPointsTask += 1;
1119   TH1F *fHistSSDChargekeVLayer5 = new TH1F("fHistSSDChargekeVLayer5",
1120                                            "Charge - Layer 5;q [keV];Entries;",
1121                                            100,0.,300.);
1122   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargekeVLayer5, 
1123                                            fGenRecPointsOffset + 24, !expert, image);
1124   fSSDhRecPointsTask += 1;
1125   TH1F *fHistSSDChargekeVLayer6 = new TH1F("fHistSSDChargekeVLayer6",
1126                                            "Charge - Layer 6;q [keV];Entries;",
1127                                            100,0.,300.);
1128   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargekeVLayer6, 
1129                                            fGenRecPointsOffset + 25, !expert, image);
1130   fSSDhRecPointsTask += 1;
1131   TH1F *fHistSSDChargePSideLayer5 = new TH1F("fHistSSDChargePSideLayer5",
1132                                              "Charge P- Layer 5;q_{P} [keV];Entries;",
1133                                              100,0.,300.);
1134   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargePSideLayer5,
1135                                            fGenRecPointsOffset + 26, expert, !image);
1136   fSSDhRecPointsTask += 1;
1137   TH1F *fHistSSDChargePSideLayer6 = new TH1F("fHistSSDChargePSideLayer6",
1138                                              "Charge P- Layer 6;q_{P} [keV];Entries;",
1139                                              100,0.,300.);
1140   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargePSideLayer6,
1141                                            fGenRecPointsOffset + 27, expert, !image);
1142   fSSDhRecPointsTask += 1;
1143   TH1F *fHistSSDChargeNSideLayer5 = new TH1F("fHistSSDChargeNSideLayer5",
1144                                              "Charge N- Layer 5;q_{N} [keV];Entries;",
1145                                              100,0.,300.);
1146   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeNSideLayer5,
1147                                            fGenRecPointsOffset + 28, expert, !image);
1148   fSSDhRecPointsTask += 1;
1149   TH1F *fHistSSDChargeNSideLayer6 = new TH1F("fHistSSDChargeNSideLayer6",
1150                                              "Charge N- Layer 6;q_{N} [keV];Entries;",
1151                                              100,0.,300.);
1152   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeNSideLayer6,
1153                                            fGenRecPointsOffset + 29, expert, !image);
1154   fSSDhRecPointsTask += 1;
1155   TH1F *fHistSSDChargeRatio2Layer5 = new TH1F("fHistSSDChargeRatio2Layer5",
1156                                               "Charge Ratio qN/qP - Layer 5;q_{N}/q_{P};Entries;",
1157                                               100,0,2);
1158   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeRatio2Layer5,
1159                                            fGenRecPointsOffset + 30, expert, !image);
1160   fSSDhRecPointsTask += 1;
1161   TH1F *fHistSSDChargeRatio2Layer6 = new TH1F("fHistSSDChargeRatio2Layer6",
1162                                               "Charge Ratio qN/qP - Layer 6;q_{N}/q_{P};Entries;",
1163                                               100,0,2);
1164   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeRatio2Layer6,
1165                                            fGenRecPointsOffset + 31, expert, !image);
1166   fSSDhRecPointsTask += 1;
1167   TH2F *fHistSSDChargePNSideLayer5 = new TH2F("fHistSSDChargePNSideLayer5",
1168                                               "Charge correlation - Layer 5;q_{P} [keV];q_{N} [keV]",
1169                                               100,0.,300.,
1170                                               100,0.,300.);
1171   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargePNSideLayer5,
1172                                            fGenRecPointsOffset + 32, expert, !image);
1173   fSSDhRecPointsTask += 1;
1174   TH2F *fHistSSDChargePNSideLayer6 = new TH2F("fHistSSDChargePNSideLayer6",
1175                                               "Charge correlation - Layer 6;q_{P} [keV];q_{N} [keV]",
1176                                               100,0.,300.,
1177                                               100,0.,300.);
1178   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargePNSideLayer6,
1179                                            fGenRecPointsOffset + 33, expert, !image);
1180   fSSDhRecPointsTask += 1;
1181   TH2F *fHistSSDChargeMapLayer5 = new TH2F("fHistSSDChargeMapLayer5",
1182                                            "Charge map;N_{modules};N_{Ladders}",
1183                                            fgkSSDMODULESPERLADDERLAYER5,
1184                                            -0.5,fgkSSDMODULESPERLADDERLAYER5+0.5,
1185                                            3*fgkSSDLADDERSLAYER5,
1186                                            -0.5,fgkSSDLADDERSLAYER5+0.5);
1187   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeMapLayer5, 
1188                                            fGenRecPointsOffset + 34, expert, !image);
1189   fSSDhRecPointsTask += 1;
1190   TH2F *fHistSSDChargeMapLayer6 = new TH2F("fHistSSDChargeMapLayer6",
1191                                            "Charge map;N_{modules};N_{Ladders}",
1192                                            fgkSSDMODULESPERLADDERLAYER6,
1193                                            -0.5,fgkSSDMODULESPERLADDERLAYER6+0.5,
1194                                            3*fgkSSDLADDERSLAYER6,
1195                                            -0.5,fgkSSDLADDERSLAYER6+0.5);
1196   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDChargeMapLayer6, 
1197                                            fGenRecPointsOffset + 35, expert, !image);
1198   fSSDhRecPointsTask += 1;
1199   TH2F *fHistSSDClusterMapLayer5 = new TH2F("fHistSSDClusterMapLayer5",
1200                                             "Layer 5;N_{module};N_{ladder}",
1201                                             22,1,23,
1202                                             34,500,534);
1203   fHistSSDClusterMapLayer5->GetXaxis()->SetTitleColor(1);
1204   fHistSSDClusterMapLayer5->SetStats(kFALSE);
1205   fHistSSDClusterMapLayer5->GetYaxis()->SetTitleOffset(1.8);
1206   fHistSSDClusterMapLayer5->GetXaxis()->SetNdivisions(22);
1207   fHistSSDClusterMapLayer5->GetYaxis()->SetNdivisions(34);
1208   fHistSSDClusterMapLayer5->GetXaxis()->SetLabelSize(0.03);
1209   fHistSSDClusterMapLayer5->GetYaxis()->SetLabelSize(0.03);
1210   fHistSSDClusterMapLayer5->GetZaxis()->SetTitleOffset(1.4);
1211   fHistSSDClusterMapLayer5->GetZaxis()->SetTitle("N_{clusters}");
1212   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDClusterMapLayer5,
1213                                            fGenRecPointsOffset + 36, !expert, image);
1214   fSSDhRecPointsTask += 1;
1215   TH2F *fHistSSDClusterMapLayer6 = new TH2F("fHistSSDClusterMapLayer6",
1216                                             "Layer 6;N_{module};N_{ladder}",
1217                                             25,1,26,
1218                                             38,600,638);
1219   fHistSSDClusterMapLayer6->GetXaxis()->SetTitleColor(1);
1220   fHistSSDClusterMapLayer6->SetStats(kFALSE);
1221   fHistSSDClusterMapLayer6->GetYaxis()->SetTitleOffset(1.8);
1222   fHistSSDClusterMapLayer6->GetXaxis()->SetNdivisions(25);
1223   fHistSSDClusterMapLayer6->GetYaxis()->SetNdivisions(38);
1224   fHistSSDClusterMapLayer6->GetXaxis()->SetLabelSize(0.03);
1225   fHistSSDClusterMapLayer6->GetYaxis()->SetLabelSize(0.03);
1226   fHistSSDClusterMapLayer6->GetZaxis()->SetTitleOffset(1.4);
1227   fHistSSDClusterMapLayer6->GetZaxis()->SetTitle("N_{clusters}");
1228   rv = fAliITSQADataMakerRec->Add2RecPointsList(fHistSSDClusterMapLayer6,
1229                                            fGenRecPointsOffset + 37, !expert, image);
1230   fSSDhRecPointsTask += 1;
1231   //printf ("%d SSD Recs histograms booked\n",fSSDhRecPointsTask);
1232   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD Recs histograms booked\n",fSSDhRecPointsTask));
1233   return rv ; 
1234 }
1235
1236 //____________________________________________________________________________ 
1237 Int_t AliITSQASSDDataMakerRec::MakeRecPoints(TTree *clustersTree)
1238 {
1239   // Fill QA for recpoints - SSD -
1240   //printf("*-*-*-*-*-*-*---*-*-*-------*-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints called \n");
1241   Int_t rv = 0 ; 
1242   Int_t gLayer = 0, gLadder = 0, gModule = 0;
1243   Int_t lLadderLocationY = 0;
1244   TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
1245   if (!branchRecP) { 
1246     AliError("can't get the branch with the ITS clusters !");
1247     return rv;
1248   }
1249
1250   // Check id histograms already created for this Event Specie
1251   if ( ! fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset) )
1252     rv = InitRecPoints() ;
1253   
1254   static TClonesArray statRecpoints("AliITSRecPoint");
1255   TClonesArray *recpoints = &statRecpoints;
1256   branchRecP->SetAddress(&recpoints);
1257   Int_t nClustersLayer5 = 0, nClustersLayer6 = 0;
1258   Int_t npoints = 0;      
1259   Float_t cluglo[3]={0.,0.,0.}; 
1260   //printf("*-*-*-*-*-*-*---*-*-*-------*-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints STEP1 \n");
1261   for(Int_t module = 0; module < clustersTree->GetEntries(); module++){
1262     branchRecP->GetEvent(module);
1263     npoints += recpoints->GetEntries();
1264     AliITSgeomTGeo::GetModuleId(module,gLayer,gLadder,gModule);
1265     //printf("SSDDataMAkerRec:::::::::::::::::::::::gLayer ========== %d \n\n",gLayer);
1266     lLadderLocationY = 3*gLadder;
1267     ////printf("*-*-*-*-*-*-*---*-*-*-------*-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints inside loop \n");
1268     for(Int_t j = 0;j < recpoints->GetEntries(); j++){
1269       ////printf("*-*-*-*-*-*-*---*-*-*-------*-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints inside loop 2\n");
1270       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);    
1271       Int_t layer = recp->GetLayer();
1272       //printf("SSDDataMAkerRec:::::::::::::::::::::::layer ========== %d \n\n",layer);
1273       recp->GetGlobalXYZ(cluglo);
1274       Float_t radius = TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); 
1275       Float_t phi = TMath::ATan2(cluglo[1],cluglo[0]);
1276       Float_t theta = TMath::ATan2(radius,cluglo[2]);
1277       Double_t chargeRatio = recp->GetChargeRatio();
1278       Double_t clusterCharge = recp->GetQ();
1279       Double_t chargePSide = clusterCharge*(1. + chargeRatio);
1280       Double_t chargeNSide = clusterCharge*(1. - chargeRatio);
1281       if(layer == 4) {
1282         //printf("-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints Filling 4 called \n");
1283         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 20)->Fill(recp->GetType());
1284
1285         if(recp->GetType() != 1) continue;
1286         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 0)->Fill(module);
1287         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 4)->Fill(recp->GetDetLocalX());
1288         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 6)->Fill(recp->GetDetLocalZ());
1289         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 8)->Fill(cluglo[0]);
1290         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 10)->Fill(cluglo[1]);
1291         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 12)->Fill(cluglo[2]);
1292         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 14)->Fill(phi);
1293         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 16)->Fill(theta);
1294         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 18)->Fill(radius);
1295         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 22)->Fill(recp->GetChargeRatio());
1296         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 24)->Fill(recp->GetQ());
1297         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 26)->Fill(chargePSide);
1298         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 28)->Fill(chargeNSide);
1299         if(chargePSide != 0.) fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 30)->Fill(chargeNSide/chargePSide);
1300         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 32)->Fill(chargePSide,chargeNSide);
1301         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 34)->SetBinContent(gModule,lLadderLocationY,recp->GetQ());
1302         ((TH2F *)fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 36))->Fill(gModule,499+gLadder,1);
1303         nClustersLayer5 += 1;
1304       }//layer 5 histograms
1305       if(layer == 5) {
1306         //printf("-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints Filling 5 called \n");
1307         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 21)->Fill(recp->GetType());
1308
1309         if(recp->GetType() != 1) continue;
1310         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 1)->Fill(module);
1311         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 5)->Fill(recp->GetDetLocalX());
1312         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 7)->Fill(recp->GetDetLocalZ());
1313         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 9)->Fill(cluglo[0]);
1314         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 11)->Fill(cluglo[1]);
1315         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 13)->Fill(cluglo[2]);
1316         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 15)->Fill(phi);
1317         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 17)->Fill(theta);
1318         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 19)->Fill(radius);
1319         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 23)->Fill(recp->GetChargeRatio());
1320         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 25)->Fill(recp->GetQ());
1321         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 27)->Fill(chargePSide);
1322         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 29)->Fill(chargeNSide);
1323         if(chargePSide != 0.) fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 31)->Fill(chargeNSide/chargePSide);
1324         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 33)->Fill(chargePSide,chargeNSide);
1325         fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 35)->SetBinContent(gModule,lLadderLocationY,recp->GetQ());
1326         ((TH2F *)fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 37))->Fill(gModule,599+gLadder,1);
1327         nClustersLayer6 += 1;
1328       }//layer 6 histograms
1329     }//rec. points loop
1330   }//module loop
1331   
1332   //printf("-*-*-*-*-*-***************AliITSQASSDataMakerRec::MakeRecpoints Filling called \n");
1333   fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 2)->Fill(nClustersLayer5);
1334   fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset + 3)->Fill(nClustersLayer6);
1335
1336   statRecpoints.Clear();
1337   return rv ; 
1338 }
1339
1340 //____________________________________________________________________________ 
1341 Int_t AliITSQASSDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task) {
1342   // Returns offset number according to the specified task 
1343   Int_t offset=0;
1344   if( task == AliQAv1::kRAWS ) {
1345     offset=fGenRawsOffset;  
1346   }
1347   else if( task == AliQAv1::kDIGITSR ) {
1348     offset=fGenDigitsOffset;   
1349   }
1350   else if( task == AliQAv1::kRECPOINTS ) {
1351     offset=fGenRecPointsOffset;   
1352   }
1353   else {
1354     AliWarning("No task has been selected. Offset set to zero.\n");
1355   }
1356
1357   return offset;
1358 }
1359
1360 //_______________________________________________________________
1361
1362 void AliITSQASSDDataMakerRec::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset) {
1363   // Returns offset number according to the specified task
1364   if( task == AliQAv1::kRAWS ) {
1365     fGenRawsOffset=offset;
1366   }
1367   else if( task == AliQAv1::kDIGITSR ) {
1368     fGenDigitsOffset=offset;
1369   }
1370   else if( task == AliQAv1::kRECPOINTS ) {
1371     fGenRecPointsOffset=offset;
1372   }
1373   else {
1374     AliInfo("No task has been selected. Offset set to zero.\n");
1375   }
1376 }
1377
1378 //____________________________________________________________________________ 
1379 Int_t AliITSQASSDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task) {
1380   // Returns the number of histograms associated to the specified task
1381   Int_t histotot=0;
1382
1383   if( task == AliQAv1::kRAWS ) {
1384     histotot=fSSDhRawsTask;  
1385   }
1386   else if( task == AliQAv1::kDIGITSR ) {
1387     histotot=fSSDhDigitsTask;
1388   }
1389   else if( task == AliQAv1::kRECPOINTS ){
1390     histotot=fSSDhRecPointsTask;   
1391   }
1392   else { 
1393     AliWarning("No task has been selected. TaskHisto set to zero.\n");
1394   }
1395
1396   return histotot;
1397 }