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