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