]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDetTypeRec.cxx
updated
[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,"ZS")) 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   // dead pixel are not used for local reconstruction
399   AliCDBEntry *entrySPD = AliCDBManager::Instance()->Get("ITS/Calib/SPDNoisy");
400   AliCDBEntry *entrySDD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSDD");
401  
402  //  AliCDBEntry *entrySSD = AliCDBManager::Instance()->Get("ITS/Calib/CalibSSD");
403     AliCDBEntry *entryNoiseSSD = AliCDBManager::Instance()->Get("ITS/Calib/NoiseSSD");
404     AliCDBEntry *entryPedestalSSD = AliCDBManager::Instance()->Get("ITS/Calib/PedestalSSD");
405     AliCDBEntry *entryGainSSD = AliCDBManager::Instance()->Get("ITS/Calib/GainSSD");
406     AliCDBEntry *entryBadChannelsSSD = AliCDBManager::Instance()->Get("ITS/Calib/BadChannelsSSD");
407
408   AliCDBEntry *entry2SPD = AliCDBManager::Instance()->Get("ITS/Calib/RespSPD");
409   AliCDBEntry *entry2SDD = AliCDBManager::Instance()->Get("ITS/Calib/RespSDD");
410   AliCDBEntry *entry2SSD = AliCDBManager::Instance()->Get("ITS/Calib/RespSSD");
411   AliCDBEntry *drSpSDD = AliCDBManager::Instance()->Get("ITS/Calib/DriftSpeedSDD");
412   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
413   AliCDBEntry *mapASDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsAnodeSDD");
414   AliCDBEntry *mapTSDD = AliCDBManager::Instance()->Get("ITS/Calib/MapsTimeSDD");
415
416   if(!entrySPD || !entrySDD || !entryNoiseSSD || !entryGainSSD || 
417      !entryPedestalSSD || !entryBadChannelsSSD || 
418      !entry2SPD || !entry2SDD || !entry2SSD || !drSpSDD || !ddlMapSDD || !mapASDD || !mapTSDD){
419     AliFatal("Calibration object retrieval failed! ");
420     return kFALSE;
421   }     
422
423   TObjArray *calSPD = (TObjArray *)entrySPD->GetObject();
424   if(!cacheStatus)entrySPD->SetObject(NULL);
425   entrySPD->SetOwner(kTRUE);
426  
427   AliITSresponseSPD *pSPD = (AliITSresponseSPD*)entry2SPD->GetObject();
428   if(!cacheStatus)entry2SPD->SetObject(NULL);
429   entry2SPD->SetOwner(kTRUE);
430     
431   TObjArray *calSDD = (TObjArray *)entrySDD->GetObject();
432   if(!cacheStatus)entrySDD->SetObject(NULL);
433   entrySDD->SetOwner(kTRUE);
434  
435   AliITSresponseSDD *pSDD = (AliITSresponseSDD*)entry2SDD->GetObject();
436   if(!cacheStatus)entry2SDD->SetObject(NULL);
437   entry2SDD->SetOwner(kTRUE);
438
439   TObjArray *drSp = (TObjArray *)drSpSDD->GetObject();
440   if(!cacheStatus)drSpSDD->SetObject(NULL);
441   drSpSDD->SetOwner(kTRUE);
442
443   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
444   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
445   ddlMapSDD->SetOwner(kTRUE);
446
447   TObjArray *mapAn = (TObjArray *)mapASDD->GetObject();
448   if(!cacheStatus)mapASDD->SetObject(NULL);
449   mapASDD->SetOwner(kTRUE);
450
451   TObjArray *mapT = (TObjArray *)mapTSDD->GetObject();
452   if(!cacheStatus)mapTSDD->SetObject(NULL);
453   mapTSDD->SetOwner(kTRUE);
454
455   TObjArray *noiseSSD = (TObjArray *)entryNoiseSSD->GetObject();
456   if(!cacheStatus)entryNoiseSSD->SetObject(NULL);
457   entryNoiseSSD->SetOwner(kTRUE);
458
459   TObjArray *pedestalSSD = (TObjArray *)entryPedestalSSD->GetObject();
460   if(!cacheStatus)entryPedestalSSD->SetObject(NULL);
461   entryPedestalSSD->SetOwner(kTRUE);
462
463   TObjArray *gainSSD = (TObjArray *)entryGainSSD->GetObject();
464   if(!cacheStatus)entryGainSSD->SetObject(NULL);
465   entryGainSSD->SetOwner(kTRUE);
466
467   TObjArray *badchannelsSSD = (TObjArray *)entryBadChannelsSSD->GetObject();
468   if(!cacheStatus)entryBadChannelsSSD->SetObject(NULL);
469   entryBadChannelsSSD->SetOwner(kTRUE);
470
471   AliITSresponseSSD *pSSD = (AliITSresponseSSD*)entry2SSD->GetObject();
472   if(!cacheStatus)entry2SSD->SetObject(NULL);
473   entry2SSD->SetOwner(kTRUE);
474
475   // DB entries are deleted. In this way metadeta objects are deleted as well
476   if(!cacheStatus){
477     delete entrySPD;
478     delete entrySDD;
479     delete entryNoiseSSD;
480     delete entryPedestalSSD;
481     delete entryGainSSD;
482     delete entryBadChannelsSSD;
483     delete entry2SPD;
484     delete entry2SDD;
485     delete entry2SSD;
486     delete mapASDD;
487     delete mapTSDD;
488     delete drSpSDD;
489     delete ddlMapSDD;
490   }
491
492   if ((!pSPD)||(!pSDD)||(!pSSD) || (!calSPD) || (!calSDD) || (!drSp) || (!ddlsdd)
493       || (!mapAn) || (!mapT) || (!noiseSSD)|| (!gainSSD)|| (!badchannelsSSD)) {
494     AliWarning("Can not get calibration from calibration database !");
495     return kFALSE;
496   }
497
498   fNMod[0] = calSPD->GetEntries();
499   fNMod[1] = calSDD->GetEntries();
500   fNMod[2] = noiseSSD->GetEntries();
501   AliInfo(Form("%i SPD, %i SDD and %i SSD in calibration database",
502                fNMod[0], fNMod[1], fNMod[2]));
503   AliITSCalibration* cal;
504   for (Int_t i=0; i<fNMod[0]; i++) {
505     cal = (AliITSCalibration*) calSPD->At(i);
506     cal->SetResponse((AliITSresponse*)pSPD);
507     SetCalibrationModel(i, cal);
508   }
509
510   fDDLMapSDD->SetDDLMap(ddlsdd);
511   for(Int_t iddl=0; iddl<AliITSDDLModuleMapSDD::GetNDDLs(); iddl++){
512     for(Int_t icar=0; icar<AliITSDDLModuleMapSDD::GetNModPerDDL();icar++){
513       Int_t iMod=fDDLMapSDD->GetModuleNumber(iddl,icar);
514       if(iMod==-1) continue;
515       Int_t i=iMod - fgkDefaultNModulesSPD;
516       cal = (AliITSCalibration*) calSDD->At(i);
517       cal->SetResponse((AliITSresponse*)pSDD);
518       Int_t i0=2*i;
519       Int_t i1=1+2*i;
520       AliITSDriftSpeedArraySDD* arr0 = (AliITSDriftSpeedArraySDD*) drSp->At(i0);
521       AliITSMapSDD* ma0 = (AliITSMapSDD*)mapAn->At(i0);
522       AliITSMapSDD* mt0 = (AliITSMapSDD*)mapT->At(i0);
523       AliITSDriftSpeedArraySDD* arr1 = (AliITSDriftSpeedArraySDD*) drSp->At(i1);
524       AliITSMapSDD* ma1 = (AliITSMapSDD*)mapAn->At(i1);
525       AliITSMapSDD* mt1 = (AliITSMapSDD*)mapT->At(i1);
526       cal->SetDriftSpeed(0,arr0);
527       cal->SetDriftSpeed(1,arr1);
528       cal->SetMapA(0,ma0);
529       cal->SetMapA(1,ma1);
530       cal->SetMapT(0,mt0);
531       cal->SetMapT(1,mt1);
532       SetCalibrationModel(iMod, cal);
533     }
534   }
535   for (Int_t i=0; i<fNMod[2]; i++) {
536
537     AliITSCalibrationSSD *calibSSD = new AliITSCalibrationSSD();
538     calibSSD->SetResponse((AliITSresponse*)pSSD);
539     
540     AliITSNoiseSSD *noise = (AliITSNoiseSSD*) (noiseSSD->At(i));
541     calibSSD->SetNoise(noise);
542     AliITSPedestalSSD *pedestal = (AliITSPedestalSSD*) (pedestalSSD->At(i));
543     calibSSD->SetPedestal(pedestal);
544     AliITSGainSSD *gain = (AliITSGainSSD*) (gainSSD->At(i));
545     calibSSD->SetGain(gain);
546     AliITSBadChannelsSSD *bad = (AliITSBadChannelsSSD*) (badchannelsSSD->At(i));
547     calibSSD->SetBadChannels(bad);
548
549     Int_t iMod = i + fgkDefaultNModulesSPD + fgkDefaultNModulesSDD;
550     SetCalibrationModel(iMod, calibSSD);
551  }
552
553   return kTRUE;
554 }
555
556
557 //________________________________________________________________
558 void AliITSDetTypeRec::SetDefaultClusterFinders(){
559   
560   //set defaults for standard cluster finder
561
562   if(!GetITSgeom()){
563     Warning("SetDefaults","null pointer to AliITSgeom!");
564     return;
565   }
566
567   AliITSClusterFinder *clf; 
568
569   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
570     //SPD
571     if(dettype==0){
572       if(!GetReconstructionModel(dettype)){
573         TClonesArray *dig0 = DigitsAddress(0);
574         TClonesArray *rec0 = ClustersAddress(0);
575         clf = new AliITSClusterFinderSPD(this,dig0,rec0);
576         SetReconstructionModel(dettype,clf);
577
578       }
579     }
580    
581     //SDD
582     if(dettype==1){
583       if(!GetReconstructionModel(dettype)){
584         TClonesArray *dig1 = DigitsAddress(1);
585         TClonesArray *rec1 = ClustersAddress(1);
586         clf = new AliITSClusterFinderSDD(this,dig1,rec1);
587         SetReconstructionModel(dettype,clf);
588       }
589
590     }
591     //SSD
592     if(dettype==2){
593       if(!GetReconstructionModel(dettype)){
594         TClonesArray* dig2 = DigitsAddress(2);
595         clf = new AliITSClusterFinderSSD(this,dig2);
596         SetReconstructionModel(dettype,clf);
597       }
598     }
599
600  }
601  
602   
603 }
604
605 //________________________________________________________________
606 void AliITSDetTypeRec::SetDefaultClusterFindersV2(Bool_t rawdata){
607
608   //Set defaults for cluster finder V2
609
610   if(!GetITSgeom()){
611     Warning("SetDefaults","Null pointer to AliITSgeom !");
612     return;
613   }
614
615   AliITSClusterFinder *clf; 
616
617   for(Int_t dettype=0;dettype<fgkNdettypes;dettype++){
618     //SPD
619     if(dettype==0){
620       if(!GetReconstructionModel(dettype)){
621         clf = new AliITSClusterFinderV2SPD(this);
622         clf->InitGeometry();
623         if(!rawdata) clf->SetDigits(DigitsAddress(0));
624         SetReconstructionModel(dettype,clf);
625
626       }
627     }
628     //SDD
629     if(dettype==1){
630       if(!GetReconstructionModel(dettype)){
631         clf = new AliITSClusterFinderV2SDD(this);
632         clf->InitGeometry();
633         if(!rawdata) clf->SetDigits(DigitsAddress(1));
634         SetReconstructionModel(dettype,clf);
635       }
636
637     }
638
639     //SSD
640     if(dettype==2){
641       if(!GetReconstructionModel(dettype)){
642         clf = new AliITSClusterFinderV2SSD(this);
643         clf->InitGeometry();
644         if(!rawdata) clf->SetDigits(DigitsAddress(2));
645         SetReconstructionModel(dettype,clf);
646       }
647     }
648
649  }
650    
651 }
652 //______________________________________________________________________
653 void AliITSDetTypeRec::MakeBranch(TTree* tree, Option_t* option){
654
655   //Creates branches for clusters and recpoints
656   Bool_t cR = (strstr(option,"R")!=0);
657   Bool_t cRF = (strstr(option,"RF")!=0);
658   
659   if(cRF)cR = kFALSE;
660
661   if(cR) MakeBranchR(tree);
662   if(cRF) MakeBranchRF(tree);
663
664 }
665
666 //___________________________________________________________________
667 void AliITSDetTypeRec::AddCluster(Int_t id, AliITSRawCluster *c){
668
669   // Adds a raw cluster to the list
670   TClonesArray &lc = *((TClonesArray*)fCtype->At(id));  
671   switch(id){
672   case 0:
673     new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
674     break;
675   case 1:
676     new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
677     break;
678   case 2:
679     new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
680     break;
681   } 
682 }
683 //___________________________________________________________________
684 void AliITSDetTypeRec::ResetDigits(){
685   // Reset number of digits and the digits array for the ITS detector.
686   
687   if(!fDigits) return;
688   for(Int_t i=0;i<fgkNdettypes;i++){
689     ResetDigits(i);
690   }
691 }
692 //___________________________________________________________________
693 void AliITSDetTypeRec::ResetDigits(Int_t branch){
694   // Reset number of digits and the digits array for this branch.
695   
696   if(fDigits->At(branch)) ((TClonesArray*)fDigits->At(branch))->Clear();
697   if(fNdtype) fNdtype[branch]=0;
698
699 }
700
701 //__________________________________________________________________
702 void AliITSDetTypeRec::ResetClusters(){
703
704   //Resets number of clusters and the cluster array 
705   for(Int_t i=0;i<fgkNdettypes;i++){
706     ResetClusters(i);
707   }
708 }
709
710 //__________________________________________________________________
711 void AliITSDetTypeRec::ResetClusters(Int_t i){
712
713   //Resets number of clusters and the cluster array for this branch
714
715   if (fCtype->At(i))    ((TClonesArray*)fCtype->At(i))->Clear();
716   if (fNctype)  fNctype[i]=0;
717 }
718 //__________________________________________________________________
719 void AliITSDetTypeRec::MakeBranchR(TTree *treeR, Option_t *opt){
720
721   //Creates tree branches for recpoints
722   // Inputs:
723   //      cont char *file  File name where RecPoints branch is to be written
724   //                       to. If blank it write the SDigits to the same
725   //                       file in which the Hits were found.
726
727   Int_t buffsz = 4000;
728   char branchname[30];
729
730   // only one branch for rec points for all detector types
731   Bool_t oFast= (strstr(opt,"Fast")!=0);
732   
733   Char_t detname[10] = "ITS";
734  
735   
736   if(oFast){
737     sprintf(branchname,"%sRecPointsF",detname);
738   } else {
739     sprintf(branchname,"%sRecPoints",detname);
740   }
741   
742   if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
743   if (treeR)
744     MakeBranchInTree(treeR,branchname,0,&fRecPoints,buffsz,99);
745 }
746 //______________________________________________________________________
747 void AliITSDetTypeRec::SetTreeAddressR(TTree *treeR){
748     // Set branch address for the Reconstructed points Trees.
749     // Inputs:
750     //      TTree *treeR   Tree containing the RecPoints.
751     // Outputs:
752     //      none.
753     // Return:
754
755     char branchname[30];
756     Char_t namedet[10]="ITS";
757
758     if(!treeR) return;
759     if(fRecPoints==0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
760     TBranch *branch;
761     sprintf(branchname,"%sRecPoints",namedet);
762     branch = treeR->GetBranch(branchname);
763     if (branch) {
764       branch->SetAddress(&fRecPoints);
765     }else {
766       sprintf(branchname,"%sRecPointsF",namedet);
767       branch = treeR->GetBranch(branchname);
768       if (branch) {
769         branch->SetAddress(&fRecPoints);
770       }
771
772     }
773 }
774 //____________________________________________________________________
775 void AliITSDetTypeRec::AddRecPoint(const AliITSRecPoint &r){
776     // Add a reconstructed space point to the list
777     // Inputs:
778     //      const AliITSRecPoint &r RecPoint class to be added to the tree
779     //                              of reconstructed points TreeR.
780     // Outputs:
781     //      none.
782     // Return:
783     //      none.
784
785     TClonesArray &lrecp = *fRecPoints;
786     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
787 }
788
789 //______________________________________________________________________
790 void AliITSDetTypeRec::DigitsToRecPoints(TTree *treeD,TTree *treeR,Int_t lastentry,Option_t *opt, Bool_t v2){
791   // cluster finding and reconstruction of space points
792   // the condition below will disappear when the geom class will be
793   // initialized for all versions - for the moment it is only for v5 !
794   // 7 is the SDD beam test version
795   // Inputs:
796   //      TTree *treeD     Digits tree
797   //      TTree *treeR     Clusters tree
798   //      Int_t lastentry  Offset for module when not all of the modules
799   //                       are processed.
800   //      Option_t *opt    String indicating which ITS sub-detectors should
801   //                       be processed. If ="All" then all of the ITS
802   //                       sub detectors are processed.
803
804   const char *all = strstr(opt,"All");
805   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
806                         strstr(opt,"SSD")};
807   if(!v2) {
808     SetDefaultClusterFinders();
809     AliInfo("Original cluster finder has been selected\n");
810   }
811   else   { 
812     SetDefaultClusterFindersV2();
813     AliInfo("V2 cluster finder has been selected \n");
814   }
815
816   AliITSClusterFinder *rec     = 0;
817   Int_t id,module,first=0;
818   for(module=0;module<GetITSgeom()->GetIndexMax();module++){
819       id       = GetITSgeom()->GetModuleType(module);
820       if (!all && !det[id]) continue;
821       if(det[id]) first = GetITSgeom()->GetStartDet(id);
822       rec = (AliITSClusterFinder*)GetReconstructionModel(id);
823       TClonesArray *itsDigits  = DigitsAddress(id);
824       if (!rec)
825           AliFatal("The reconstruction class was not instanciated!");
826       ResetDigits();  // MvL: Not sure we neeed this when rereading anyways
827       if (all) {
828           treeD->GetEvent(lastentry+module);
829       }else {
830           treeD->GetEvent(lastentry+(module-first));
831       }
832       Int_t ndigits = itsDigits->GetEntriesFast();
833       if(ndigits>0){
834         rec->SetDetTypeRec(this);
835         rec->SetDigits(DigitsAddress(id));
836         //      rec->SetClusters(ClustersAddress(id));
837         rec->FindRawClusters(module);
838       } // end if
839       treeR->Fill();
840       ResetRecPoints();
841       ResetClusters();
842   } 
843 }
844 //______________________________________________________________________
845 void AliITSDetTypeRec::DigitsToRecPoints(AliRawReader* rawReader,TTree *treeR,Option_t *opt){
846   // cluster finding and reconstruction of space points
847   // the condition below will disappear when the geom class will be
848   // initialized for all versions - for the moment it is only for v5 !
849   // 7 is the SDD beam test version
850   // Inputs:
851   //      AliRawReader *rawReader  Pointer to the raw-data reader
852   //      TTree *treeR             Clusters tree
853   // Outputs:
854   //      none.
855   // Return:
856   //      none.
857   const char *all = strstr(opt,"All");
858   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
859                         strstr(opt,"SSD")};
860   AliITSClusterFinderV2 *rec     = 0;
861   Int_t id=0;
862
863   TClonesArray *array=new TClonesArray("AliITSRecPoint",1000);
864   TBranch *branch = treeR->Branch("ITSRecPoints",&array);
865   delete array;
866  
867   TClonesArray** clusters = new TClonesArray*[GetITSgeom()->GetIndexMax()]; 
868   for (Int_t iModule = 0; iModule < GetITSgeom()->GetIndexMax(); iModule++) {
869     clusters[iModule] = NULL;
870   }
871   for(id=0;id<3;id++){
872     if (!all && !det[id]) continue;
873     rec = (AliITSClusterFinderV2*)GetReconstructionModel(id);
874     if (!rec)
875       AliFatal("The reconstruction class was not instanciated");
876     rec->SetDetTypeRec(this);
877     rec->RawdataToClusters(rawReader,clusters);    
878   } 
879   Int_t nClusters =0;
880   TClonesArray *emptyArray=new TClonesArray("AliITSRecPoint");
881   for(Int_t iModule=0;iModule<GetITSgeom()->GetIndexMax();iModule++){
882     id = GetITSgeom()->GetModuleType(iModule);
883     if (!all && !det[id]) continue;
884     array = clusters[iModule];
885     if(!array){
886       AliDebug(1,Form("data for module %d missing!",iModule));
887       array = emptyArray;
888     }
889     branch->SetAddress(&array);
890     treeR->Fill();
891     nClusters+=array->GetEntriesFast();
892
893     if (array != emptyArray) {
894       array->Delete();
895       delete array;
896     }
897   }
898   delete emptyArray;
899
900   delete[] clusters;
901   Info("DigitsToRecPoints", "total number of found recpoints in ITS: %d\n", 
902        nClusters);
903   
904 }
905
906