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