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