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