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