]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeRec.cxx
ZDC automatic scripts updates (Marco Leoncino) + updates in QA config
[u/mrichter/AliRoot.git] / ITS / AliITSDetTypeRec.cxx
1 /***************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Conributors 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
17
18 ////////////////////////////////////////////////////////////////////////
19 // This class defines the "Standard" reconstruction for the ITS       // 
20 // detector.                                                          //
21 //                                                                    //
22 ////////////////////////////////////////////////////////////////////////
23 #include "TObjArray.h"
24 #include "TTree.h"
25
26 #include "AliCDBManager.h"
27 #include "AliCDBEntry.h"
28 #include "AliITSClusterFinder.h"
29 #include "AliITSClusterFinderV2SPD.h"
30 #include "AliITSClusterFinderV2SDD.h"
31 #include "AliITSClusterFinderSDDfast.h"
32 #include "AliITSClusterFinderV2SSD.h"
33 #include "AliITSDetTypeRec.h"
34 #include "AliITSDDLModuleMapSDD.h"
35 #include "AliITSRecPoint.h"
36 #include "AliITSRecPointContainer.h"
37 #include "AliITSCalibrationSDD.h"
38 #include "AliITSMapSDD.h"
39 #include "AliITSCalibrationSSD.h"
40 #include "AliITSNoiseSSDv2.h"
41 #include "AliITSGainSSDv2.h"
42 #include "AliITSBadChannelsSSDv2.h"
43 #include "AliITSNoiseSSD.h"
44 #include "AliITSGainSSD.h"
45 #include "AliITSBadChannelsSSD.h"
46 #include "AliITSresponseSDD.h"
47 #include "AliITSsegmentationSPD.h"
48 #include "AliITSsegmentationSDD.h"
49 #include "AliITSsegmentationSSD.h"
50 #include "AliLog.h"
51 #include "AliITSRawStreamSPD.h"
52 #include "AliITSTriggerConditions.h"
53 #include "AliITSFOSignalsSPD.h"
54 #include "AliRunLoader.h"
55 #include "AliDataLoader.h"
56 #include "AliITSLoader.h"
57
58
59 class AliITSDriftSpeedArraySDD;
60 class AliITSCorrMapSDD;
61 class AliITSRecoParam;
62
63 const Int_t AliITSDetTypeRec::fgkNdettypes = 3;
64 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSPD =  240;
65 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSDD =  260;
66 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSSD = 1698;
67
68 ClassImp(AliITSDetTypeRec)
69
70 //________________________________________________________________
71 AliITSDetTypeRec::AliITSDetTypeRec(): TObject(),
72 fNMod(0),
73 fITSgeom(0),
74 fReconstruction(0),
75 fSegmentation(0),
76 fCalibration(0),
77 fSSDCalibration(0),
78 fSPDDead(0),
79 fSPDSparseDead(0),
80 fTriggerConditions(0),
81 fDigits(0),
82 fFOSignals(0),
83 fDDLMapSDD(0),
84 fRespSDD(0),
85 fAveGainSDD(0),
86 fRecPoints(0),
87 fNRecPoints(0),
88 fFirstcall(kTRUE),
89 fLoadOnlySPDCalib(0),
90 fFastOrFiredMap(1200){
91     // Standard Constructor
92     // Inputs:
93     //    none.
94     // Outputs:
95     //    none.
96     // Return:
97     //   
98
99   fReconstruction = new TObjArray(fgkNdettypes);
100   fDigits = new TObjArray(fgkNdettypes);
101   for(Int_t i=0; i<3; i++){
102     fkDigClassName[i]=0;
103   }
104   fSSDCalibration=new AliITSCalibrationSSD();
105   fNMod = new Int_t [fgkNdettypes];
106   fNMod[0] = fgkDefaultNModulesSPD;
107   fNMod[1] = fgkDefaultNModulesSDD;
108   fNMod[2] = fgkDefaultNModulesSSD;
109   fNRecPoints = 0;
110   
111   
112 }
113
114 //______________________________________________________________________
115 AliITSDetTypeRec::AliITSDetTypeRec(const AliITSDetTypeRec & rec):TObject(rec),
116 fNMod(rec.fNMod),
117 fITSgeom(rec.fITSgeom),
118 fReconstruction(rec.fReconstruction),
119 fSegmentation(rec.fSegmentation),
120 fCalibration(rec.fCalibration),
121 fSSDCalibration(rec.fSSDCalibration),
122 fSPDDead(rec.fSPDDead),
123 fSPDSparseDead(rec.fSPDSparseDead),
124 fTriggerConditions(rec.fTriggerConditions),
125 fDigits(rec.fDigits),
126 fFOSignals(rec.fFOSignals),
127 fDDLMapSDD(rec.fDDLMapSDD),
128 fRespSDD(rec.fRespSDD),
129 fAveGainSDD(rec.fAveGainSDD),
130 fRecPoints(rec.fRecPoints),
131 fNRecPoints(rec.fNRecPoints),
132 fFirstcall(rec.fFirstcall),
133 fLoadOnlySPDCalib(rec.fLoadOnlySPDCalib),
134 fFastOrFiredMap(rec.fFastOrFiredMap){
135
136   // Copy constructor. 
137   for(Int_t i=0; i<3; i++){    
138     fkDigClassName[i]=rec.fkDigClassName[i];  // NB only copies Char_t*, so not so safe, but this code should never be reached anyways
139   }
140 }
141 //______________________________________________________________________
142 AliITSDetTypeRec& AliITSDetTypeRec::operator=(const AliITSDetTypeRec& source){
143     // Assignment operator. 
144     this->~AliITSDetTypeRec();
145     new(this) AliITSDetTypeRec(source);
146     return *this;
147
148 }
149
150 //_____________________________________________________________________
151 AliITSDetTypeRec::~AliITSDetTypeRec(){
152   //Destructor
153
154   if(fReconstruction){
155     fReconstruction->Delete();
156     delete fReconstruction;
157   }
158   if(fSegmentation){
159     fSegmentation->Delete();
160     delete fSegmentation;
161   }
162   if(fCalibration){
163     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
164       fCalibration->Delete();
165       delete fCalibration;
166       if(fRespSDD) delete fRespSDD;
167       if(fDDLMapSDD) delete fDDLMapSDD;
168    }
169   }
170   if(fSSDCalibration){
171     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
172       delete fSSDCalibration;
173     }
174   }
175    if(fSPDDead){
176     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
177       fSPDDead->Delete();
178       delete fSPDDead;
179     }
180   } 
181      if(fSPDSparseDead){
182     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
183       fSPDSparseDead->Delete();
184       delete fSPDSparseDead;
185     }
186   } 
187   if(fTriggerConditions){
188     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
189       fTriggerConditions->Delete();
190       delete fTriggerConditions;
191     }
192   } 
193   if(fDigits){
194     fDigits->Delete();
195     delete fDigits;
196   }
197   if(fRecPoints){
198     fRecPoints->Delete();
199     delete fRecPoints;
200   }
201   delete [] fNMod;
202   
203   if (fITSgeom) delete fITSgeom;
204  
205 }
206
207 //___________________________________________________________________
208 void AliITSDetTypeRec::SetReconstructionModel(Int_t dettype,AliITSClusterFinder *clf){
209
210   //Set reconstruction model for detector type
211
212   if(fReconstruction==0) fReconstruction = new TObjArray(fgkNdettypes);
213   if(fReconstruction->At(dettype)!=0) delete fReconstruction->At(dettype);
214   fReconstruction->AddAt(clf,dettype);
215 }
216 //______________________________________________________________________
217 AliITSClusterFinder* AliITSDetTypeRec::GetReconstructionModel(Int_t dettype) const{
218
219   //Get reconstruction model for detector type
220   if(fReconstruction==0)  {
221     Warning("GetReconstructionModel","fReconstruction is 0!");
222     return 0;     
223   }
224   return (AliITSClusterFinder*)fReconstruction->At(dettype);
225 }
226
227 //______________________________________________________________________
228 void AliITSDetTypeRec::SetSegmentationModel(Int_t dettype,AliITSsegmentation *seg){
229    
230   //Set segmentation model for detector type
231   
232   if(fSegmentation==0) fSegmentation = new TObjArray(fgkNdettypes);
233   if(fSegmentation->At(dettype)!=0) delete fSegmentation->At(dettype);
234   fSegmentation->AddAt(seg,dettype);
235
236 }
237 //______________________________________________________________________
238 AliITSsegmentation* AliITSDetTypeRec::GetSegmentationModel(Int_t dettype) const {
239
240   //Get segmentation model for detector type
241    
242    if(fSegmentation==0) {
243      Warning("GetSegmentationModel","fSegmentation is 0!");
244      return 0; 
245    } 
246    return (AliITSsegmentation*)fSegmentation->At(dettype);
247
248 }
249 //_______________________________________________________________________
250 void AliITSDetTypeRec::SetCalibrationModel(Int_t iMod, AliITSCalibration *cal){
251
252   //Set calibration (response) for the module iMod of type dettype
253   if (fCalibration==0) {
254     fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
255     fCalibration->SetOwner(kTRUE);
256     fCalibration->Clear();
257   }
258
259   if (fCalibration->At(iMod) != 0)
260     delete (AliITSCalibration*) fCalibration->At(iMod);
261   fCalibration->AddAt(cal,iMod);
262
263 }
264 //_______________________________________________________________________
265 void AliITSDetTypeRec::SetSPDDeadModel(Int_t iMod, AliITSCalibration *cal){
266
267   //Set dead pixel info for the SPD module iMod
268   if (fSPDDead==0) {
269     fSPDDead = new TObjArray(fgkDefaultNModulesSPD);
270     fSPDDead->SetOwner(kTRUE);
271     fSPDDead->Clear();
272   }
273
274   if (fSPDDead->At(iMod) != 0)
275     delete (AliITSCalibration*) fSPDDead->At(iMod);
276   fSPDDead->AddAt(cal,iMod);
277 }
278 //_______________________________________________________________________
279 void AliITSDetTypeRec::SetSPDSparseDeadModel(Int_t iMod, AliITSCalibration *cal){
280
281   //Set dead pixel info for the SPD ACTIVE module iMod
282   if (fSPDSparseDead==0) {
283     fSPDSparseDead = new TObjArray(fgkDefaultNModulesSPD);
284     fSPDSparseDead->SetOwner(kTRUE);
285     fSPDSparseDead->Clear();
286   }
287
288   if (fSPDSparseDead->At(iMod) != 0)
289     delete (AliITSCalibration*) fSPDSparseDead->At(iMod);
290   fSPDSparseDead->AddAt(cal,iMod);
291 }
292 //_______________________________________________________________________
293 AliITSCalibration* AliITSDetTypeRec::GetCalibrationModel(Int_t iMod) const {
294   
295   //Get calibration model for module type
296   
297   if(fCalibration==0) {
298     Warning("GetalibrationModel","fCalibration is 0!");
299     return 0; 
300   }  
301
302   if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
303     return (AliITSCalibration*)fCalibration->At(iMod);
304   }else{
305     Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
306     fSSDCalibration->SetModule(i);
307     return (AliITSCalibration*)fSSDCalibration;
308   }
309
310 }
311 //_______________________________________________________________________
312 AliITSCalibration* AliITSDetTypeRec::GetSPDDeadModel(Int_t iMod) const {
313   
314   //Get SPD dead for module iMod
315   
316   if(fSPDDead==0) {
317     AliWarning("fSPDDead is 0!");
318     return 0; 
319   }  
320   return (AliITSCalibration*)fSPDDead->At(iMod);
321 }
322 //_______________________________________________________________________
323 AliITSCalibration* AliITSDetTypeRec::GetSPDSparseDeadModel(Int_t iMod) const {
324   
325   //Get SPD dead for module iMod
326   
327   if(fSPDSparseDead==0) {
328     AliWarning("fSPDSparseDead is 0!");
329     return 0; 
330   }  
331   return (AliITSCalibration*)fSPDSparseDead->At(iMod);
332 }
333 //_______________________________________________________________________
334 AliITSTriggerConditions* AliITSDetTypeRec::GetTriggerConditions() const {
335   //Get Pixel Trigger Conditions
336   if (fTriggerConditions==0) {
337     AliWarning("fTriggerConditions is 0!");
338   }
339   return fTriggerConditions;
340 }
341 //______________________________________________________________________
342 void AliITSDetTypeRec::SetTreeAddressD(TTree* const treeD){
343     // Set branch address for the tree of digits.
344
345   const char *det[4] = {"SPD","SDD","SSD","ITS"};
346   TBranch *branch;
347   const Char_t* digclass;
348   Int_t i;
349   char branchname[30];
350
351   if(!treeD) return;
352   if (fDigits == 0x0) {
353     fDigits = new TObjArray(fgkNdettypes);
354   }
355   else {
356     ResetDigits();
357   }
358   for (i=0; i<fgkNdettypes; i++) {
359     digclass = GetDigitClassName(i);
360     fDigits->AddAt(new TClonesArray(digclass,1000),i); 
361     if (fgkNdettypes==3) snprintf(branchname,29,"%sDigits%s",det[3],det[i]);
362     else  snprintf(branchname,29,"%sDigits%d",det[3],i+1);
363     branch = treeD->GetBranch(branchname);
364     if (branch) branch->SetAddress(&((*fDigits)[i]));
365   } 
366
367 }
368
369 //_______________________________________________________________________
370 TBranch* AliITSDetTypeRec::MakeBranchInTree(TTree* const tree, 
371                            const char* name, const char *classname, 
372                            void* address,Int_t size,Int_t splitlevel)
373
374 //
375 // Makes branch in given tree and diverts them to a separate file
376 // 
377 //
378 //
379     
380   if (tree == 0x0) {
381     Error("MakeBranchInTree","Making Branch %s Tree is NULL",name);
382     return 0x0;
383   }
384   TBranch *branch = tree->GetBranch(name);
385   if (branch) {  
386     return branch;
387   }
388   if (classname){
389     branch = tree->Branch(name,classname,address,size,splitlevel);
390   }
391   else {
392     branch = tree->Bronch(name, "TClonesArray", address, size, splitlevel);
393   }
394   
395   return branch;
396 }
397
398 //____________________________________________________________________
399 void AliITSDetTypeRec::SetDefaults(){
400   
401   //Set defaults for segmentation and response
402
403   if(!GetITSgeom()){
404     Warning("SetDefaults","null pointer to AliITSgeomGeom !");
405     return;
406   }
407
408   AliITSsegmentation* seg;
409   if(!GetCalibration()) {AliFatal("Exit");exit(0);}  
410
411   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
412     if(dettype==0){
413       seg = new AliITSsegmentationSPD();
414       SetSegmentationModel(dettype,seg);
415       SetDigitClassName(dettype,"AliITSdigitSPD");
416     }
417     if(dettype==1){
418       seg = new AliITSsegmentationSDD();
419       if(fLoadOnlySPDCalib==kFALSE){
420         AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
421         if(cal->IsAMAt20MHz()){ 
422           seg->SetPadSize(seg->Dpz(0),20.);
423           seg->SetNPads(seg->Npz()/2,128);
424         }
425       }
426       SetSegmentationModel(dettype,seg);
427       SetDigitClassName(dettype,"AliITSdigitSDD");
428     }
429     if(dettype==2){
430       AliITSsegmentationSSD* seg2 = new AliITSsegmentationSSD();
431       SetSegmentationModel(dettype,seg2);
432       SetDigitClassName(dettype,"AliITSdigitSSD");
433     }
434   }
435 }
436 //______________________________________________________________________
437 Bool_t AliITSDetTypeRec::GetCalibration() {
438   // Get Default calibration if a storage is not defined.
439
440   if(!fFirstcall){
441     AliITSCalibration* cal = GetCalibrationModel(0);
442     if(cal)return kTRUE;
443   }else {
444     fFirstcall = kFALSE;
445   }
446
447   //  SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
448   //  Int_t run=GetRunNumber();
449
450   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
451   if (fCalibration==0) {
452     fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
453     fCalibration->SetOwner(!cacheStatus);
454     fCalibration->Clear();
455   }
456     
457   Bool_t retCode=GetCalibrationSPD(cacheStatus);
458   if(retCode==kFALSE) return kFALSE;
459
460   if(fLoadOnlySPDCalib==kFALSE){
461     retCode=GetCalibrationSDD(cacheStatus);
462     if(retCode==kFALSE) return kFALSE;
463     retCode=GetCalibrationSSD(cacheStatus);
464     if(retCode==kFALSE) return kFALSE;
465   }
466
467   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
468                fNMod[0], fNMod[1], fNMod[2]));
469   return kTRUE;
470 }
471 //______________________________________________________________________
472 Bool_t AliITSDetTypeRec::GetCalibrationSPD(Bool_t cacheStatus) {
473   // Get SPD calibration objects from OCDB
474   // dead pixel are not used for local reconstruction
475
476  
477   AliCDBEntry *noisySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy");
478   AliCDBEntry *deadSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead");
479   AliCDBEntry *deadSparseSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDSparseDead");
480   AliCDBEntry *pitCond = AliCDBManager::Instance()->Get("TRIGGER/SPD/PITConditions");
481   if(!noisySPD || !deadSPD || !pitCond ){
482     AliFatal("SPD Calibration object retrieval failed! ");
483     return kFALSE;
484   }
485
486   TObjArray *calNoisySPD = (TObjArray*) noisySPD->GetObject();
487   if (!cacheStatus) noisySPD->SetObject(NULL);
488   noisySPD->SetOwner(kTRUE);
489  
490   TObjArray *calDeadSPD = (TObjArray*) deadSPD->GetObject();
491   if (!cacheStatus) deadSPD->SetObject(NULL);
492   deadSPD->SetOwner(kTRUE);
493
494   TObjArray *calSparseDeadSPD = (TObjArray*) deadSparseSPD->GetObject();
495   if (!cacheStatus) deadSparseSPD->SetObject(NULL);
496   deadSparseSPD->SetOwner(kTRUE);
497
498   
499   AliITSTriggerConditions *calPitCond = (AliITSTriggerConditions*) pitCond->GetObject();
500   if (!cacheStatus) pitCond->SetObject(NULL);
501   pitCond->SetOwner(kTRUE);
502
503   if(!cacheStatus){
504     delete noisySPD;
505     delete deadSPD;
506     delete deadSparseSPD;
507     delete pitCond;
508   }
509   if ((!calNoisySPD) || (!calDeadSPD) || (!calSparseDeadSPD) || (!calPitCond)){ 
510     AliWarning("Can not get SPD calibration from calibration database !");
511     return kFALSE;
512   }
513   fNMod[0] = calNoisySPD->GetEntries();
514
515   AliITSCalibration* cal;
516   for (Int_t i=0; i<fNMod[0]; i++) {
517     cal = (AliITSCalibration*) calNoisySPD->At(i);
518     SetCalibrationModel(i, cal);
519     cal = (AliITSCalibration*) calDeadSPD->At(i);
520     SetSPDDeadModel(i, cal);
521     cal = (AliITSCalibration*) calSparseDeadSPD->At(i);
522     SetSPDSparseDeadModel(i, cal);
523   }
524   fTriggerConditions = calPitCond;
525
526   return kTRUE;
527 }
528
529 //______________________________________________________________________
530 Bool_t AliITSDetTypeRec::GetCalibrationSDD(Bool_t cacheStatus) {
531   // Get SDD calibration objects from OCDB
532
533   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
534   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD");
535   AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD");
536   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
537   //   AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD");
538   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD");
539
540   if(!entrySDD || !entry2SDD || !drSpSDD || !ddlMapSDD || !mapTSDD ){
541     AliFatal("SDD Calibration object retrieval failed! ");
542     return kFALSE;
543   }     
544
545
546     
547   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
548   if(!cacheStatus)entrySDD->SetObject(NULL);
549   entrySDD->SetOwner(kTRUE);
550  
551   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
552   if(!cacheStatus)entry2SDD->SetObject(NULL);
553   entry2SDD->SetOwner(kTRUE);
554
555   TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
556   if(!cacheStatus)drSpSDD->SetObject(NULL);
557   drSpSDD->SetOwner(kTRUE);
558
559   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
560   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
561   ddlMapSDD->SetOwner(kTRUE);
562
563 //   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
564 //   if(!cacheStatus)mapASDD->SetObject(NULL);
565 //   mapASDD->SetOwner(kTRUE);
566
567   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
568   if(!cacheStatus)mapTSDD->SetObject(NULL);
569   mapTSDD->SetOwner(kTRUE);
570
571
572   // DB entries are deleted. In this way metadeta objects are deleted as well
573   if(!cacheStatus){
574     delete entrySDD;
575     delete entry2SDD;
576     //delete mapASDD;
577     delete mapTSDD;
578     delete drSpSDD;
579     delete ddlMapSDD;
580   }
581
582   if ((!pSDD)||(!calSDD) || (!drSp) || (!ddlsdd) || (!mapT) ){
583     AliWarning("Can not get SDD calibration from calibration database !");
584     return kFALSE;
585   }
586
587   fNMod[1] = calSDD->GetEntries();
588
589   fDDLMapSDD=ddlsdd;
590   fRespSDD=pSDD;
591   AliITSCalibration* cal;
592   Float_t avegain=0.;
593   Float_t nGdAnodes=0;
594   Bool_t oldMapFormat=kFALSE;
595   TObject* objmap=(TObject*)mapT->At(0);
596   TString cname(objmap->ClassName());
597   if(cname.CompareTo("AliITSMapSDD")==0){ 
598     oldMapFormat=kTRUE;
599     AliInfo("SDD Maps converted to new format");
600   }
601   for(Int_t iddl=0; iddl<AliITSDDLModuleMapSDD::GetNDDLs(); iddl++){
602     for(Int_t icar=0; icar<AliITSDDLModuleMapSDD::GetNModPerDDL();icar++){
603       Int_t iMod=fDDLMapSDD->GetModuleNumber(iddl,icar);
604       if(iMod==-1) continue;
605       Int_t i=iMod - fgkDefaultNModulesSPD;
606       cal = (AliITSCalibration*) calSDD->At(i);
607       Int_t i0=2*i;
608       Int_t i1=1+2*i;
609       for(Int_t iAnode=0;iAnode< ((AliITSCalibrationSDD*)cal)->NOfAnodes(); iAnode++){
610         if(((AliITSCalibrationSDD*)cal)->IsBadChannel(iAnode)) continue;
611         avegain+= ((AliITSCalibrationSDD*)cal)->GetChannelGain(iAnode);
612         nGdAnodes++;
613       }
614       AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
615       AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
616
617       AliITSCorrMapSDD* mt0 = 0;
618       AliITSCorrMapSDD* mt1 = 0;
619       if(oldMapFormat){ 
620         AliITSMapSDD* oldmap0=(AliITSMapSDD*)mapT->At(i0);
621         AliITSMapSDD* oldmap1=(AliITSMapSDD*)mapT->At(i1);
622         mt0=oldmap0->ConvertToNewFormat();
623         mt1=oldmap1->ConvertToNewFormat();
624       }else{
625         mt0=(AliITSCorrMapSDD*)mapT->At(i0);
626         mt1=(AliITSCorrMapSDD*)mapT->At(i1);
627       }
628       cal->SetDriftSpeed(0,arr0);
629       cal->SetDriftSpeed(1,arr1);
630       cal->SetMapT(0,mt0);
631       cal->SetMapT(1,mt1);
632       SetCalibrationModel(iMod, cal);
633     }
634   }
635   if(nGdAnodes) fAveGainSDD=avegain/nGdAnodes;
636   return kTRUE;
637 }
638
639
640 //______________________________________________________________________
641 Bool_t AliITSDetTypeRec::GetCalibrationSSD(Bool_t cacheStatus) {
642   // Get SSD calibration objects from OCDB
643   //  AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD");
644
645   AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
646   AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
647   AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
648
649   if(!entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD){
650     AliFatal("SSD Calibration object retrieval failed! ");
651     return kFALSE;
652   }     
653
654   TObject *emptyssd = 0; TString ssdobjectname;
655   AliITSNoiseSSDv2 *noiseSSD = NULL; 
656   emptyssd = (TObject *)entryNoiseSSD->GetObject();
657   ssdobjectname = emptyssd->GetName();
658   if(ssdobjectname=="TObjArray") {
659     TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
660     noiseSSD = new AliITSNoiseSSDv2(); 
661     ReadOldSSDNoise(noiseSSDOld, noiseSSD);
662   }
663   else if(ssdobjectname=="AliITSNoiseSSDv2")
664     noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
665   if(!cacheStatus)entryNoiseSSD->SetObject(NULL);
666   entryNoiseSSD->SetOwner(kTRUE);
667
668   AliITSGainSSDv2 *gainSSD = NULL;;
669   emptyssd = (TObject *)entryGainSSD->GetObject();
670   ssdobjectname = emptyssd->GetName();
671   if(ssdobjectname=="Gain") {
672     TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
673     gainSSD = new AliITSGainSSDv2();
674     ReadOldSSDGain(gainSSDOld, gainSSD);
675   }
676   else if(ssdobjectname=="AliITSGainSSDv2")
677     gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
678   if(!cacheStatus)entryGainSSD->SetObject(NULL);
679   entryGainSSD->SetOwner(kTRUE);
680
681   AliITSBadChannelsSSDv2 *badChannelsSSD = NULL;
682   emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
683   ssdobjectname = emptyssd->GetName();
684   if(ssdobjectname=="TObjArray") {
685     TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
686     badChannelsSSD = new AliITSBadChannelsSSDv2();
687     ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
688   }
689   else if(ssdobjectname=="AliITSBadChannelsSSDv2")
690     badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
691   if(!cacheStatus)entryBadChannelsSSD->SetObject(NULL);
692   entryBadChannelsSSD->SetOwner(kTRUE);
693
694   // DB entries are deleted. In this way metadeta objects are deleted as well
695   if(!cacheStatus){
696     delete entryNoiseSSD;
697     delete entryGainSSD;
698     delete entryBadChannelsSSD;
699   }
700
701   if ((!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
702     AliWarning("Can not get SSD calibration from calibration database !");
703     return kFALSE;
704   }
705
706   fSSDCalibration->SetNoise(noiseSSD);
707   fSSDCalibration->SetGain(gainSSD);
708   fSSDCalibration->SetBadChannels(badChannelsSSD);
709   //fSSDCalibration->FillBadChipMap();
710
711   return kTRUE;
712 }
713
714 //________________________________________________________________
715 void AliITSDetTypeRec::SetDefaultClusterFindersV2(Bool_t rawdata, Bool_t fastSDD){
716
717   //Set defaults for cluster finder V2
718
719   if(!GetITSgeom()){
720     Warning("SetDefaults","Null pointer to AliITSgeom !");
721     return;
722   }
723
724   AliITSClusterFinder *clf; 
725
726   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
727     //SPD
728     if(dettype==0){
729       if(!GetReconstructionModel(dettype)){
730         clf = new AliITSClusterFinderV2SPD(this);
731         clf->InitGeometry();
732         if(!rawdata) clf->SetDigits(DigitsAddress(0));
733         SetReconstructionModel(dettype,clf);
734
735       }
736     }
737     //SDD
738     if(dettype==1){
739       if(!GetReconstructionModel(dettype)){
740         if(fastSDD){
741           clf = new AliITSClusterFinderSDDfast(this);
742         }
743         else {
744           clf = new AliITSClusterFinderV2SDD(this);
745         }
746         clf->InitGeometry();
747         if(!rawdata) clf->SetDigits(DigitsAddress(1));
748         SetReconstructionModel(dettype,clf);
749       }
750
751     }
752
753     //SSD
754     if(dettype==2){
755       if(!GetReconstructionModel(dettype)){
756         clf = new AliITSClusterFinderV2SSD(this);
757         clf->InitGeometry();
758         if(!rawdata) clf->SetDigits(DigitsAddress(2));
759         SetReconstructionModel(dettype,clf);
760       }
761     }
762
763  }
764    
765 }
766 //______________________________________________________________________
767 void AliITSDetTypeRec::MakeBranch(TTree* tree, Option_t* option){
768
769   //Creates branches for clusters and recpoints
770   Bool_t cR = (strstr(option,"R")!=0);
771   Bool_t cRF = (strstr(option,"RF")!=0);
772   
773   if(cRF)cR = kFALSE;
774
775   if(cR) MakeBranchR(tree);
776   if(cRF) MakeBranchRF(tree);
777
778 }
779
780 //___________________________________________________________________
781 void AliITSDetTypeRec::ResetDigits(){
782   // Reset number of digits and the digits array for the ITS detector.
783   
784   if(!fDigits) return;
785   for(Int_t i=0;i<fgkNdettypes;i++){
786     ResetDigits(i);
787   }
788 }
789 //___________________________________________________________________
790 void AliITSDetTypeRec::ResetDigits(Int_t branch){
791   // Reset number of digits and the digits array for this branch.
792   
793   if(fDigits->At(branch)) ((TClonesArray*)fDigits->At(branch))->Clear();
794
795 }
796 //__________________________________________________________________
797 void AliITSDetTypeRec::MakeBranchR(TTree *treeR, Option_t *opt){
798
799   //Creates tree branches for recpoints
800   // Inputs:
801   //      cont char *file  File name where RecPoints branch is to be written
802   //                       to. If blank it write the SDigits to the same
803   //                       file in which the Hits were found.
804
805   Int_t buffsz = 4000;
806   char branchname[30];
807
808   // only one branch for rec points for all detector types
809   Bool_t oFast= (strstr(opt,"Fast")!=0);
810   
811   Char_t detname[10] = "ITS";
812  
813   
814   if(oFast){
815     snprintf(branchname,29,"%sRecPointsF",detname);
816   } else {
817     snprintf(branchname,29,"%sRecPoints",detname);
818   }
819   
820   if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
821   if (treeR)
822     MakeBranchInTree(treeR,branchname,0,&fRecPoints,buffsz,99);
823 }
824 //______________________________________________________________________
825 void AliITSDetTypeRec::SetTreeAddressR(TTree* const treeR){
826     // Set branch address for the Reconstructed points Trees.
827     // Inputs:
828     //      TTree *treeR   Tree containing the RecPoints.
829     // Outputs:
830     //      none.
831     // Return:
832
833    char branchname[30];
834    Char_t namedet[10]="ITS";
835
836    if(!treeR) return;
837    if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
838    TBranch *branch;
839    snprintf(branchname,29,"%sRecPoints",namedet);
840    branch = treeR->GetBranch(branchname);
841    if (branch) {
842       branch->SetAddress(&fRecPoints);
843     } 
844     else {
845       snprintf(branchname,29,"%sRecPointsF",namedet);
846       branch = treeR->GetBranch(branchname);
847       if (branch) {
848         branch->SetAddress(&fRecPoints);
849       }
850    }
851 }
852 //____________________________________________________________________
853 void AliITSDetTypeRec::AddRecPoint(const AliITSRecPoint &r){
854     // Add a reconstructed space point to the list
855     // Inputs:
856     //      const AliITSRecPoint &r RecPoint class to be added to the tree
857     //                              of reconstructed points TreeR.
858     // Outputs:
859     //      none.
860     // Return:
861     //      none.
862     TClonesArray &lrecp = *fRecPoints;
863     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
864 }
865
866 //______________________________________________________________________
867 void AliITSDetTypeRec::DigitsToRecPoints(TTree *treeD,TTree *treeR,Int_t lastentry,Option_t *opt, Int_t optCluFind){
868   // cluster finding and reconstruction of space points
869   // the condition below will disappear when the geom class will be
870   // initialized for all versions - for the moment it is only for v5 !
871   // 7 is the SDD beam test version
872   // Inputs:
873   //      TTree *treeD     Digits tree
874   //      TTree *treeR     Clusters tree
875   //      Int_t lastentry  Offset for module when not all of the modules
876   //                       are processed.
877   //      Option_t *opt    String indicating which ITS sub-detectors should
878   //                       be processed. If ="All" then all of the ITS
879   //                       sub detectors are processed.
880
881   const char *all = strstr(opt,"All");
882   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
883                         strstr(opt,"SSD")};
884   if(optCluFind==0){
885     SetDefaultClusterFindersV2();
886     AliDebug(1,"V2 cluster finder has been selected \n");
887   }else{
888     SetDefaultClusterFindersV2(kFALSE,kTRUE);
889     AliDebug(1,"SPD and SSD V2 Cluster Finder - SDD fast Cluster Finder \n");    
890   }
891
892   
893   // Reset Fast-OR fired map
894   ResetFastOrFiredMap();
895
896   if (all || det[0]) { // SPD present
897     // Get the FO signals for this event
898     AliRunLoader* runLoader = AliRunLoader::Instance();
899     AliITSLoader* itsLoader = (AliITSLoader*) runLoader->GetLoader("ITSLoader");
900     if (!itsLoader) {
901       AliError("ITS loader is NULL.");
902     }
903    else {
904       fFOSignals = (AliITSFOSignalsSPD*)itsLoader->TreeD()->GetUserInfo()->FindObject("AliITSFOSignalsSPD");
905       if(!fFOSignals) AliError("FO signals not retrieved");
906      }
907
908   }
909
910   
911   AliITSClusterFinder *rec     = 0;
912   Int_t id,module,first=0;
913   for(module=0;module<GetITSgeom()->GetIndexMax();module++){
914       id       = GetITSgeom()->GetModuleType(module);
915       if (!all && !det[id]) continue;
916       if(det[id]) first = GetITSgeom()->GetStartDet(id);
917       rec = (AliITSClusterFinder*)GetReconstructionModel(id);
918       TClonesArray *itsDigits  = DigitsAddress(id);
919       if (!rec){
920         AliFatal("The reconstruction class was not instanciated!");
921         return;
922       }
923       ResetDigits();  // MvL: Not sure we neeed this when rereading anyways
924       if (all) {
925           treeD->GetEvent(lastentry+module);
926         }
927     else {
928       treeD->GetEvent(lastentry+(module-first));
929     }
930     Int_t ndigits = itsDigits->GetEntriesFast();
931     if (ndigits>0 || id==0) { // for SPD we always want to call FindRawClusters (to process FO signals)
932       rec->SetDetTypeRec(this);
933       rec->SetDigits(DigitsAddress(id));
934       //        rec->SetClusters(ClustersAddress(id));
935       rec->FindRawClusters(module);
936     } // end if
937     treeR->Fill();
938     ResetRecPoints();
939   }
940   
941    // Remove PIT in-active chips from Fast-OR fired map
942   if (all || det[0]) { // SPD present
943     RemoveFastOrFiredInActive();
944     // here removing bits which have no associated clusters 
945     RemoveFastOrFiredFromDead(GetFiredChipMap(treeR));  
946   }
947
948   AliITSRecPointContainer* rpcont = AliITSRecPointContainer::Instance();
949   Int_t nClu[6];
950   nClu[0]=rpcont->GetNClustersInLayer(1,treeR);
951   for(Int_t iLay=2; iLay<=6; iLay++) nClu[iLay-1]=rpcont->GetNClustersInLayerFast(iLay);
952   AliInfo(Form("Number of RecPoints in ITS Layers = %d %d %d %d %d %d",
953                nClu[0],nClu[1],nClu[2],nClu[3],nClu[4],nClu[5]));
954 }
955 //______________________________________________________________________
956 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader,TTree *treeR,Option_t *opt){
957   // cluster finding and reconstruction of space points
958   // the condition below will disappear when the geom class will be
959   // initialized for all versions - for the moment it is only for v5 !
960   // 7 is the SDD beam test version
961   // Inputs:
962   //      AliRawReader *rawReader  Pointer to the raw-data reader
963   //      TTree *treeR             Clusters tree
964   // Outputs:
965   //      none.
966   // Return:
967   //      none.
968   const char *all = strstr(opt,"All");
969   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
970                         strstr(opt,"SSD")};
971   
972   Int_t id=0;
973   AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
974   rpc->FullReset();
975   TClonesArray* array = rpc->UncheckedGetClusters(0);
976   TBranch *branch = treeR->Branch("ITSRecPoints",&array);
977   DigitsToRecPoints(rawReader,opt); 
978
979   Int_t nClusters =0;
980   for(Int_t iModule=0;iModule<GetITSgeom()->GetIndexMax();iModule++){
981     id = GetITSgeom()->GetModuleType(iModule);
982     if (!all && !det[id]) continue;
983     array = rpc->UncheckedGetClusters(iModule);
984     if(!array){
985       AliDebug(1,Form("data for module %d missing!",iModule));
986     }
987     branch->SetAddress(&array);
988     treeR->Fill();
989     nClusters+=array->GetEntriesFast();
990   }
991
992   rpc->FullReset();
993
994   AliITSRecPointContainer* rpcont = AliITSRecPointContainer::Instance();
995   Int_t nClu[6];
996   nClu[0]=rpcont->GetNClustersInLayer(1,treeR);
997   for(Int_t iLay=2; iLay<=6; iLay++) nClu[iLay-1]=rpcont->GetNClustersInLayerFast(iLay);
998   AliInfo(Form("Number of RecPoints in ITS Layers = %d %d %d %d %d %d, Total = %d",
999                nClu[0],nClu[1],nClu[2],nClu[3],nClu[4],nClu[5],nClusters));
1000 }
1001 //______________________________________________________________________
1002 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader,Option_t *opt){
1003   // cluster finding and reconstruction of space points
1004   // the condition below will disappear when the geom class will be
1005   // initialized for all versions - for the moment it is only for v5 !
1006   // 7 is the SDD beam test version
1007   // Inputs:
1008   //      AliRawReader *rawReader  Pointer to the raw-data reader
1009   // Outputs:
1010   //      none.
1011   // Return:
1012   //      none.
1013   const char *all = strstr(opt,"All");
1014   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1015                         strstr(opt,"SSD")};
1016   
1017   // Reset Fast-OR fired map
1018   ResetFastOrFiredMap();
1019   
1020   AliITSClusterFinder *rec     = 0;
1021   Int_t id=0;
1022
1023   for(id=0;id<3;id++){
1024     if (!all && !det[id]) continue;
1025     rec = (AliITSClusterFinder*)GetReconstructionModel(id);
1026     if (!rec){
1027       AliFatal("The reconstruction class was not instantiated");
1028       return;
1029     }
1030     rec->SetDetTypeRec(this);
1031     rec->RawdataToClusters(rawReader);    
1032   } 
1033    
1034   // Remove PIT in-active chips from Fast-OR fired map
1035   if (all || det[0]) { // SPD present
1036     RemoveFastOrFiredInActive();
1037     // here removing bits which have no associated clusters 
1038     RemoveFastOrFiredFromDead(GetFiredChipMap());
1039    
1040   }  
1041 }
1042 //______________________________________________________________________
1043 void AliITSDetTypeRec::ReadOldSSDNoise(const TObjArray *array, 
1044                                        AliITSNoiseSSDv2 *noiseSSD) {
1045   //Reads the old SSD calibration object and converts it to the new format
1046   const Int_t fgkSSDSTRIPSPERMODULE = 1536;
1047   const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
1048
1049   Int_t gNMod = array->GetEntries();
1050   AliInfo("Converting old calibration object for noise...\n");
1051
1052   //NOISE
1053   Double_t noise = 0.0;
1054   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1055     AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
1056     for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
1057       noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
1058       if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
1059         noiseSSD->AddNoiseP(iModule,iStrip,noise);
1060       if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
1061         noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
1062     }//loop over strips
1063   }//loop over modules      
1064 }
1065
1066 //______________________________________________________________________
1067 void AliITSDetTypeRec::ReadOldSSDBadChannels(const TObjArray *array, 
1068                                              AliITSBadChannelsSSDv2 *badChannelsSSD) {
1069   //Reads the old SSD calibration object and converts it to the new format
1070   Int_t gNMod = array->GetEntries();
1071   AliInfo("Converting old calibration object for bad channels...");
1072   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1073     //for (Int_t iModule = 0; iModule < 1; iModule++) {
1074     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
1075     TArrayI arrayPSide = bad->GetBadPChannelsList();
1076     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++) 
1077       badChannelsSSD->AddBadChannelP(iModule,
1078                                      iPCounter,
1079                                      (Char_t)arrayPSide.At(iPCounter));
1080         
1081     TArrayI arrayNSide = bad->GetBadNChannelsList();
1082     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++) 
1083       badChannelsSSD->AddBadChannelN(iModule,
1084                                      iNCounter,
1085                                      (Char_t)arrayNSide.At(iNCounter));
1086     
1087   }//loop over modules      
1088 }
1089
1090 //______________________________________________________________________
1091 void AliITSDetTypeRec::ReadOldSSDGain(const TObjArray *array, 
1092                                       AliITSGainSSDv2 *gainSSD) {
1093   //Reads the old SSD calibration object and converts it to the new format
1094
1095   Int_t gNMod = array->GetEntries();
1096   AliInfo("Converting old calibration object for gain...\n");
1097
1098   //GAIN
1099   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1100     AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
1101     TArrayF arrayPSide = gainModule->GetGainP();
1102     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
1103       gainSSD->AddGainP(iModule,
1104                         iPCounter,
1105                         arrayPSide.At(iPCounter));
1106     TArrayF arrayNSide = gainModule->GetGainN();
1107     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
1108       gainSSD->AddGainN(iModule,
1109                         iNCounter,
1110                         arrayNSide.At(iNCounter));
1111   }//loop over modules 
1112 }
1113 //______________________________________________________________________
1114 void AliITSDetTypeRec::RemoveFastOrFiredInActive() {
1115   // Removes the chips that were in-active in the pixel trigger (from fast-or fired map)
1116
1117   if (fTriggerConditions==NULL) {
1118     AliError("Pixel trigger conditions are missing.");
1119     return;
1120   }
1121   Int_t eq   = -1;
1122   Int_t hs   = -1;
1123   Int_t chip = -1;
1124   while (fTriggerConditions->GetNextInActiveChip(eq,hs,chip)) {
1125     UInt_t chipKey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1126     fFastOrFiredMap.SetBitNumber(chipKey,kFALSE);
1127   }
1128 }
1129 //______________________________________________________________________
1130 TBits AliITSDetTypeRec::GetFiredChipMap() const {
1131   
1132   //
1133   // TBits of the fired chips  
1134   //
1135  
1136   AliITSRecPointContainer* rpc = AliITSRecPointContainer::Instance();
1137
1138   TBits isfiredchip(1200);
1139   
1140    AliITSsegmentationSPD *segSPD = (AliITSsegmentationSPD*)GetSegmentationModel(0);
1141    if(!segSPD) {
1142     AliError("no segmentation model for SPD available, the fired chip map is empty. Exiting"); 
1143     return isfiredchip;
1144    }
1145    
1146   
1147   for(Int_t imod =0; imod < fgkDefaultNModulesSPD; imod++){
1148     TClonesArray *array = rpc->UncheckedGetClusters(imod);
1149     if(!array) continue;
1150     Int_t nCluster = array->GetEntriesFast();
1151  
1152     while(nCluster--) {
1153       AliITSRecPoint* cluster = (AliITSRecPoint*)array->UncheckedAt(nCluster);
1154      if (cluster->GetLayer()>1)continue;
1155      Float_t local[3]={-1,-1};
1156      local[1]=cluster->GetDetLocalX();
1157      local[0]=cluster->GetDetLocalZ();
1158      
1159      Int_t eq = AliITSRawStreamSPD::GetOnlineEqIdFromOffline(imod);
1160      Int_t hs = AliITSRawStreamSPD::GetOnlineHSFromOffline(imod);
1161      Int_t row, col;
1162      segSPD->LocalToDet(0.5,local[0],row,col);
1163      Int_t chip = AliITSRawStreamSPD::GetOnlineChipFromOffline(imod,col);
1164      Int_t chipkey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1165      isfiredchip.SetBitNumber(chipkey,kTRUE);
1166     }
1167     
1168   } 
1169  
1170   return isfiredchip;
1171   
1172 }
1173 //______________________________________________________________________
1174 TBits AliITSDetTypeRec::GetFiredChipMap(TTree *treeR) const{
1175   //
1176   // TBits of the fired chips  
1177   //
1178   TBits isfiredchip(1200);
1179   
1180   if(!treeR) {
1181      AliError("no treeR. fired chip map stays empty. Exiting.");
1182      return isfiredchip;
1183    }
1184    
1185   AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
1186   TClonesArray *recpoints = NULL;
1187   rpcont->FetchClusters(0,treeR);
1188   if(!rpcont->GetStatusOK() || !rpcont->IsSPDActive()){
1189     AliError("no clusters. fired chip map stays empty. Exiting.");
1190      return isfiredchip;
1191   }
1192   
1193    AliITSsegmentationSPD *segSPD = (AliITSsegmentationSPD*)GetSegmentationModel(0);
1194       
1195    for(Int_t imod =0; imod < fgkDefaultNModulesSPD; imod++){
1196     recpoints = rpcont->UncheckedGetClusters(imod);
1197     Int_t nCluster = recpoints->GetEntriesFast();
1198     
1199     // loop over clusters
1200     while(nCluster--) {
1201       AliITSRecPoint* cluster = (AliITSRecPoint*)recpoints->UncheckedAt(nCluster);
1202       if (cluster->GetLayer()>1)continue;
1203       Float_t local[3]={-1,-1};
1204       local[1]=cluster->GetDetLocalX();
1205       local[0]=cluster->GetDetLocalZ();
1206       
1207       Int_t eq = AliITSRawStreamSPD::GetOnlineEqIdFromOffline(imod);
1208       Int_t hs = AliITSRawStreamSPD::GetOnlineHSFromOffline(imod);
1209       Int_t row, col;
1210       segSPD->LocalToDet(0.5,local[0],row,col);
1211       Int_t chip = AliITSRawStreamSPD::GetOnlineChipFromOffline(imod,col);
1212       Int_t chipkey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1213       isfiredchip.SetBitNumber(chipkey,kTRUE);
1214     }
1215   }
1216  
1217   return isfiredchip;
1218 }
1219 //______________________________________________________________________
1220 void  AliITSDetTypeRec::RemoveFastOrFiredFromDead(TBits firedchipmap){
1221   //
1222   // resetting of the fast-or bit on cluster basis. 
1223   // fast-or bits can be remnant from SPD ideal simulation (no dead channels)
1224   //
1225   
1226   for(Int_t chipKey=0; chipKey<1200; chipKey++){
1227     // FO masked chips have been previously removed  
1228    if(!fFastOrFiredMap.TestBitNumber(chipKey)) continue; 
1229    if(!firedchipmap.TestBitNumber(chipKey))  {
1230     fFastOrFiredMap.SetBitNumber(chipKey,kFALSE);
1231     AliDebug(2,Form("removing bit in key %i \n ",chipKey));
1232   }
1233  }
1234    
1235 }
1236 //______________________________________________________________________
1237 void AliITSDetTypeRec::SetFastOrFiredMapOnline(UInt_t eq, UInt_t hs, UInt_t chip) {
1238   // Set fast-or fired map for this chip
1239   Int_t chipKey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1240   return SetFastOrFiredMap(chipKey);
1241 }
1242