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