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