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