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