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