]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeRec.cxx
Most of these updates are related to the offline software needed for the pixel trigge...
[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 "AliITSCalibrationSDD.h"
38 #include "AliITSHLTforSDD.h"
39 #include "AliITSCalibrationSSD.h"
40 #include "AliITSNoiseSSDv2.h"
41 #include "AliITSGainSSDv2.h"
42 #include "AliITSBadChannelsSSDv2.h"
43 #include "AliITSNoiseSSD.h"
44 #include "AliITSGainSSD.h"
45 #include "AliITSBadChannelsSSD.h"
46 #include "AliITSresponseSDD.h"
47 #include "AliITSsegmentationSPD.h"
48 #include "AliITSsegmentationSDD.h"
49 #include "AliITSsegmentationSSD.h"
50 #include "AliLog.h"
51 #include "AliITSRawStreamSPD.h"
52 #include "AliITSTriggerConditions.h"
53 #include "AliITSFOSignalsSPD.h"
54 #include "AliRunLoader.h"
55 #include "AliDataLoader.h"
56 #include "AliITSLoader.h"
57
58 class AliITSDriftSpeedArraySDD;
59 class AliITSMapSDD;
60 class AliITSRecoParam;
61
62 const Int_t AliITSDetTypeRec::fgkNdettypes = 3;
63 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSPD =  240;
64 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSDD =  260;
65 const Int_t AliITSDetTypeRec::fgkDefaultNModulesSSD = 1698;
66
67 ClassImp(AliITSDetTypeRec)
68
69 //________________________________________________________________
70 AliITSDetTypeRec::AliITSDetTypeRec(): TObject(),
71 fNMod(0),
72 fITSgeom(0),
73 fReconstruction(0),
74 fSegmentation(0),
75 fCalibration(0),
76 fSSDCalibration(0),
77 fSPDDead(0),
78 fTriggerConditions(0),
79 fDigits(0),
80 fFOSignals(0),
81 fDDLMapSDD(0),
82 fRespSDD(0),
83 fAveGainSDD(0),
84 fIsHLTmodeC(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 fIsHLTmodeC(rec.fIsHLTmodeC),
130 fRecPoints(rec.fRecPoints),
131 fNRecPoints(rec.fNRecPoints),
132 fFirstcall(rec.fFirstcall),
133 fLoadOnlySPDCalib(rec.fLoadOnlySPDCalib),
134 fFastOrFiredMap(rec.fFastOrFiredMap){
135
136   // Copy constructor. 
137
138 }
139 //______________________________________________________________________
140 AliITSDetTypeRec& AliITSDetTypeRec::operator=(const AliITSDetTypeRec& source){
141     // Assignment operator. 
142     this->~AliITSDetTypeRec();
143     new(this) AliITSDetTypeRec(source);
144     return *this;
145
146 }
147
148 //_____________________________________________________________________
149 AliITSDetTypeRec::~AliITSDetTypeRec(){
150   //Destructor
151
152   if(fReconstruction){
153     fReconstruction->Delete();
154     delete fReconstruction;
155     fReconstruction = 0;
156   }
157   if(fSegmentation){
158     fSegmentation->Delete();
159     delete fSegmentation;
160     fSegmentation = 0;
161   }
162   if(fCalibration){
163     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
164       fCalibration->Delete();
165       delete fCalibration;
166       fCalibration = 0;
167       if(fRespSDD) delete fRespSDD;
168       if(fDDLMapSDD) delete fDDLMapSDD;
169    }
170   }
171   if(fSSDCalibration) delete fSSDCalibration;
172    if(fSPDDead){
173     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
174       fSPDDead->Delete();
175       delete fSPDDead;
176       fSPDDead = 0;
177     }
178   } 
179   if(fTriggerConditions){
180     if(!(AliCDBManager::Instance()->GetCacheFlag())) {
181       fTriggerConditions->Delete();
182       delete fTriggerConditions;
183       fTriggerConditions = 0;
184     }
185   } 
186   if(fDigits){
187     fDigits->Delete();
188     delete fDigits;
189     fDigits=0;
190   }
191   if(fRecPoints){
192     fRecPoints->Delete();
193     delete fRecPoints;
194     fRecPoints=0;
195   }
196   delete [] fNMod;
197   
198   if (fITSgeom) delete fITSgeom;
199  
200 }
201
202 //___________________________________________________________________
203 void AliITSDetTypeRec::SetReconstructionModel(Int_t dettype,AliITSClusterFinder *clf){
204
205   //Set reconstruction model for detector type
206
207   if(fReconstruction==0) fReconstruction = new TObjArray(fgkNdettypes);
208   if(fReconstruction->At(dettype)!=0) delete fReconstruction->At(dettype);
209   fReconstruction->AddAt(clf,dettype);
210 }
211 //______________________________________________________________________
212 AliITSClusterFinder* AliITSDetTypeRec::GetReconstructionModel(Int_t dettype) const{
213
214   //Get reconstruction model for detector type
215   if(fReconstruction==0)  {
216     Warning("GetReconstructionModel","fReconstruction is 0!");
217     return 0;     
218   }
219   return (AliITSClusterFinder*)fReconstruction->At(dettype);
220 }
221
222 //______________________________________________________________________
223 void AliITSDetTypeRec::SetSegmentationModel(Int_t dettype,AliITSsegmentation *seg){
224    
225   //Set segmentation model for detector type
226   
227   if(fSegmentation==0) fSegmentation = new TObjArray(fgkNdettypes);
228   if(fSegmentation->At(dettype)!=0) delete fSegmentation->At(dettype);
229   fSegmentation->AddAt(seg,dettype);
230
231 }
232 //______________________________________________________________________
233 AliITSsegmentation* AliITSDetTypeRec::GetSegmentationModel(Int_t dettype) const {
234
235   //Get segmentation model for detector type
236    
237    if(fSegmentation==0) {
238      Warning("GetSegmentationModel","fSegmentation is 0!");
239      return 0; 
240    } 
241    return (AliITSsegmentation*)fSegmentation->At(dettype);
242
243 }
244 //_______________________________________________________________________
245 void AliITSDetTypeRec::SetCalibrationModel(Int_t iMod, AliITSCalibration *cal){
246
247   //Set calibration (response) for the module iMod of type dettype
248   if (fCalibration==0) {
249     fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
250     fCalibration->SetOwner(kTRUE);
251     fCalibration->Clear();
252   }
253
254   if (fCalibration->At(iMod) != 0)
255     delete (AliITSCalibration*) fCalibration->At(iMod);
256   fCalibration->AddAt(cal,iMod);
257
258 }
259 //_______________________________________________________________________
260 void AliITSDetTypeRec::SetSPDDeadModel(Int_t iMod, AliITSCalibration *cal){
261
262   //Set dead pixel info for the SPD module iMod
263   if (fSPDDead==0) {
264     fSPDDead = new TObjArray(fgkDefaultNModulesSPD);
265     fSPDDead->SetOwner(kTRUE);
266     fSPDDead->Clear();
267   }
268
269   if (fSPDDead->At(iMod) != 0)
270     delete (AliITSCalibration*) fSPDDead->At(iMod);
271   fSPDDead->AddAt(cal,iMod);
272 }
273 //_______________________________________________________________________
274 AliITSCalibration* AliITSDetTypeRec::GetCalibrationModel(Int_t iMod) const {
275   
276   //Get calibration model for module type
277   
278   if(fCalibration==0) {
279     Warning("GetalibrationModel","fCalibration is 0!");
280     return 0; 
281   }  
282
283   if(iMod<fgkDefaultNModulesSPD+fgkDefaultNModulesSDD){
284     return (AliITSCalibration*)fCalibration->At(iMod);
285   }else{
286     Int_t i=iMod-(fgkDefaultNModulesSPD+fgkDefaultNModulesSDD);
287     fSSDCalibration->SetModule(i);
288     return (AliITSCalibration*)fSSDCalibration;
289   }
290
291 }
292 //_______________________________________________________________________
293 AliITSCalibration* AliITSDetTypeRec::GetSPDDeadModel(Int_t iMod) const {
294   
295   //Get SPD dead for module iMod
296   
297   if(fSPDDead==0) {
298     AliWarning("fSPDDead is 0!");
299     return 0; 
300   }  
301   return (AliITSCalibration*)fSPDDead->At(iMod);
302 }
303 //_______________________________________________________________________
304 AliITSTriggerConditions* AliITSDetTypeRec::GetTriggerConditions() const {
305   //Get Pixel Trigger Conditions
306   if (fTriggerConditions==0) {
307     AliWarning("fTriggerConditions is 0!");
308   }
309   return fTriggerConditions;
310 }
311 //______________________________________________________________________
312 void AliITSDetTypeRec::SetTreeAddressD(TTree* const treeD){
313     // Set branch address for the tree of digits.
314
315   const char *det[4] = {"SPD","SDD","SSD","ITS"};
316   TBranch *branch;
317   const Char_t* digclass;
318   Int_t i;
319   char branchname[30];
320
321   if(!treeD) return;
322   if (fDigits == 0x0) fDigits = new TObjArray(fgkNdettypes);
323   for (i=0; i<fgkNdettypes; i++) {
324     digclass = GetDigitClassName(i);
325     if(!(fDigits->At(i))) {
326       fDigits->AddAt(new TClonesArray(digclass,1000),i);
327     }
328     else{
329       ResetDigits(i);
330     } 
331     if (fgkNdettypes==3) sprintf(branchname,"%sDigits%s",det[3],det[i]);
332     else  sprintf(branchname,"%sDigits%d",det[3],i+1);
333     if (fDigits) {
334       branch = treeD->GetBranch(branchname);
335       if (branch) branch->SetAddress(&((*fDigits)[i]));
336     }
337   } 
338
339 }
340
341 //_______________________________________________________________________
342 TBranch* AliITSDetTypeRec::MakeBranchInTree(TTree* const tree, 
343                            const char* name, const char *classname, 
344                            void* address,Int_t size,Int_t splitlevel)
345
346 //
347 // Makes branch in given tree and diverts them to a separate file
348 // 
349 //
350 //
351     
352   if (tree == 0x0) {
353     Error("MakeBranchInTree","Making Branch %s Tree is NULL",name);
354     return 0x0;
355   }
356   TBranch *branch = tree->GetBranch(name);
357   if (branch) {  
358     return branch;
359   }
360   if (classname){
361     branch = tree->Branch(name,classname,address,size,splitlevel);
362   }
363   else {
364     branch = tree->Bronch(name, "TClonesArray", address, size, splitlevel);
365   }
366   
367   return branch;
368 }
369
370 //____________________________________________________________________
371 void AliITSDetTypeRec::SetDefaults(){
372   
373   //Set defaults for segmentation and response
374
375   if(!GetITSgeom()){
376     Warning("SetDefaults","null pointer to AliITSgeomGeom !");
377     return;
378   }
379
380   AliITSsegmentation* seg;
381   if(!GetCalibration()) {AliFatal("Exit");exit(0);}  
382
383   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
384     if(dettype==0){
385       seg = new AliITSsegmentationSPD();
386       SetSegmentationModel(dettype,seg);
387       SetDigitClassName(dettype,"AliITSdigitSPD");
388     }
389     if(fLoadOnlySPDCalib==kFALSE){
390       if(dettype==1){
391         seg = new AliITSsegmentationSDD();
392         AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)GetCalibrationModel(fgkDefaultNModulesSPD+1);
393         if(cal->IsAMAt20MHz()){ 
394           seg->SetPadSize(seg->Dpz(0),20.);
395           seg->SetNPads(seg->Npz()/2,128);
396         }
397         SetSegmentationModel(dettype,seg);
398         SetDigitClassName(dettype,"AliITSdigitSDD");
399       }
400     }
401     if(dettype==2){
402       AliITSsegmentationSSD* seg2 = new AliITSsegmentationSSD();
403       SetSegmentationModel(dettype,seg2);
404       SetDigitClassName(dettype,"AliITSdigitSSD");
405     }
406   }
407 }
408 //______________________________________________________________________
409 Bool_t AliITSDetTypeRec::GetCalibration() {
410   // Get Default calibration if a storage is not defined.
411
412   if(!fFirstcall){
413     AliITSCalibration* cal = GetCalibrationModel(0);
414     if(cal)return kTRUE;
415   }else {
416     fFirstcall = kFALSE;
417   }
418
419   //  SetRunNumber((Int_t)AliCDBManager::Instance()->GetRun());
420   //  Int_t run=GetRunNumber();
421
422   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
423   if (fCalibration==0) {
424     fCalibration = new TObjArray(GetITSgeom()->GetIndexMax());
425     fCalibration->SetOwner(!cacheStatus);
426     fCalibration->Clear();
427   }
428     
429   Bool_t retCode=GetCalibrationSPD(cacheStatus);
430   if(retCode==kFALSE) return kFALSE;
431
432   if(fLoadOnlySPDCalib==kFALSE){
433     retCode=GetCalibrationSDD(cacheStatus);
434     if(retCode==kFALSE) return kFALSE;
435     retCode=GetCalibrationSSD(cacheStatus);
436     if(retCode==kFALSE) return kFALSE;
437   }
438
439   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
440                fNMod[0], fNMod[1], fNMod[2]));
441   return kTRUE;
442 }
443 //______________________________________________________________________
444 Bool_t AliITSDetTypeRec::GetCalibrationSPD(Bool_t cacheStatus) {
445   // Get SPD calibration objects from OCDB
446   // dead pixel are not used for local reconstruction
447
448  
449   AliCDBEntry *noisySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy");
450   AliCDBEntry *deadSPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDDead");
451   AliCDBEntry *pitCond = AliCDBManager::Instance()->Get("ITS/Calib/PITConditions");
452   if(!noisySPD || !deadSPD || !pitCond ){
453     AliFatal("SPD Calibration object retrieval failed! ");
454     return kFALSE;
455   }
456
457   TObjArray *calNoisySPD = (TObjArray*) noisySPD->GetObject();
458   if (!cacheStatus) noisySPD->SetObject(NULL);
459   noisySPD->SetOwner(kTRUE);
460  
461   TObjArray *calDeadSPD = (TObjArray*) deadSPD->GetObject();
462   if (!cacheStatus) deadSPD->SetObject(NULL);
463   deadSPD->SetOwner(kTRUE);
464
465   AliITSTriggerConditions *calPitCond = (AliITSTriggerConditions*) pitCond->GetObject();
466   if (!cacheStatus) pitCond->SetObject(NULL);
467   pitCond->SetOwner(kTRUE);
468
469   if(!cacheStatus){
470     delete noisySPD;
471     delete deadSPD;
472     delete pitCond;
473   }
474   if ((!calNoisySPD) || (!calDeadSPD) || (!calPitCond)){ 
475     AliWarning("Can not get SPD calibration from calibration database !");
476     return kFALSE;
477   }
478
479   fNMod[0] = calNoisySPD->GetEntries();
480
481   AliITSCalibration* cal;
482   for (Int_t i=0; i<fNMod[0]; i++) {
483     cal = (AliITSCalibration*) calNoisySPD->At(i);
484     SetCalibrationModel(i, cal);
485     cal = (AliITSCalibration*) calDeadSPD->At(i);
486     SetSPDDeadModel(i, cal);
487   }
488   fTriggerConditions = calPitCond;
489
490   return kTRUE;
491 }
492
493 //______________________________________________________________________
494 Bool_t AliITSDetTypeRec::GetCalibrationSDD(Bool_t cacheStatus) {
495   // Get SDD calibration objects from OCDB
496
497   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
498   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD");
499   AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD");
500   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
501   AliCDBEntry *hltforSDD = AliCDBManager::Instance()->Get("ITS/Calib/HLTforSDD");
502   //   AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD");
503   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD");
504
505   if(!entrySDD || !entry2SDD || !drSpSDD || !ddlMapSDD || !hltforSDD || !mapTSDD ){
506     AliFatal("SDD Calibration object retrieval failed! ");
507     return kFALSE;
508   }     
509
510
511     
512   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
513   if(!cacheStatus)entrySDD->SetObject(NULL);
514   entrySDD->SetOwner(kTRUE);
515  
516   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
517   if(!cacheStatus)entry2SDD->SetObject(NULL);
518   entry2SDD->SetOwner(kTRUE);
519
520   TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
521   if(!cacheStatus)drSpSDD->SetObject(NULL);
522   drSpSDD->SetOwner(kTRUE);
523
524   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
525   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
526   ddlMapSDD->SetOwner(kTRUE);
527
528   AliITSHLTforSDD* hltsdd=(AliITSHLTforSDD*)hltforSDD->GetObject();
529   if(!cacheStatus)hltforSDD->SetObject(NULL);
530   hltforSDD->SetOwner(kTRUE);
531   
532 //   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
533 //   if(!cacheStatus)mapASDD->SetObject(NULL);
534 //   mapASDD->SetOwner(kTRUE);
535
536   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
537   if(!cacheStatus)mapTSDD->SetObject(NULL);
538   mapTSDD->SetOwner(kTRUE);
539
540
541   // DB entries are deleted. In this way metadeta objects are deleted as well
542   if(!cacheStatus){
543     delete entrySDD;
544     delete entry2SDD;
545     //delete mapASDD;
546     delete hltforSDD;
547     delete mapTSDD;
548     delete drSpSDD;
549     delete ddlMapSDD;
550   }
551
552   if ((!pSDD)||(!calSDD) || (!drSp) || (!ddlsdd) || (!hltsdd) || (!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   fIsHLTmodeC=hltsdd->IsHLTmodeC();
562   AliITSCalibration* cal;
563   Float_t avegain=0.;
564   Float_t nGdAnodes=0;
565   for(Int_t iddl=0; iddl<AliITSDDLModuleMapSDD::GetNDDLs(); iddl++){
566     for(Int_t icar=0; icar<AliITSDDLModuleMapSDD::GetNModPerDDL();icar++){
567       Int_t iMod=fDDLMapSDD->GetModuleNumber(iddl,icar);
568       if(iMod==-1) continue;
569       Int_t i=iMod - fgkDefaultNModulesSPD;
570       cal = (AliITSCalibration*) calSDD->At(i);
571       Int_t i0=2*i;
572       Int_t i1=1+2*i;
573       for(Int_t iAnode=0;iAnode< ((AliITSCalibrationSDD*)cal)->NOfAnodes(); iAnode++){
574         if(((AliITSCalibrationSDD*)cal)->IsBadChannel(iAnode)) continue;
575         avegain+= ((AliITSCalibrationSDD*)cal)->GetChannelGain(iAnode);
576         nGdAnodes++;
577       }
578       AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
579       //      AliITSMapSDD* ma0 = (AliITSMapSDD*)mapAn->At(i0);
580       AliITSMapSDD* mt0 = (AliITSMapSDD*)mapT->At(i0);
581       AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
582       //      AliITSMapSDD* ma1 = (AliITSMapSDD*)mapAn->At(i1);
583       AliITSMapSDD* mt1 = (AliITSMapSDD*)mapT->At(i1);
584       cal->SetDriftSpeed(0,arr0);
585       cal->SetDriftSpeed(1,arr1);
586 //       cal->SetMapA(0,ma0);
587 //       cal->SetMapA(1,ma1);
588       cal->SetMapT(0,mt0);
589       cal->SetMapT(1,mt1);
590       SetCalibrationModel(iMod, cal);
591     }
592   }
593   if(nGdAnodes) fAveGainSDD=avegain/nGdAnodes;
594   return kTRUE;
595 }
596
597
598 //______________________________________________________________________
599 Bool_t AliITSDetTypeRec::GetCalibrationSSD(Bool_t cacheStatus) {
600   // Get SSD calibration objects from OCDB
601   //  AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD");
602
603   AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
604   AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
605   AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
606
607   if(!entryNoiseSSD || !entryGainSSD || !entryBadChannelsSSD){
608     AliFatal("SSD Calibration object retrieval failed! ");
609     return kFALSE;
610   }     
611
612   TObject *emptyssd = 0; TString ssdobjectname = 0;
613   AliITSNoiseSSDv2 *noiseSSD = new AliITSNoiseSSDv2();  
614   emptyssd = (TObject *)entryNoiseSSD->GetObject();
615   ssdobjectname = emptyssd->GetName();
616   if(ssdobjectname=="TObjArray") {
617     TObjArray *noiseSSDOld = (TObjArray *)entryNoiseSSD->GetObject();
618     ReadOldSSDNoise(noiseSSDOld, noiseSSD);
619   }
620   else if(ssdobjectname=="AliITSNoiseSSDv2")
621     noiseSSD = (AliITSNoiseSSDv2 *)entryNoiseSSD->GetObject();
622   if(!cacheStatus)entryNoiseSSD->SetObject(NULL);
623   entryNoiseSSD->SetOwner(kTRUE);
624
625   AliITSGainSSDv2 *gainSSD = new AliITSGainSSDv2();
626   emptyssd = (TObject *)entryGainSSD->GetObject();
627   ssdobjectname = emptyssd->GetName();
628   if(ssdobjectname=="Gain") {
629     TObjArray *gainSSDOld = (TObjArray *)entryGainSSD->GetObject();
630     ReadOldSSDGain(gainSSDOld, gainSSD);
631   }
632   else if(ssdobjectname=="AliITSGainSSDv2")
633     gainSSD = (AliITSGainSSDv2 *)entryGainSSD->GetObject();
634   if(!cacheStatus)entryGainSSD->SetObject(NULL);
635   entryGainSSD->SetOwner(kTRUE);
636
637   AliITSBadChannelsSSDv2 *badChannelsSSD = new AliITSBadChannelsSSDv2();
638   emptyssd = (TObject *)entryBadChannelsSSD->GetObject();
639   ssdobjectname = emptyssd->GetName();
640   if(ssdobjectname=="TObjArray") {
641     TObjArray *badChannelsSSDOld = (TObjArray *)entryBadChannelsSSD->GetObject();
642     ReadOldSSDBadChannels(badChannelsSSDOld, badChannelsSSD);
643   }
644   else if(ssdobjectname=="AliITSBadChannelsSSDv2")
645     badChannelsSSD = (AliITSBadChannelsSSDv2*)entryBadChannelsSSD->GetObject();
646   if(!cacheStatus)entryBadChannelsSSD->SetObject(NULL);
647   entryBadChannelsSSD->SetOwner(kTRUE);
648
649   // DB entries are deleted. In this way metadeta objects are deleted as well
650   if(!cacheStatus){
651     delete entryNoiseSSD;
652     delete entryGainSSD;
653     delete entryBadChannelsSSD;
654   }
655
656   if ((!noiseSSD)|| (!gainSSD)|| (!badChannelsSSD)) {
657     AliWarning("Can not get SSD calibration from calibration database !");
658     return kFALSE;
659   }
660
661   fSSDCalibration->SetNoise(noiseSSD);
662   fSSDCalibration->SetGain(gainSSD);
663   fSSDCalibration->SetBadChannels(badChannelsSSD);
664   //fSSDCalibration->FillBadChipMap();
665
666   return kTRUE;
667 }
668
669 //________________________________________________________________
670 void AliITSDetTypeRec::SetDefaultClusterFindersV2(Bool_t rawdata){
671
672   //Set defaults for cluster finder V2
673
674   if(!GetITSgeom()){
675     Warning("SetDefaults","Null pointer to AliITSgeom !");
676     return;
677   }
678
679   AliITSClusterFinder *clf; 
680
681   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
682     //SPD
683     if(dettype==0){
684       if(!GetReconstructionModel(dettype)){
685         clf = new AliITSClusterFinderV2SPD(this);
686         clf->InitGeometry();
687         if(!rawdata) clf->SetDigits(DigitsAddress(0));
688         SetReconstructionModel(dettype,clf);
689
690       }
691     }
692     //SDD
693     if(dettype==1){
694       if(!GetReconstructionModel(dettype)){
695         clf = new AliITSClusterFinderV2SDD(this);
696         clf->InitGeometry();
697         if(!rawdata) clf->SetDigits(DigitsAddress(1));
698         SetReconstructionModel(dettype,clf);
699       }
700
701     }
702
703     //SSD
704     if(dettype==2){
705       if(!GetReconstructionModel(dettype)){
706         clf = new AliITSClusterFinderV2SSD(this);
707         clf->InitGeometry();
708         if(!rawdata) clf->SetDigits(DigitsAddress(2));
709         SetReconstructionModel(dettype,clf);
710       }
711     }
712
713  }
714    
715 }
716 //______________________________________________________________________
717 void AliITSDetTypeRec::MakeBranch(TTree* tree, Option_t* option){
718
719   //Creates branches for clusters and recpoints
720   Bool_t cR = (strstr(option,"R")!=0);
721   Bool_t cRF = (strstr(option,"RF")!=0);
722   
723   if(cRF)cR = kFALSE;
724
725   if(cR) MakeBranchR(tree);
726   if(cRF) MakeBranchRF(tree);
727
728 }
729
730 //___________________________________________________________________
731 void AliITSDetTypeRec::ResetDigits(){
732   // Reset number of digits and the digits array for the ITS detector.
733   
734   if(!fDigits) return;
735   for(Int_t i=0;i<fgkNdettypes;i++){
736     ResetDigits(i);
737   }
738 }
739 //___________________________________________________________________
740 void AliITSDetTypeRec::ResetDigits(Int_t branch){
741   // Reset number of digits and the digits array for this branch.
742   
743   if(fDigits->At(branch)) ((TClonesArray*)fDigits->At(branch))->Clear();
744
745 }
746 //__________________________________________________________________
747 void AliITSDetTypeRec::MakeBranchR(TTree *treeR, Option_t *opt){
748
749   //Creates tree branches for recpoints
750   // Inputs:
751   //      cont char *file  File name where RecPoints branch is to be written
752   //                       to. If blank it write the SDigits to the same
753   //                       file in which the Hits were found.
754
755   Int_t buffsz = 4000;
756   char branchname[30];
757
758   // only one branch for rec points for all detector types
759   Bool_t oFast= (strstr(opt,"Fast")!=0);
760   
761   Char_t detname[10] = "ITS";
762  
763   
764   if(oFast){
765     sprintf(branchname,"%sRecPointsF",detname);
766   } else {
767     sprintf(branchname,"%sRecPoints",detname);
768   }
769   
770   if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
771   if (treeR)
772     MakeBranchInTree(treeR,branchname,0,&fRecPoints,buffsz,99);
773 }
774 //______________________________________________________________________
775 void AliITSDetTypeRec::SetTreeAddressR(TTree* const treeR){
776     // Set branch address for the Reconstructed points Trees.
777     // Inputs:
778     //      TTree *treeR   Tree containing the RecPoints.
779     // Outputs:
780     //      none.
781     // Return:
782
783    char branchname[30];
784    Char_t namedet[10]="ITS";
785
786    if(!treeR) return;
787    if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
788    TBranch *branch;
789    sprintf(branchname,"%sRecPoints",namedet);
790    branch = treeR->GetBranch(branchname);
791    if (branch) {
792       branch->SetAddress(&fRecPoints);
793     } 
794     else {
795       sprintf(branchname,"%sRecPointsF",namedet);
796       branch = treeR->GetBranch(branchname);
797       if (branch) {
798         branch->SetAddress(&fRecPoints);
799       }
800    }
801 }
802 //____________________________________________________________________
803 void AliITSDetTypeRec::AddRecPoint(const AliITSRecPoint &r){
804     // Add a reconstructed space point to the list
805     // Inputs:
806     //      const AliITSRecPoint &r RecPoint class to be added to the tree
807     //                              of reconstructed points TreeR.
808     // Outputs:
809     //      none.
810     // Return:
811     //      none.
812
813     TClonesArray &lrecp = *fRecPoints;
814     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
815 }
816
817 //______________________________________________________________________
818 void AliITSDetTypeRec::DigitsToRecPoints(TTree *treeD,TTree *treeR,Int_t lastentry,Option_t *opt, Int_t optCluFind){
819   // cluster finding and reconstruction of space points
820   // the condition below will disappear when the geom class will be
821   // initialized for all versions - for the moment it is only for v5 !
822   // 7 is the SDD beam test version
823   // Inputs:
824   //      TTree *treeD     Digits tree
825   //      TTree *treeR     Clusters tree
826   //      Int_t lastentry  Offset for module when not all of the modules
827   //                       are processed.
828   //      Option_t *opt    String indicating which ITS sub-detectors should
829   //                       be processed. If ="All" then all of the ITS
830   //                       sub detectors are processed.
831
832   const char *all = strstr(opt,"All");
833   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
834                         strstr(opt,"SSD")};
835   if(optCluFind==0){
836     SetDefaultClusterFindersV2();
837     AliInfo("V2 cluster finder has been selected \n");
838   }else{
839     SetDefaultClusterFindersV2();
840     AliInfo("Cluster Finder Option not implemented, V2 cluster finder will be used \n");    
841   }
842
843   
844   // Reset Fast-OR fired map
845   ResetFastOrFiredMap();
846
847   if (all || det[0]) { // SPD present
848     // Get the FO signals for this event
849     AliRunLoader* runLoader = AliRunLoader::Instance();
850     AliITSLoader* itsLoader = (AliITSLoader*) runLoader->GetLoader("ITSLoader");
851     if (!itsLoader) {
852       AliError("ITS loader is NULL.");
853     }
854     else {
855       AliBaseLoader* foLoader = itsLoader->GetFOSignalsLoader();
856       if (!foLoader) {
857         AliError("FO signals base loader not retrieved.");
858       }
859       else {
860         foLoader->Load();
861         fFOSignals = (AliITSFOSignalsSPD*) foLoader->Get();
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   }
899 }
900 //______________________________________________________________________
901 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader,TTree *treeR,Option_t *opt){
902   // cluster finding and reconstruction of space points
903   // the condition below will disappear when the geom class will be
904   // initialized for all versions - for the moment it is only for v5 !
905   // 7 is the SDD beam test version
906   // Inputs:
907   //      AliRawReader *rawReader  Pointer to the raw-data reader
908   //      TTree *treeR             Clusters tree
909   // Outputs:
910   //      none.
911   // Return:
912   //      none.
913   const char *all = strstr(opt,"All");
914   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
915                         strstr(opt,"SSD")};
916   
917   // Reset Fast-OR fired map
918   ResetFastOrFiredMap();
919   
920   AliITSClusterFinder *rec     = 0;
921   Int_t id=0;
922
923   TClonesArray *array=new TClonesArray("AliITSRecPoint",1000);
924   TBranch *branch = treeR->Branch("ITSRecPoints",&array);
925   delete array;
926  
927   TClonesArray** clusters = new TClonesArray*[GetITSgeom()->GetIndexMax()]; 
928   for (Int_t iModule = 0; iModule < GetITSgeom()->GetIndexMax(); iModule++) {
929     clusters[iModule] = NULL;
930   }
931   for(id=0;id<3;id++){
932     if (!all && !det[id]) continue;
933     rec = (AliITSClusterFinder*)GetReconstructionModel(id);
934     if (!rec)
935       AliFatal("The reconstruction class was not instantiated");
936     rec->SetDetTypeRec(this);
937     rec->RawdataToClusters(rawReader,clusters);    
938   } 
939   Int_t nClusters =0;
940   TClonesArray *emptyArray=new TClonesArray("AliITSRecPoint");
941   for(Int_t iModule=0;iModule<GetITSgeom()->GetIndexMax();iModule++){
942     id = GetITSgeom()->GetModuleType(iModule);
943     if (!all && !det[id]) continue;
944     array = clusters[iModule];
945     if(!array){
946       AliDebug(1,Form("data for module %d missing!",iModule));
947       array = emptyArray;
948     }
949     branch->SetAddress(&array);
950     treeR->Fill();
951     nClusters+=array->GetEntriesFast();
952
953     if (array != emptyArray) {
954       array->Delete();
955       delete array;
956     }
957   }
958   delete emptyArray;
959
960   delete[] clusters;
961   Info("DigitsToRecPoints", "total number of found recpoints in ITS: %d\n", 
962        nClusters);
963   
964   // Remove PIT in-active chips from Fast-OR fired map
965   if (all || det[0]) { // SPD present
966     RemoveFastOrFiredInActive();
967   }  
968 }
969
970 //______________________________________________________________________
971 void AliITSDetTypeRec::ReadOldSSDNoise(const TObjArray *array, 
972                                        AliITSNoiseSSDv2 *noiseSSD) {
973   //Reads the old SSD calibration object and converts it to the new format
974   const Int_t fgkSSDSTRIPSPERMODULE = 1536;
975   const Int_t fgkSSDPSIDESTRIPSPERMODULE = 768;
976
977   Int_t gNMod = array->GetEntries();
978   cout<<"Converting old calibration object for noise..."<<endl;
979
980   //NOISE
981   Double_t noise = 0.0;
982   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
983     AliITSNoiseSSD *noiseModule = (AliITSNoiseSSD*) (array->At(iModule));
984     for(Int_t iStrip = 0; iStrip < fgkSSDSTRIPSPERMODULE; iStrip++) {
985       noise = (iStrip < fgkSSDPSIDESTRIPSPERMODULE) ? noiseModule->GetNoiseP(iStrip) : noiseModule->GetNoiseN(1535 - iStrip);
986       if(iStrip < fgkSSDPSIDESTRIPSPERMODULE)
987         noiseSSD->AddNoiseP(iModule,iStrip,noise);
988       if(iStrip >= fgkSSDPSIDESTRIPSPERMODULE)
989         noiseSSD->AddNoiseN(iModule,1535 - iStrip,noise);
990     }//loop over strips
991   }//loop over modules      
992 }
993
994 //______________________________________________________________________
995 void AliITSDetTypeRec::ReadOldSSDBadChannels(const TObjArray *array, 
996                                              AliITSBadChannelsSSDv2 *badChannelsSSD) {
997   //Reads the old SSD calibration object and converts it to the new format
998   Int_t gNMod = array->GetEntries();
999   cout<<"Converting old calibration object for bad channels..."<<endl;
1000   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1001     //for (Int_t iModule = 0; iModule < 1; iModule++) {
1002     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (array->At(iModule));
1003     TArrayI arrayPSide = bad->GetBadPChannelsList();
1004     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++) 
1005       badChannelsSSD->AddBadChannelP(iModule,
1006                                      iPCounter,
1007                                      (Char_t)arrayPSide.At(iPCounter));
1008         
1009     TArrayI arrayNSide = bad->GetBadNChannelsList();
1010     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++) 
1011       badChannelsSSD->AddBadChannelN(iModule,
1012                                      iNCounter,
1013                                      (Char_t)arrayNSide.At(iNCounter));
1014     
1015   }//loop over modules      
1016 }
1017
1018 //______________________________________________________________________
1019 void AliITSDetTypeRec::ReadOldSSDGain(const TObjArray *array, 
1020                                       AliITSGainSSDv2 *gainSSD) {
1021   //Reads the old SSD calibration object and converts it to the new format
1022
1023   Int_t gNMod = array->GetEntries();
1024   cout<<"Converting old calibration object for gain..."<<endl;
1025
1026   //GAIN
1027   for (Int_t iModule = 0; iModule < gNMod; iModule++) {
1028     AliITSGainSSD *gainModule = (AliITSGainSSD*) (array->At(iModule));
1029     TArrayF arrayPSide = gainModule->GetGainP();
1030     for(Int_t iPCounter = 0; iPCounter < arrayPSide.GetSize(); iPCounter++)
1031       gainSSD->AddGainP(iModule,
1032                         iPCounter,
1033                         arrayPSide.At(iPCounter));
1034     TArrayF arrayNSide = gainModule->GetGainN();
1035     for(Int_t iNCounter = 0; iNCounter < arrayNSide.GetSize(); iNCounter++)
1036       gainSSD->AddGainN(iModule,
1037                         iNCounter,
1038                         arrayNSide.At(iNCounter));
1039   }//loop over modules 
1040 }
1041 //______________________________________________________________________
1042 void AliITSDetTypeRec::RemoveFastOrFiredInActive() {
1043   // Removes the chips that were in-active in the pixel trigger (from fast-or fired map)
1044   if (fTriggerConditions==NULL) {
1045     AliError("Pixel trigger conditions are missing.");
1046     return;
1047   }
1048   Int_t eq   = -1;
1049   Int_t hs   = -1;
1050   Int_t chip = -1;
1051   while (fTriggerConditions->GetNextInActiveChip(eq,hs,chip)) {
1052     UInt_t chipKey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1053     fFastOrFiredMap.SetBitNumber(chipKey,kFALSE);
1054   }
1055 }
1056 //______________________________________________________________________
1057 void AliITSDetTypeRec::SetFastOrFiredMapOnline(UInt_t eq, UInt_t hs, UInt_t chip) {
1058   // Set fast-or fired map for this chip
1059   Int_t chipKey = AliITSRawStreamSPD::GetOfflineChipKeyFromOnline(eq,hs,chip);
1060   return SetFastOrFiredMap(chipKey);
1061 }
1062
1063