]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerRec.cxx
Fixes for bug #61900: Restore capaibility to chooses ITS PID method (Marco)
[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 
25 //  INFN Torino
26
27 // --- ROOT system ---
28
29 #include <TProfile2D.h>
30 #include <TH2D.h>
31 #include <TH1F.h>
32 #include <TBranch.h>
33 #include <TTree.h>
34 #include <TGaxis.h>
35 #include <TMath.h>
36 #include <TF1.h>
37 #include <TDirectory.h>
38 #include <TSystem.h>
39 // --- Standard library ---
40
41 // --- AliRoot header files ---
42 #include "AliITSQASDDDataMakerRec.h"
43 #include "AliLog.h"
44 #include "AliQAv1.h"
45 #include "AliQAChecker.h"
46 #include "AliRawReader.h"
47 #include "AliITSRawStream.h"
48 #include "AliITSRawStreamSDD.h"
49 #include "AliITSRawStreamSDDCompressed.h"
50 #include "AliITSDetTypeRec.h"
51 #include "AliITSdigit.h"
52 #include "AliITSRecPoint.h"
53 #include "AliITSRecPointContainer.h"
54 #include "AliITSgeomTGeo.h"
55 #include "AliCDBManager.h"
56 #include "AliCDBStorage.h"
57 #include "AliCDBEntry.h"
58 #include "Riostream.h"
59 #include "AliITSdigitSDD.h"
60 #include "AliITS.h"
61 #include "AliRunLoader.h"
62 #include "AliITSLoader.h"
63 #include "AliITSDetTypeRec.h"
64
65
66
67 ClassImp(AliITSQASDDDataMakerRec)
68
69 //____________________________________________________________________________ 
70 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) :
71 TObject(),
72 fAliITSQADataMakerRec(aliITSQADataMakerRec),
73 fkOnline(kMode),
74 fLDC(ldc),
75 fSDDhRawsTask(0),
76 fSDDhDigitsTask(0),
77 fSDDhRecPointsTask(0),
78 fGenRawsOffset(0),
79 fGenDigitsOffset(0),
80 fGenRecPointsOffset(0),
81 fTimeBinSize(1),
82 fDDLModuleMap(0),
83 fGoodAnodes(0),
84 fBadAnodes(0),
85 fGoodAnodesCurrent(0),
86 fBadAnodesCurrent(0)
87 {
88   //ctor used to discriminate OnLine-Offline analysis
89   if(fLDC < 0 || fLDC > 6) {
90         AliError("Error: LDC number out of range; return\n");
91   }
92         fGenRawsOffset = new Int_t[AliRecoParam::kNSpecies];
93         fGenRecPointsOffset = new Int_t[AliRecoParam::kNSpecies];
94         fGenDigitsOffset = new Int_t[AliRecoParam::kNSpecies];
95         for(Int_t i=0; i<AliRecoParam::kNSpecies; i++) {
96                 fGenRawsOffset[i] = 0;
97                 fGenRecPointsOffset[i] = 0;
98                 fGenDigitsOffset[i]=0;
99         }
100 }
101
102 //____________________________________________________________________________ 
103 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) :
104 TObject(),
105 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
106 fkOnline(qadm.fkOnline),
107 fLDC(qadm.fLDC),
108 fSDDhRawsTask(qadm.fSDDhRawsTask),
109 fSDDhDigitsTask(qadm.fSDDhDigitsTask),
110 fSDDhRecPointsTask(qadm.fSDDhRecPointsTask),
111 fGenRawsOffset(qadm.fGenRawsOffset),
112 fGenDigitsOffset(qadm.fGenDigitsOffset),
113 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
114 fTimeBinSize(1),
115 fDDLModuleMap(0),
116 fGoodAnodes(qadm.fGoodAnodes),
117 fBadAnodes(qadm.fBadAnodes),
118 fGoodAnodesCurrent(qadm.fGoodAnodesCurrent),
119 fBadAnodesCurrent(qadm.fBadAnodesCurrent)
120 {
121   //copy ctor 
122   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; 
123   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
124   fDDLModuleMap=NULL;
125 }
126
127 //____________________________________________________________________________ 
128 AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){
129   // destructor
130   //if(fDDLModuleMap) delete fDDLModuleMap;
131 }
132 //__________________________________________________________________
133 AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac )
134 {
135   // Equal operator.
136   this->~AliITSQASDDDataMakerRec();
137   new(this) AliITSQASDDDataMakerRec(qac);
138   return *this;
139 }
140
141 //____________________________________________________________________________ 
142 void AliITSQASDDDataMakerRec::StartOfDetectorCycle()
143 {
144   //Detector specific actions at start of cycle
145   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SDD Cycle\n");
146 }
147
148 //____________________________________________________________________________ 
149 void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
150 {
151   
152   // launch the QA checking
153         if(fkOnline) {
154                 AnalyseBNG(); // Analyse Baseline, Noise, Gain
155                 AnalyseINJ(); // Analyse Injectors
156         }
157   
158         AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); 
159 }
160
161 //____________________________________________________________________________ 
162 Int_t AliITSQASDDDataMakerRec::InitRaws()
163
164   // Initialization for RAW data - SDD -
165   const Bool_t expert   = kTRUE ; 
166   const Bool_t saveCorr = kTRUE ; 
167   const Bool_t image    = kTRUE ; 
168
169
170
171   Int_t rv = 0 ; 
172   /*
173   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
174   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
175   if(!ddlMapSDD)
176     {
177       AliError("Calibration object retrieval failed! SDD will not be processed");
178       fDDLModuleMap = NULL;
179       return rv;
180     }
181   fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
182   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
183   ddlMapSDD->SetOwner(kTRUE);
184   if(!cacheStatus)
185     {
186       delete ddlMapSDD;
187     }
188   */
189   Int_t lay, lad, det;
190   Int_t indexlast = 0;
191   Int_t index1 = 0;
192
193   if(fkOnline) 
194     {
195       AliInfo("Book Online Histograms for SDD\n");
196     }
197   else 
198     {
199       AliInfo("Book Offline Histograms for SDD\n ");
200     }
201   TH1D *h0 = new TH1D("SDDModPattern","HW Modules pattern",fgknSDDmodules,239.5,499.5); //0
202   h0->GetXaxis()->SetTitle("Module Number");
203   h0->GetYaxis()->SetTitle("Counts");
204   rv = fAliITSQADataMakerRec->Add2RawsList(h0,0+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
205   fSDDhRawsTask++;
206   
207   //zPhi distribution using ladder and modules numbers
208   TH2D *hphil3 = new TH2D("SDDphizL3","SDD #varphiz Layer3 ",12,0.5,6.5,14,0.5,14.5);
209   hphil3->GetXaxis()->SetTitle("z[#Module L3 ]");
210   hphil3->GetYaxis()->SetTitle("#varphi[#Ladder L3]");
211   rv = fAliITSQADataMakerRec->Add2RawsList(hphil3,1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, saveCorr); 
212   fSDDhRawsTask++;
213   
214   TH2D *hphil4 = new TH2D("SDDphizL4","SDD #varphiz Layer4 ",16,0.5,8.5,22,0.5,22.5); 
215   hphil4->GetXaxis()->SetTitle("z[#Module L4]");
216   hphil4->GetYaxis()->SetTitle("#varphi[#Ladder L4]");
217    rv = fAliITSQADataMakerRec->Add2RawsList(hphil4,2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, saveCorr); 
218   fSDDhRawsTask++;
219   
220
221   if(fkOnline) 
222     {
223
224       //DDL Pattern 
225       TH2D *hddl = new TH2D("SDDDDLPattern","SDD DDL Pattern ",24,-0.5,11.5,24,-0.5,23.5); 
226       hddl->GetXaxis()->SetTitle("Channel");
227       hddl->GetYaxis()->SetTitle("#DDL");
228       rv = fAliITSQADataMakerRec->Add2RawsList(hddl,3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
229       fSDDhRawsTask++;
230       Int_t indexlast1 = 0;
231   
232       fTimeBinSize = 4;
233       indexlast = 0;
234       index1 = 0;
235       indexlast1 = fSDDhRawsTask;
236       char *hname[3];
237       for(Int_t i=0; i<3; i++) hname[i]= new char[50];
238       for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
239         for(Int_t iside=0;iside<fgknSide;iside++){
240           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
241           sprintf(hname[0],"SDDchargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside);
242           sprintf(hname[1],"SDDChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside);
243           //      sprintf(hname[2],"SDDhmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);
244           TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
245           fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin");
246           fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode");
247            rv = fAliITSQADataMakerRec->Add2RawsList(fModuleChargeMapFSE,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);        
248           fSDDhRawsTask++;
249           index1++;      
250         }
251       }
252       
253       for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
254         for(Int_t iside=0;iside<fgknSide;iside++){
255           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
256           sprintf(hname[0],"SDDchargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
257           sprintf(hname[1],"SDDChargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
258           TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
259           fModuleChargeMap->GetXaxis()->SetTitle("Time Bin");
260           fModuleChargeMap->GetYaxis()->SetTitle("Anode");
261            rv = fAliITSQADataMakerRec->Add2RawsList(fModuleChargeMap,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
262   
263           fSDDhRawsTask++;
264           index1++;      
265         }
266       }
267       
268       //Event Size 
269       TH1F *hsize = new TH1F("SDDEventSize","SDD Event Size ",1000,-0.5,999.5); 
270       hsize->GetXaxis()->SetTitle("Event Size [kB]");
271       hsize->GetYaxis()->SetTitle("Entries");
272       rv = fAliITSQADataMakerRec->Add2RawsList(hsize,indexlast1 + index1 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
273       fSDDhRawsTask++;
274           
275     }  // kONLINE
276   
277   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Raws histograms booked\n",fSDDhRawsTask));
278   return rv ; 
279 }
280
281
282 //____________________________________________________________________________
283 Int_t AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader)
284
285   // Fill QA for RAW - SDD -
286         Int_t rv = 0;
287   // Check id histograms already created for this Event Specie
288
289   if(!fDDLModuleMap){
290     CreateTheMap();
291     //AliError("SDD DDL module map not available - skipping SDD QA");
292     //return rv;
293   
294   }
295   if(rawReader->GetType() != 7) return rv;  // skips non physical triggers
296   AliDebug(AliQAv1::GetQADebugLevel(),"entering MakeRaws\n");                 
297   rawReader->Reset();       
298   rawReader->Select("ITSSDD");
299   AliITSRawStream *stream=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader);
300    stream->SetDDLModuleMap(fDDLModuleMap);
301   
302   Int_t lay, lad, det; 
303   
304   Int_t index=0;
305   if(fkOnline) {
306     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
307       for(Int_t iside=0;iside<fgknSide;iside++) {
308                 if(fSDDhRawsTask > 4 + index) fAliITSQADataMakerRec->GetRawsData(4 + index +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();   
309         // 4  because the 2D histos for single events start after the fourth position
310         index++;
311       }
312     }
313   }
314   
315   Int_t cnt = 0;
316   Int_t ildcID = -1;
317   Int_t iddl = -1;
318   Int_t isddmod = -1;
319   Int_t coord1, coord2, signal, moduleSDD, activeModule, index1; 
320   //if(fkOnline)
321   //{
322       Int_t prevDDLID = -1;
323       UInt_t size = 0;
324       int totalddl=static_cast<int>(fDDLModuleMap->GetNDDLs());
325       Bool_t *ddldata=new Bool_t[totalddl];
326       for(Int_t jddl=0;jddl<totalddl;jddl++){ddldata[jddl]=kFALSE;}
327       //}
328   while(stream->Next()) {
329     ildcID = rawReader->GetLDCId();
330     iddl = rawReader->GetDDLID();// - fgkDDLIDshift;
331
332
333     isddmod = fDDLModuleMap->GetModuleNumber(iddl,stream->GetCarlosId());
334     if(isddmod==-1){
335       AliDebug(AliQAv1::GetQADebugLevel(),Form("Found module with iddl: %d, stream->GetCarlosId: %d \n",iddl,stream->GetCarlosId()));
336       continue;
337     }
338     if(stream->IsCompletedModule()) {
339       AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedModule == KTRUE\n"));
340       continue;
341     } 
342     if(stream->IsCompletedDDL()) {
343
344       if(fkOnline){
345         if ((rawReader->GetDDLID() != prevDDLID)&&(ddldata[iddl])==kFALSE){
346           size += rawReader->GetDataSize();//in bytes
347           prevDDLID = rawReader->GetDDLID();
348           ddldata[iddl]=kTRUE;
349         }
350       }
351       AliDebug(AliQAv1::GetQADebugLevel(),Form("IsCompletedDDL == KTRUE\n"));
352       continue;
353     } 
354     
355     coord1 = stream->GetCoord1();
356     coord2 = stream->GetCoord2();
357     signal = stream->GetSignal();
358     
359     moduleSDD = isddmod - fgkmodoffset;
360     
361     if(isddmod <fgkmodoffset|| isddmod>fgknSDDmodules+fgkmodoffset-1) {
362       AliDebug(AliQAv1::GetQADebugLevel(),Form( "Module SDD = %d, resetting it to 1 \n",isddmod));
363       isddmod = 1;
364     }
365     
366     AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);
367     Short_t iside = stream->GetChannel();
368
369     //printf(" \n%i %i %i %i \n ",lay, lad, det,iside );
370     fAliITSQADataMakerRec->GetRawsData( 0 + fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()] )->Fill(isddmod);   
371     if(lay==3)    fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad); 
372     if(lay==4) { 
373       fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(det+0.5*iside-0.5,lad);}  
374  
375     if(fkOnline) {
376       fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill((stream->GetCarlosId())+0.5*iside -0.5,iddl);
377
378       activeModule = moduleSDD;
379       index1 = activeModule * 2 + iside;
380       
381       if(index1<0){
382         AliDebug(AliQAv1::GetQADebugLevel(),Form("Wrong index number %d - patched to 0\n",index1));
383         index1 = 0;
384       }      
385
386       if(fSDDhRawsTask > 4 + index1) {                                  
387         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(4 + index1 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(coord2, coord1, signal);     
388         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(4 + index1 + 260*2 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(coord2, coord1, signal); 
389       }
390
391     }//online
392     cnt++;
393     if(!(cnt%10000)) AliDebug(AliQAv1::GetQADebugLevel(),Form(" %d raw digits read",cnt));
394   }//end next()
395   if(fkOnline)
396     {
397       ((TH1F*)(fAliITSQADataMakerRec->GetRawsData(4 + 260*4 +fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])))->Fill(size/1024.);//KB
398     }
399   AliDebug(AliQAv1::GetQADebugLevel(),Form("Event completed, %d raw digits read",cnt)); 
400   delete stream;
401   stream = NULL; 
402
403 //      if(fkOnline) {
404 //              AnalyseBNG(); // Analyse Baseline, Noise, Gain
405 //              AnalyseINJ(); // Analyse Injectors
406 //      }
407
408   delete []ddldata;
409   ddldata=NULL;
410   return rv ; 
411 }
412
413 //____________________________________________________________________________ 
414 Int_t AliITSQASDDDataMakerRec::InitDigits()
415
416
417
418   // Initialization for DIGIT data - SDD -  
419   const Bool_t expert   = kTRUE ; 
420   const Bool_t image    = kTRUE ;
421   Int_t rv = 0 ; 
422 //  fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
423   //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List
424   TH1F* h0=new TH1F("SDD DIGITS Module Pattern","SDD DIGITS Module Pattern",260,239.5,499.5);       //hmod
425   h0->GetXaxis()->SetTitle("SDD Module Number");
426   h0->GetYaxis()->SetTitle("# DIGITS");
427   rv = fAliITSQADataMakerRec->Add2DigitsList(h0,fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
428   fSDDhDigitsTask ++;
429   // printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
430   TH1F* h1=new TH1F("SDD Anode Distribution","DIGITS Anode Distribution",512,-0.5,511.5);      //hanocc
431   h1->GetXaxis()->SetTitle("Anode Number");
432   h1->GetYaxis()->SetTitle("# DIGITS");
433   rv = fAliITSQADataMakerRec->Add2DigitsList(h1,1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
434   fSDDhDigitsTask ++;
435   //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
436   TH1F* h2=new TH1F("SDD Tbin Distribution","DIGITS Tbin Distribution",256,-0.5,255.5);      //htbocc
437   h2->GetXaxis()->SetTitle("Tbin Number");
438   h2->GetYaxis()->SetTitle("# DIGITS");
439   rv = fAliITSQADataMakerRec->Add2DigitsList(h2,2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
440   fSDDhDigitsTask ++;
441   //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
442   TH1F* h3=new TH1F("SDD ADC Counts Distribution","DIGITS ADC Counts Distribution",200,0.,1024.);          //hsig
443   h3->GetXaxis()->SetTitle("ADC Value");
444   h3->GetYaxis()->SetTitle("# DIGITS");
445   rv = fAliITSQADataMakerRec->Add2DigitsList(h3,3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
446   fSDDhDigitsTask ++;
447   //printf("Add %s \t the task offset is %i\n",fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->GetName() , fSDDhDigitsTask );
448   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Digits histograms booked\n",fSDDhDigitsTask));
449   return rv ; 
450 }
451
452 //____________________________________________________________________________
453 Int_t AliITSQASDDDataMakerRec::MakeDigits(TTree * digits)
454
455
456   // Fill QA for DIGIT - SDD -
457   //AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");
458   //fITS->SetTreeAddress();
459   //TClonesArray *iITSdigits  = fITS->DigitsAddress(1);
460
461
462   Int_t rv = 0 ; 
463
464   TBranch *branchD = digits->GetBranch("ITSDigitsSDD");
465
466   if (!branchD) {
467     AliError("can't get the branch with the ITS SDD digits !");
468     return rv ;
469   }
470   // Check id histograms already created for this Event Specie
471 //  if ( ! fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset) )
472 //    rv = InitDigits() ;
473   
474   static TClonesArray statDigits("AliITSdigitSDD");
475   TClonesArray *iITSdigits = &statDigits;
476   branchD->SetAddress(&iITSdigits);
477
478   for(Int_t i=0; i<260; i++){
479     Int_t nmod=i+240;
480     digits->GetEvent(nmod);
481     Int_t ndigits = iITSdigits->GetEntries();
482     fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nmod,ndigits);
483
484     for (Int_t idig=0; idig<ndigits; idig++) {
485       AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
486       Int_t iz=dig->GetCoord1();  // cell number z
487       Int_t ix=dig->GetCoord2();  // cell number x
488       Int_t sig=dig->GetSignal();
489       fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(iz);
490       fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(ix);
491       fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(sig);
492     }
493   }
494   return rv ; 
495 }
496
497 //____________________________________________________________________________ 
498 Int_t AliITSQASDDDataMakerRec::InitRecPoints()
499 {
500
501         //AliInfo("Initialize SDD recpoints histos\n");
502   // Initialization for RECPOINTS - SDD -
503   const Bool_t expert   = kTRUE ; 
504   const Bool_t image    = kTRUE ; 
505   Int_t rv = 0 ; 
506 //  fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
507
508   Int_t  nOnline=1;
509   Int_t  nOnline2=1;
510   Int_t  nOnline3=1; 
511   Int_t  nOnline4=1;
512   Int_t nOnline5=1;
513   if(fkOnline)
514     {
515       nOnline=4;
516       nOnline2=28;
517       nOnline3=64;
518       nOnline4=14;
519     }
520
521   //AliInfo(Form("fAliITSQADataMakerRec->GetEventSpecie() %d\n",fAliITSQADataMakerRec->GetEventSpecie()));
522   //AliInfo(Form("fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()] %d\n",fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]));
523   TH1F *h0 = new TH1F("SDDLay3TotCh","Layer 3 total charge",1000/nOnline,-0.5, 499.5); //position number 0
524   h0->GetXaxis()->SetTitle("ADC value");
525   h0->GetYaxis()->SetTitle("Entries");
526   rv = fAliITSQADataMakerRec->Add2RecPointsList(h0, 0 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
527   fSDDhRecPointsTask++;
528  
529   TH1F *h1 = new TH1F("SDDLay4TotCh","Layer 4 total charge",1000/nOnline,-0.5, 499.5);//position number 1
530   h1->GetXaxis()->SetTitle("ADC value");
531   h1->GetYaxis()->SetTitle("Entries");
532   rv = fAliITSQADataMakerRec->Add2RecPointsList(h1, 1 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
533   fSDDhRecPointsTask++;
534
535   char hisnam[50];
536   TH2F *h2 = new TH2F("SDDGlobalCoordDistribYX","YX Global Coord Distrib",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 2
537   h2->GetYaxis()->SetTitle("Y[cm]");
538   h2->GetXaxis()->SetTitle("X[cm]");
539   rv = fAliITSQADataMakerRec->Add2RecPointsList(h2,2+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image
540   fSDDhRecPointsTask++;
541
542   TH2F *h3 = new TH2F("SDDGlobalCoordDistribRZ","RZ Global Coord Distrib",6400/nOnline3,-32,32,1400/nOnline4,12,26);//position number 3
543   h3->GetYaxis()->SetTitle("R[cm]");
544   h3->GetXaxis()->SetTitle("Z[cm]");
545   rv = fAliITSQADataMakerRec->Add2RecPointsList(h3,3+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert image
546   fSDDhRecPointsTask++;
547   
548   TH2F *h4 = new TH2F("SDDGlobalCoordDistribL3PHIZ","#varphi Z Global Coord Distrib L3",6400/nOnline3,-32,32,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 4
549   h4->GetYaxis()->SetTitle("#phi[rad]");
550   h4->GetXaxis()->SetTitle("Z[cm]");
551   rv = fAliITSQADataMakerRec->Add2RecPointsList(h4,4+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
552   fSDDhRecPointsTask++;
553
554   TH2F *h5 = new TH2F("SDDGlobalCoordDistribL4PHIZ","#varphi Z Global Coord Distrib L4",6400/nOnline3,-32,32,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 5
555   h5->GetYaxis()->SetTitle("#phi[rad]");
556   h5->GetXaxis()->SetTitle("Z[cm]");
557   rv = fAliITSQADataMakerRec->Add2RecPointsList(h5,5+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);//NON expert image
558   fSDDhRecPointsTask++;
559   
560   TH1F *h6 = new TH1F("SDDModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); //position number 6
561   h6->GetXaxis()->SetTitle("Module number"); //spd offset = 240
562   h6->GetYaxis()->SetTitle("Entries");
563   rv = fAliITSQADataMakerRec->Add2RecPointsList(h6,6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
564   fSDDhRecPointsTask++;
565
566   TH1F *h7 = new TH1F("SDDLadPatternL3RP","Ladder pattern L3 RP",14,0.5,14.5);  //position number 7
567   h7->GetXaxis()->SetTitle("Ladder #, Layer 3");
568   h7->GetYaxis()->SetTitle("Entries");
569   rv = fAliITSQADataMakerRec->Add2RecPointsList(h7,7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
570   fSDDhRecPointsTask++;
571
572   TH1F *h8 = new TH1F("SDDLadPatternL4RP","Ladder pattern L4 RP",22,0.5,22.5); //position number 8
573   h8->GetXaxis()->SetTitle("Ladder #, Layer 4");
574   h8->GetYaxis()->SetTitle("Entries");
575   rv = fAliITSQADataMakerRec->Add2RecPointsList(h8,8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);//  expert  NO image
576   fSDDhRecPointsTask++;
577
578   TH2F *h9 = new TH2F("SDDLocalCoordDistrib","Local Coord Distrib",1000/nOnline,-4,4,1000/nOnline,-4,4);//position number 9
579   h9->GetXaxis()->SetTitle("X local coord, drift, cm");
580   h9->GetYaxis()->SetTitle("Z local coord, anode, cm");
581   rv = fAliITSQADataMakerRec->Add2RecPointsList(h9,9 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);//  expert  NO image
582   fSDDhRecPointsTask++;
583   
584   //AliInfo("Create SDD recpoints histos\n");
585   
586   TH1F *h10 = new TH1F("SDDrdistrib_Layer3" ,"SDD r distribution Layer3" ,100,14.,18.);//position number 10 (L3)
587   h10->GetXaxis()->SetTitle("r[cm]");
588   h10->GetXaxis()->CenterTitle();
589   h10->GetYaxis()->SetTitle("Entries");
590   rv = fAliITSQADataMakerRec->Add2RecPointsList(h10,10 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
591   fSDDhRecPointsTask++;
592   
593   TH1F *h11 = new TH1F("SDDrdistrib_Layer4" ,"SDD r distribution Layer4" ,100,22.,26.);// and position number 11 (L4)
594   h11->GetXaxis()->SetTitle("r[cm]");
595   h11->GetXaxis()->CenterTitle();
596   h11->GetYaxis()->SetTitle("Entries");
597   rv = fAliITSQADataMakerRec->Add2RecPointsList(h11,11 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
598   fSDDhRecPointsTask++;
599   
600   for(Int_t iLay=0; iLay<=1; iLay++){
601     sprintf(hisnam,"SDDphidistrib_Layer%d",iLay+3);
602     TH1F *h12 = new TH1F(hisnam,hisnam,180,-TMath::Pi(),TMath::Pi());//position number 12 (L3) and position number 13 (L4)
603     h12->GetXaxis()->SetTitle("#varphi[rad]");
604     h12->GetXaxis()->CenterTitle();
605     h12->GetYaxis()->SetTitle("Entries");
606     rv = fAliITSQADataMakerRec->Add2RecPointsList(h12,iLay+12+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
607     fSDDhRecPointsTask++;
608   }
609   
610   for(Int_t iLay=0; iLay<=1; iLay++){
611     sprintf(hisnam,"SDDdrifttime_Layer%d",iLay+3);
612     TH1F *h14 = new TH1F(hisnam,hisnam,200/nOnline5,-0.5,9999.5);//position number 14 (L3) and position number 15 (L4)
613     h14->GetXaxis()->SetTitle("drift time[#mus]");
614     h14->GetXaxis()->CenterTitle();
615     h14->GetYaxis()->SetTitle("Entries");
616     rv = fAliITSQADataMakerRec->Add2RecPointsList(h14,iLay+14+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);// NON expert  image
617     fSDDhRecPointsTask++;
618   }
619   
620   if(fkOnline)
621     {
622       TH2F *h16 = new TH2F("SDDGlobalCoordDistribYXFSE","YX Global Coord Distrib FSE",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 16
623       h16->GetYaxis()->SetTitle("Y[cm]");
624       h16->GetXaxis()->SetTitle("X[cm]");
625       rv = fAliITSQADataMakerRec->Add2RecPointsList(h16,16+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
626       fSDDhRecPointsTask++;
627       
628       TH2F *h17 = new TH2F("SDDGlobalCoordDistribRZFSE","RZ Global Coord Distrib FSE",Int_t(6400/nOnline3),-32,32,1400/nOnline4,12,26);//position number 17
629       h17->GetYaxis()->SetTitle("R[cm]");
630       h17->GetXaxis()->SetTitle("Z[cm]");
631       rv = fAliITSQADataMakerRec->Add2RecPointsList(h17,17+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);// expert NO image
632       fSDDhRecPointsTask++;
633       
634     }//online
635   
636   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Recs histograms booked\n",fSDDhRecPointsTask));
637   
638   return rv ; 
639 }
640
641 //____________________________________________________________________________ 
642 Int_t AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)
643 {
644  // Fill QA for RecPoints - SDD -
645   Int_t rv = 0 ; 
646
647   Int_t lay, lad, det; 
648   //AliInfo("get the branch with the ITS clusters !\n");
649   AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
650   TClonesArray *recpoints=NULL; 
651   if(fkOnline){
652     recpoints = rpcont->FetchClusters(0,clustersTree,fAliITSQADataMakerRec->GetEventNumber());
653   }
654   else{
655     recpoints = rpcont->FetchClusters(0,clustersTree); 
656   }
657   if(!rpcont->GetStatusOK() || !rpcont->IsSDDActive()){
658     AliError("can't get SDD clusters !");
659     return rv;
660   }
661
662   Int_t npoints = 0;      
663   Float_t cluglo[3]={0.,0.,0.}; 
664   if(fkOnline)
665     {
666       for(Int_t i=16;i<18;i++)
667         {
668           fAliITSQADataMakerRec->GetRecPointsData(i+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Reset();
669         }
670     }
671   Int_t firMod=AliITSgeomTGeo::GetModuleIndex(3,1,1);
672   Int_t lasMod=AliITSgeomTGeo::GetModuleIndex(5,1,1);
673   for(Int_t module=firMod; module<lasMod;module++){
674     //AliInfo(Form("Module %d\n",module));
675     recpoints = rpcont->UncheckedGetClusters(module);
676     npoints += recpoints->GetEntries();
677     //AliInfo(Form("modnumb %d, npoints %d, total points %d\n",module, recpoints->GetEntries(),npoints));
678     //AliITSgeomTGeo::GetModuleId(module, lay, lad, det);
679     //AliInfo(Form("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det));
680     
681     //AliInfo(Form("modnumb %d, entries %d\n",module, recpoints->GetEntries()));
682     for(Int_t j=0;j<recpoints->GetEntries();j++){
683       //AliInfo(Form("modnumb %d, entry %d \n",module, j));
684       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j); 
685       Int_t index = recp->GetDetectorIndex();
686       lay=recp->GetLayer();
687       Int_t modnumb=index+AliITSgeomTGeo::GetModuleIndex(lay+1,1,1);
688       //printf("modnumb  %i\n",modnumb);  
689       AliITSgeomTGeo::GetModuleId(modnumb, lay, lad, det);  
690       fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(modnumb);//modpatternrp
691       recp->GetGlobalXYZ(cluglo);
692       Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); 
693       Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);
694       Float_t drifttime=recp->GetDriftTime();
695       fAliITSQADataMakerRec->GetRecPointsData(9 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());//local distribution
696       if(lay==3||lay==4)fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[0],cluglo[1]);//global distribution YX
697       fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],rad);//global distribution rz
698       if(fkOnline) {
699         fAliITSQADataMakerRec->GetRecPointsData(16 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[0],cluglo[1]);//global distribution YX FSE
700         fAliITSQADataMakerRec->GetRecPointsData(17 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],rad);//global distribution rz FSE
701       }
702       if(recp->GetLayer() == 2) {
703         fAliITSQADataMakerRec->GetRecPointsData(0  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetQ()) ;//total charge of layer 3
704         fAliITSQADataMakerRec->GetRecPointsData(7  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lad);//lad pattern layer 3
705         fAliITSQADataMakerRec->GetRecPointsData(10 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);//r distribution layer 3
706         fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);// phi distribution layer 3
707         fAliITSQADataMakerRec->GetRecPointsData(4  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],phi);// zphi distribution layer
708         fAliITSQADataMakerRec->GetRecPointsData(14  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(drifttime);// time distribution layer 3
709       } else if(recp->GetLayer() == 3) {
710         fAliITSQADataMakerRec->GetRecPointsData(1  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(recp->GetQ()) ;//total charge layer 4
711         fAliITSQADataMakerRec->GetRecPointsData(8  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lad);//ladpatternlayer4
712         fAliITSQADataMakerRec->GetRecPointsData(11 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);//r distribution
713         fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);//phi distribution
714         fAliITSQADataMakerRec->GetRecPointsData(5  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluglo[2],phi);// zphi distribution layer 4
715         fAliITSQADataMakerRec->GetRecPointsData(15  +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(drifttime);// time distribution layer 4
716       }
717     }
718   }
719   return rv ; 
720 }
721
722 //_______________________________________________________________
723
724 Int_t AliITSQASDDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task, Int_t specie)
725 {
726   Int_t offset=0;
727   if( task == AliQAv1::kRAWS )
728     {
729       offset=fGenRawsOffset[specie];  
730     }
731   else if(task == AliQAv1::kDIGITSR )
732     {
733       offset=fGenDigitsOffset[specie];
734     }
735   else if( task == AliQAv1::kRECPOINTS )
736     {
737       offset=fGenRecPointsOffset[specie];   
738     }
739   return offset;
740 }
741
742 //_______________________________________________________________
743
744 void AliITSQASDDDataMakerRec::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset, Int_t specie) {
745   // Returns offset number according to the specified task
746   if( task == AliQAv1::kRAWS ) {
747     fGenRawsOffset[specie]=offset;
748   }
749   else if( task == AliQAv1::kDIGITSR ) {
750     fGenDigitsOffset[specie]=offset;
751   }
752   else if( task == AliQAv1::kRECPOINTS ) {
753     fGenRecPointsOffset[specie]=offset;
754   }
755 }
756
757 //_______________________________________________________________
758
759 Int_t AliITSQASDDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task)
760 {
761
762   Int_t histotot=0;
763
764   if( task == AliQAv1::kRAWS )
765     {
766       histotot=fSDDhRawsTask ;  
767     }
768   else if(task == AliQAv1::kDIGITSR)
769     {
770       histotot=fSDDhDigitsTask;
771     }
772   else if( task == AliQAv1::kRECPOINTS )
773     {
774       histotot=fSDDhRecPointsTask;   
775     }
776   else {
777     AliInfo("No task has been selected. TaskHisto set to zero.\n");
778   }
779   return histotot;
780 }
781
782 //_______________________________________________________________
783
784 void AliITSQASDDDataMakerRec::AnalyseBNG()
785 {
786
787 // get file time for Previous test
788         AliInfo("AnalyseBNG\n");
789         Int_t bngtimeBasPrevious; 
790         FILE *filepinPreviousBas = fopen( "SDDbasHistos.time", "r" );
791         if(filepinPreviousBas) {
792           fscanf(filepinPreviousBas,"%d",&bngtimeBasPrevious);
793           fclose (filepinPreviousBas);
794         } else 
795           bngtimeBasPrevious = 0;
796         Int_t bngtimeBasCurrent = 0; 
797         gSystem->Exec("perl -e '@d=localtime ((stat(shift))[9]); printf \"%02d%02d%02d%02d%02d\n\", $d[5]-100,$d[4]+1,$d[3],$d[2],$d[1]'  SDDbasHistos.root > SDDbasHistos.time");
798         FILE *filepinBas = fopen( "SDDbasHistos.time", "r" );
799         fscanf(filepinBas,"%d",&bngtimeBasCurrent);
800         if(bngtimeBasCurrent>bngtimeBasPrevious )AliInfo(Form("bngtimeBasCurrent %d, bngtimeBasPrevious %d\n",bngtimeBasCurrent,bngtimeBasPrevious));
801        
802         Bool_t kAnalyseBas = kTRUE;
803         if(bngtimeBasCurrent <= bngtimeBasPrevious) kAnalyseBas = kFALSE;
804         if(kAnalyseBas) {
805                 // new bng file found
806                 bngtimeBasPrevious = bngtimeBasCurrent;
807                 Bool_t kFilesExist = kTRUE;
808                 TFile basFile("SDDbasHistos.root");
809                 if(basFile.IsZombie()) kFilesExist = kFALSE;
810                 if(kFilesExist) {
811                   AnodeStatus();
812                   AnalyseHistos(1); // Baseline
813                   AnalyseHistos(2); // Uncorrected Noise
814                   AnalyseHistos(3); // Common Mode Noise
815                   AnalyseHistos(4); // Corrected Noise
816                   gSystem->Exec("cp SDDbasHistos.root SDDbasHistosPrevious.root");
817                 } else {
818                         AliInfo("file SDDbasHistos.root not found \n");
819                 }
820         }
821         fclose (filepinBas);
822 //      delete filepinBas;
823
824         Int_t bngtimeGainPrevious; 
825         FILE *filepinPrevious = fopen( "SDDgainHistos.time", "r" );
826         if(filepinPrevious) {
827           fscanf(filepinPrevious,"%d",&bngtimeGainPrevious);
828           fclose (filepinPrevious);
829         } else 
830           bngtimeGainPrevious = 0;
831         Int_t bngtimeGainCurrent = 0; 
832         gSystem->Exec("perl -e '@d=localtime ((stat(shift))[9]); printf \"%02d%02d%02d%02d%02d\n\", $d[5]-100,$d[4]+1,$d[3],$d[2],$d[1]'  SDDgainHistos.root > SDDgainHistos.time");
833         FILE *filepin = fopen( "SDDgainHistos.time", "r" );
834         fscanf(filepin,"%d",&bngtimeGainCurrent);
835         if(bngtimeGainCurrent>bngtimeGainPrevious )AliInfo(Form("bngtimeGainCurrent %d, bngtimeGainPrevious %d\n",bngtimeGainCurrent,bngtimeGainPrevious));
836        
837         Bool_t kAnalyse = kTRUE;
838         if(bngtimeGainCurrent <= bngtimeGainPrevious) kAnalyse = kFALSE;
839         if(kAnalyse) {
840                 // new bng file found
841                 bngtimeGainPrevious = bngtimeGainCurrent;
842                 Bool_t kFilesExist = kTRUE;
843                 TFile gainFile("SDDgainHistos.root");
844                 if(gainFile.IsZombie()) kFilesExist = kFALSE;
845                 if(kFilesExist) {
846                   AnalyseHistos(5); // Gain
847                   gSystem->Exec("cp SDDgainHistos.root SDDgainHistosPrevious.root");
848                 } else {
849                         AliInfo("file SDDgainHistos.root not found \n");
850                 }
851         }
852         fclose (filepin);
853 //      delete filepin;
854
855 }
856
857 //_______________________________________________________________
858
859 void AliITSQASDDDataMakerRec::AnalyseINJ()
860 {
861 // get file time for last test
862
863         AliInfo("AnalyseINJ\n");
864         Int_t injtimePrevious; 
865         FILE *filepinPrevious = fopen( "SDDinjectHistos.time", "r" );
866         if(filepinPrevious) {
867           fscanf(filepinPrevious,"%d",&injtimePrevious);
868           fclose (filepinPrevious);
869         } else 
870           injtimePrevious = 0;
871         Int_t injtimeCurrent = 0; 
872         gSystem->Exec("perl -e '@d=localtime ((stat(shift))[9]); printf \"%02d%02d%02d%02d%02d\n\", $d[5]-100,$d[4]+1,$d[3],$d[2],$d[1]'  SDDinjectHistos.root > SDDinjectHistos.time");
873         FILE *filepin = fopen( "SDDinjectHistos.time", "r" );
874         fscanf(filepin,"%d",&injtimeCurrent);
875         if(injtimeCurrent>injtimePrevious )AliInfo(Form("injtimeCurrent %d, injtimePrevious %d\n",injtimeCurrent,injtimePrevious));
876        
877         Bool_t kAnalyse = kTRUE;
878         if(injtimeCurrent <= injtimePrevious) kAnalyse = kFALSE;
879         if(kAnalyse) {
880                 // new inj file found
881                 injtimePrevious = injtimeCurrent;
882                 Bool_t kFilesExist = kTRUE;
883                 TFile gainFile("SDDinjectHistos.root");
884                 if(gainFile.IsZombie()) kFilesExist = kFALSE;
885                 
886                 if(kFilesExist) {
887                         AnalyseHistos(6); // Drift Speed
888                         gSystem->Exec("cp SDDinjectHistos.root SDDinjectHistosPrevious.root");
889                 } else {
890                         AliInfo("file(s) SDDinjectHistos.root not found \n");
891                 }
892         }
893         fclose (filepin);
894 //      delete filepin;
895 }
896
897 //_______________________________________________________________
898
899 void AliITSQASDDDataMakerRec::AnodeStatus()
900 {
901         char *hnamePrevious = new char[50];
902         fGoodAnodes = 0;
903
904         TFile basFilePrevious("SDDbasHistosPrevious.root");
905         if(!basFilePrevious.IsZombie()) {
906           for(Int_t ddl =0; ddl<fDDLModuleMap->GetNDDLs(); ddl++){
907                 for(Int_t crx =0; crx<fDDLModuleMap->GetNModPerDDL(); crx++){
908                         for(Int_t iside=0;iside<fgknSide;iside++){
909                                 Int_t moduleSDD = fDDLModuleMap->GetModuleNumber(ddl,crx);
910                 //AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
911                                 sprintf(hnamePrevious,"hgood%02dc%02ds%d",ddl,crx,iside);
912                         //AliInfo(Form("get histo %s\n",hnamePrevious));
913                                 TH1F *hgood = (TH1F *) basFilePrevious.Get(hnamePrevious);
914                                 if(!hgood) continue;
915                                 for(Int_t i=0; i<hgood->GetNbinsX();i++) {
916                                         fkAnodeMap[moduleSDD-fgkmodoffset][iside][i] = hgood->GetBinContent(i);
917                                         if(fkAnodeMap[moduleSDD-fgkmodoffset][iside][i]) fGoodAnodes++;
918                                 }
919                                 delete hgood;
920                         }
921                 }
922           }
923           basFilePrevious.Close();
924         }
925         fGoodAnodesCurrent = 0;
926         fBadAnodesCurrent = 0;
927         char *hname = new char[50];
928         Int_t nChangedStatus = 0;
929         Bool_t CurrentAnodeMap[fgknSDDmodules][fgknSide][fgknAnode];    
930         TFile basFile("SDDbasHistos.root");
931         if(!basFile.IsZombie()) {
932           for(Int_t ddl =0; ddl<fDDLModuleMap->GetNDDLs(); ddl++){
933                 for(Int_t crx =0; crx<fDDLModuleMap->GetNModPerDDL(); crx++){
934                         for(Int_t iside=0;iside<fgknSide;iside++){
935                                 Int_t moduleSDD = fDDLModuleMap->GetModuleNumber(ddl,crx);
936                         //AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
937                                 sprintf(hname,"hgood%02dc%02ds%d",ddl,crx,iside);
938                         //AliInfo(Form("get histo %s\n",hname));
939                                 TH1F *hgood = (TH1F *) basFile.Get(hname);
940                                 if(!hgood) continue;
941                                 for(Int_t i=0; i<hgood->GetNbinsX();i++) {
942                                         CurrentAnodeMap[moduleSDD-fgkmodoffset][iside][i] = hgood->GetBinContent(i);
943                                         if(CurrentAnodeMap[moduleSDD-fgkmodoffset][iside][i]) fGoodAnodesCurrent++;
944                                         else fBadAnodesCurrent++;
945                                         if(CurrentAnodeMap[moduleSDD-fgkmodoffset][iside][i] != fkAnodeMap[moduleSDD-fgkmodoffset][iside][i]) {
946                                                 fkAnodeMap[moduleSDD-fgkmodoffset][iside][i] = CurrentAnodeMap[moduleSDD-fgkmodoffset][iside][i];
947                                                 nChangedStatus++;
948                                                 //      AliWarning(Form("DDL %d, CRX %d, Side %d, Anode %d changed status to %d \n",ddl,crx,iside,i,fkAnodeMap[moduleSDD-fgkmodoffset][iside][i]));
949                                         }
950                                 }
951                                 delete hgood;
952                         }
953                 }
954           }
955           basFile.Close();
956         }
957
958         AliWarning(Form("Number of good anodes changed from %d to %d, that is %f %%\n",fGoodAnodes,fGoodAnodesCurrent,((Float_t) TMath::Abs(fGoodAnodes-fGoodAnodesCurrent))/(fBadAnodesCurrent+fGoodAnodesCurrent)));
959         if(fGoodAnodesCurrent != fGoodAnodes) {
960                 fGoodAnodes = fGoodAnodesCurrent;
961         }
962         AliWarning(Form("Number of bad anodes changed from %d to %d, that is %f %%\n",fBadAnodes,fBadAnodesCurrent,((Float_t) TMath::Abs(fBadAnodes-fBadAnodesCurrent))/(fBadAnodesCurrent+fGoodAnodesCurrent)));
963         if(fBadAnodesCurrent != fBadAnodes) {
964                 fBadAnodes = fBadAnodesCurrent;
965         }
966 //      delete hname;
967 }
968
969 //_______________________________________________________________
970
971 void AliITSQASDDDataMakerRec::AnalyseHistos(Int_t type)
972 {
973
974         if(type < 1 || type > 6) {
975           AliWarning(Form("Wrong type (%d), must be between 1 and 6\n",type));
976           return;
977         }
978
979         Double_t Current[fgknSDDmodules][fgknSide][fgknAnode];  
980         char *hnamePrevious = new char[50];
981         TString *filenamePrevious=NULL;
982
983                 if(type < 5) filenamePrevious = new TString("SDDbasHistosPrevious.root");
984                 else if(type == 5) filenamePrevious = new TString("SDDgainHistosPrevious.root");
985                 else if(type == 6) filenamePrevious = new TString("SDDinjectHistosPrevious.root");
986                 TFile *gainFilePrevious = new TFile(filenamePrevious->Data());
987                 if(!gainFilePrevious->IsZombie()) {
988
989                   for(Int_t ddl =0; ddl<fDDLModuleMap->GetNDDLs(); ddl++){
990                     for(Int_t crx =0; crx<fDDLModuleMap->GetNModPerDDL(); crx++){
991                                 for(Int_t iside=0;iside<fgknSide;iside++){
992                                         Int_t moduleSDD = fDDLModuleMap->GetModuleNumber(ddl,crx);
993                         //AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
994                                         if(type == 1) sprintf(hnamePrevious,"hbase%02dc%02ds%d",ddl,crx,iside);
995                                         else if(type == 2) sprintf(hnamePrevious,"hnois%02dc%02ds%d",ddl,crx,iside);
996                                         else if(type == 3) sprintf(hnamePrevious,"hcmn%02dc%02ds%d",ddl,crx,iside);
997                                         else if(type == 4) sprintf(hnamePrevious,"hcorn%02dc%02ds%d",ddl,crx,iside);
998                                         else if(type == 5) sprintf(hnamePrevious,"hgain%02dc%02ds%d",ddl,crx,iside);
999                                         else if(type == 6) sprintf(hnamePrevious,"hdrsp%02dc%02ds%d",ddl,crx,iside);
1000                                 //AliInfo(Form("get histo %s\n",hnamePrevious));
1001                                         TH1F *hhist = (TH1F *) gainFilePrevious->Get(hnamePrevious);
1002                                         if(!hhist) continue;
1003                                         for(Int_t i=0; i<hhist->GetNbinsX();i++) {
1004                                                 Current[moduleSDD-fgkmodoffset][iside][i] = hhist->GetBinContent(i);
1005                                         }
1006                                         delete hhist;
1007                                 }
1008                         }
1009                   }
1010                   gainFilePrevious->Close();
1011                   delete gainFilePrevious;
1012                 }
1013                 delete filenamePrevious;
1014
1015         Float_t xmin = 0.;
1016         Float_t xmax = 0;
1017         Int_t nbins = 1;
1018         TH1F *hDist = 0;
1019         TH1F *hDistDiff = 0;
1020         if(type == 1) {
1021                 xmin = 0.;
1022                 xmax = 500.;
1023                 nbins = (Int_t)(xmax-xmin);
1024                 hDist = new TH1F("hBaseline","Baseline",nbins,xmin,xmax);
1025                 hDistDiff = new TH1F("hBaselineDiff","Baseline Difference",200,-100.,100.);
1026         } else if(type == 2) {
1027                 xmin = 0.;
1028                 xmax = 10.;
1029                 nbins = (Int_t) (100.*(xmax-xmin));
1030                 hDist = new TH1F("hNoiseUnc","Noise (before correction)",nbins,xmin,xmax);
1031                 hDistDiff = new TH1F("hNoiseUncDiff","Noise (before correction) Difference",200,-10.,10.);
1032         } else if(type == 3) {
1033                 xmin = 0.;
1034                 xmax = 10.;
1035                 nbins = (Int_t)( 100.*(xmax-xmin));
1036                 hDist = new TH1F("hNoiseCMN","Noise (common mode)",nbins,xmin,xmax);
1037                 hDistDiff = new TH1F("hNoiseCMNDiff","Noise (common mode) Difference",200,-10.,10.);
1038         } else if(type == 4) {
1039                 xmin = 0.;
1040                 xmax = 10.;
1041                 nbins = (Int_t)(100.*(xmax-xmin));
1042                 hDist = new TH1F("hNoiseCor","Noise (after correction)",nbins,xmin,xmax);
1043                 hDistDiff = new TH1F("hNoiseCorDiff","Noise (after correction) Difference",200,-10.,10.);
1044         } else if(type == 5) {
1045                 xmin = 0.;
1046                 xmax = 5.;
1047                 nbins = (Int_t)(100.*(xmax-xmin));
1048                 hDist = new TH1F("hGain","Gain",nbins,xmin,xmax);
1049                 hDistDiff = new TH1F("hGainDiff","Gain Difference",200,-10.,10.);
1050         } else if(type == 6) {
1051                 xmin = 0.;
1052                 xmax = 10.;
1053                 nbins = (Int_t)(100.*(xmax-xmin));
1054                 hDist = new TH1F("hDriftSpeed","Drift Speed",nbins,xmin,xmax);
1055                 hDistDiff = new TH1F("hDriftSpeedDiff","Drift Speed Difference",200,-10.,10.);
1056         }
1057
1058         Float_t binw = (xmax-xmin)/nbins;
1059
1060         TString *filenamePrevious2=NULL;
1061
1062                 if(type < 5) filenamePrevious2 = new TString("SDDbasHistosPrevious.root");
1063                 else if(type == 5) filenamePrevious2 = new TString("SDDgainHistosPrevious.root");
1064                 else if(type == 6) filenamePrevious2 = new TString("SDDinjectHistosPrevious.root");
1065                 TFile *gainFile = new TFile(filenamePrevious2->Data());
1066                 if(!gainFile->IsZombie()) {
1067
1068                   for(Int_t ddl =0; ddl<fDDLModuleMap->GetNDDLs(); ddl++){
1069                         for(Int_t crx =0; crx<fDDLModuleMap->GetNModPerDDL(); crx++){
1070                                 for(Int_t iside=0;iside<fgknSide;iside++){
1071                                         Int_t moduleSDD = fDDLModuleMap->GetModuleNumber(ddl,crx);
1072                         //AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
1073                                         if(type == 1) sprintf(hnamePrevious,"hbase%02dc%02ds%d",ddl,crx,iside);
1074                                         else if(type == 2) sprintf(hnamePrevious,"hnois%02dc%02ds%d",ddl,crx,iside);
1075                                         else if(type == 3) sprintf(hnamePrevious,"hcmn%02dc%02ds%d",ddl,crx,iside);
1076                                         else if(type == 4) sprintf(hnamePrevious,"hcorn%02dc%02ds%d",ddl,crx,iside);
1077                                         else if(type == 5) sprintf(hnamePrevious,"hgain%02dc%02ds%d",ddl,crx,iside);
1078                                         else if(type == 6) sprintf(hnamePrevious,"hdrsp%02dc%02ds%d",ddl,crx,iside);
1079                                 //AliInfo(Form("get histo %s\n",hname));
1080                                         TH1F *hhist = (TH1F *) gainFile->Get(hnamePrevious);
1081                                         if(!hhist) continue;
1082                                         for(Int_t i=0; i<hhist->GetNbinsX();i++) {
1083                                                 if(!fkAnodeMap[moduleSDD-fgkmodoffset][iside][i]) continue;
1084                                                 hDist->Fill(hhist->GetBinContent(i));
1085                                                 hDistDiff->Fill(hhist->GetBinContent(i)-Current[moduleSDD-fgkmodoffset][iside][i]);
1086                                         }
1087                                         delete hhist;
1088                                 }
1089                         }
1090                   }
1091                   gainFile->Close();
1092                   delete gainFile;
1093                 }
1094                 delete filenamePrevious2;
1095
1096         TF1 ff("ff", "gaus", xmin+0.1, xmax-0.1);
1097         hDist->Fit("ff","NWR");
1098 //      hDist->Fit("gaus","","",xmin+0.1, xmax-0.1);
1099 //      Float_t ChiSquared = (Float_t) ff.GetChisquare();
1100 //      Int_t NDF = ff.GetNumberFitPoints() - ff.GetNpar();
1101         Float_t average = (Float_t) ff.GetParameter(1);
1102         Float_t sigma = (Float_t) ff.GetParameter(2);
1103 //      Float_t mean = hDist->GetMean();
1104 //      Float_t rms = hDist->GetRMS();
1105         Int_t badB = 0;
1106         for(Int_t i=0; i<hDist->GetNbinsX();i++) {
1107 //              if(type < 6) 
1108           if(TMath::Abs(i*binw-average) > 4.*sigma) badB += (Int_t)hDist->GetBinContent(i);
1109 //              else
1110 //                      if(TMath::Abs(i-mean) > 4*rms) badB += hDist->GetBinContent(i);
1111         }
1112         Double_t denomi = hDist->GetEntries();
1113         if(denomi == 0) {
1114           denomi = 1;
1115           badB = 0; 
1116         } 
1117         if(type == 1) {
1118                 AliInfo(Form("Number of anodes with baseline out of 4*sigma from average: %d, %f%%\n",badB,100.*((Float_t) badB)/denomi));
1119         } else if(type == 2) {
1120                 AliInfo(Form("Number of anodes with uncorrected noise out of 4*sigma from average: %d, %f%%\n",badB,100.*((Float_t) badB)/denomi));
1121         } else if(type == 3) {
1122                 AliInfo(Form("Number of anodes with common mode noise out of 4*sigma from average: %d, %f%%\n",badB,100.*((Float_t) badB)/denomi));
1123         } else if(type == 4) {
1124                 AliInfo(Form("Number of anodes with corrected noise out of 4*sigma from average: %d, %f%%\n",badB,100.*((Float_t) badB)/denomi));
1125         } else if(type == 5) {
1126                 AliInfo(Form("Number of anodes with gain out of 4*sigma from average: %d, %f%%\n",badB,100.*((Float_t) badB)/denomi));
1127         } else if(type == 6) {
1128                 Int_t badspeed = (Int_t)hDist->GetBinContent(1);
1129                 AliInfo(Form("Number of anodes with drift speed equal to 0: %d\n",badspeed));
1130                 AliInfo(Form("Number of anodes with drift speed out of 4*sigma from average: %d, %f%%\n",badB-badspeed,100.*((Float_t) (badB-badspeed))/(denomi-badspeed)));
1131         }
1132         
1133         TH1F *hDistHistoryCurrent = NULL;
1134         TH1F *hDistHistoryPrevious = NULL;
1135
1136         TFile *gainHistoryFile=NULL;
1137         if(type < 5) 
1138                 gainHistoryFile = new TFile("SDDbasHistory.root","UPDATE");
1139         else if(type ==5) 
1140                 gainHistoryFile = new TFile("SDDgainHistory.root","UPDATE");
1141         else if(type == 6)
1142                 gainHistoryFile = new TFile("SDDinjectHistory.root","UPDATE");
1143         hDist->Write();
1144         hDistDiff->Write();
1145
1146         //AliInfo("SDDgainHistory.root file opened\n");
1147         if(!gainHistoryFile->IsZombie()) {
1148           if(type == 1) {
1149             hDistHistoryPrevious = (TH1F *) gainHistoryFile->Get("hBaselineHistory"); 
1150           }
1151           else if(type == 2){ 
1152             hDistHistoryPrevious = (TH1F *) gainHistoryFile->Get("hNoiseUncHistory");
1153           }
1154           else if(type == 3) {
1155             hDistHistoryPrevious = (TH1F *) gainHistoryFile->Get("hNoiseCMNHistory");
1156           }
1157           else if(type == 4) {
1158             hDistHistoryPrevious = (TH1F *) gainHistoryFile->Get("hNoiseCorHistory");
1159           }
1160           else if(type == 5) {
1161             hDistHistoryPrevious = (TH1F *) gainHistoryFile->Get("hGainHistory");
1162           }
1163           else if(type == 6){
1164             hDistHistoryPrevious = (TH1F *) gainHistoryFile->Get("hDriftSpeedHistory");
1165           }
1166           Char_t newname[30];
1167           if(hDistHistoryPrevious){
1168             sprintf(newname,"%sPrev",hDistHistoryPrevious->GetName());
1169             hDistHistoryPrevious->SetName(newname);
1170             sprintf(newname,"%sPrev",hDistHistoryPrevious->GetTitle());
1171             hDistHistoryPrevious->SetTitle(newname);
1172           }
1173                 //AliInfo(Form("hDistHistoryPrevious %x\n",hDistHistoryPrevious));
1174
1175
1176                 if(!hDistHistoryPrevious) {
1177                         if(type == 1) hDistHistoryCurrent = new TH1F("hBaselineHistory","Average Baseline History",2,0,2);
1178                         else if(type == 2) hDistHistoryCurrent = new TH1F("hNoiseUncHistory","Average Uncorrected Noise History",2,0,2);
1179                         else if(type == 3) hDistHistoryCurrent = new TH1F("hNoiseCMNHistory","Average Common Mode Noise History",2,0,2);
1180                         else if(type == 4) hDistHistoryCurrent = new TH1F("hNoiseCorHistory","Average Corrected Noise History",2,0,2);
1181                         else if(type == 5) hDistHistoryCurrent = new TH1F("hGainHistory","Average Gain History",2,0,2);
1182                         else if(type == 6) hDistHistoryCurrent = new TH1F("hDriftSpeedHistory","Average Drift Speed History",2,0,2);
1183                         //AliInfo(Form("hDistHistoryCurrent 1 %x\n",hDistHistoryCurrent));
1184 //                      if(type < 6) {
1185                                 hDistHistoryCurrent->SetBinContent(1,average);
1186                                 hDistHistoryCurrent->SetBinError(1,sigma);
1187 /*
1188                         } else {
1189                                 hDistHistoryCurrent->SetBinContent(1,mean);
1190                                 hDistHistoryCurrent->SetBinError(1,rms);
1191                         }
1192 */
1193                 } else {
1194                         if(type == 1) hDistHistoryCurrent = new TH1F("hBaselineHistory","Average Baseline History",hDistHistoryPrevious->GetNbinsX()+1,0,hDistHistoryPrevious->GetNbinsX()+1);
1195                         else if(type == 2) hDistHistoryCurrent = new TH1F("hNoiseUncHistory","Average Uncorrected Noise History",hDistHistoryPrevious->GetNbinsX()+1,0,hDistHistoryPrevious->GetNbinsX()+1);
1196                         else if(type == 3) hDistHistoryCurrent = new TH1F("hNoiseCMNHistory","Average Common Mode Noise History",hDistHistoryPrevious->GetNbinsX()+1,0,hDistHistoryPrevious->GetNbinsX()+1);
1197                         else if(type == 4) hDistHistoryCurrent = new TH1F("hNoiseCorHistory","Average Corrected Noise History",hDistHistoryPrevious->GetNbinsX()+1,0,hDistHistoryPrevious->GetNbinsX()+1);
1198                         else if(type == 5) hDistHistoryCurrent = new TH1F("hGainHistory","Average Gain History",hDistHistoryPrevious->GetNbinsX()+1,0,hDistHistoryPrevious->GetNbinsX()+1);
1199                         else if(type == 6) hDistHistoryCurrent = new TH1F("hDriftSpeedHistory","Average Drift Speed History",hDistHistoryPrevious->GetNbinsX()+1,0,hDistHistoryPrevious->GetNbinsX()+1);
1200                         //AliInfo(Form("hBaselineHistory 2 %x\n",hDistHistory));
1201                         for(Int_t i=0;i<hDistHistoryPrevious->GetNbinsX();i++) {
1202                                 hDistHistoryCurrent->SetBinContent(i,hDistHistoryPrevious->GetBinContent(i));
1203                                 hDistHistoryCurrent->SetBinError(i,hDistHistoryPrevious->GetBinError(i));
1204                         }
1205                         hDistHistoryCurrent->SetBinContent(hDistHistoryPrevious->GetNbinsX(),average);
1206                         hDistHistoryCurrent->SetBinError(hDistHistoryPrevious->GetNbinsX(),sigma);
1207                 }
1208
1209         }
1210         hDistHistoryCurrent->Write();
1211         gainHistoryFile->Close();
1212         delete gainHistoryFile;
1213         //delete hname;
1214         delete hDist;
1215         delete hDistDiff;
1216
1217         //if(hDistHistoryCurrent) delete hDistHistoryCurrent;
1218         //if(hDistHistoryPrevious) delete hDistHistoryPrevious;
1219 }
1220
1221 //_______________________________________________________________
1222
1223
1224 void AliITSQASDDDataMakerRec::CreateTheMap()
1225 {
1226
1227   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
1228   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
1229   if(!ddlMapSDD)
1230     {
1231       AliError("Calibration object retrieval failed! SDD will not be processed");
1232       fDDLModuleMap = NULL;
1233       //return rv;
1234     }
1235   fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
1236   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
1237   ddlMapSDD->SetOwner(kTRUE);
1238   if(!cacheStatus)
1239     {
1240       delete ddlMapSDD;
1241     }
1242   AliInfo("DDL Map Created\n ");
1243 }