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