]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUReconstructor.cxx
04ca9db58e26d1da5b1e9348c58417e9d40accd8
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUReconstructor.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 /* $Id: AliITSUReconstructor.cxx 58442 2012-09-04 17:17:06Z masera $ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // class for ITS reconstruction                                              //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include "Riostream.h"
25 #include "AliITSUReconstructor.h"
26 #include "AliRun.h"
27 #include "AliRawReader.h"
28 #include "AliESDEvent.h"
29
30 #include "AliTracker.h"
31 #include "AliITSUTrackerGlo.h"
32
33 #include "AliITSUGeomTGeo.h"
34 #include "AliITSUSegmentationPix.h"
35 #include "AliITSUDigitPix.h"
36 #include "AliITSUClusterizer.h"
37 #include "AliITSUClusterPix.h"
38 #include "AliMagF.h"
39
40 ClassImp(AliITSUReconstructor)
41
42 //___________________________________________________________________________
43 AliITSUReconstructor::AliITSUReconstructor() 
44 :  AliReconstructor()
45   ,fGeom(0)
46   ,fClusterFinders(0)
47   ,fClusters(0)
48 {
49   // Default constructor
50
51 }
52
53 //___________________________________________________________________________
54 AliITSUReconstructor::~AliITSUReconstructor()
55 {
56   // destructor
57   //
58   if (!fGeom) return; // was not initialized
59   //
60   // same cluster finders and recpoint arrays might be attached to different layers
61   for (int i=fGeom->GetNLayers();i--;) {
62     TObject* clFinder = fClusterFinders.At(i);
63     if (clFinder) {
64       while (fClusterFinders.Remove(clFinder)) {}
65       delete clFinder;
66     }
67     //
68     delete fClusters[i];
69   }
70   delete[] fClusters;
71   //
72   delete fGeom;
73
74
75 //______________________________________________________________________
76 void AliITSUReconstructor::Init() 
77 {
78   // Initalize this constructor 
79   AliInfo("Initializing");
80   if (fGeom) AliFatal("was already done, something is wrong...");
81   //
82   fGeom = new AliITSUGeomTGeo(kTRUE,kTRUE);
83   AliITSUClusterPix::SetGeom(fGeom);
84   //  
85   AliITSUClusterizer* clusPIX = 0;
86   fClusters = new TClonesArray*[fGeom->GetNLayers()];
87   //
88   for (int ilr=fGeom->GetNLayers();ilr--;) {
89     fClusters[ilr] = 0;
90     int tpDet = fGeom->GetLayerDetTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerDetType;
91     if (tpDet == AliITSUGeomTGeo::kDetTypePix) {
92       if (!clusPIX)    clusPIX    = new AliITSUClusterizer();
93       fClusterFinders.AddAtAndExpand(clusPIX, ilr);
94       fClusters[ilr] = new TClonesArray(AliITSUClusterPix::Class());
95       //
96       // to expand the buffers to max.size
97       clusPIX->SetSegmentation((AliITSUSegmentationPix*)fGeom->GetSegmentation(ilr)); 
98       continue;
99     }
100     else {
101       AliFatal(Form("ClusterFinder for detector type %d is not defined",tpDet));
102     }
103   }
104   //
105   return;
106
107 }
108
109 //_____________________________________________________________________________
110 void AliITSUReconstructor::Reconstruct(TTree *digitsTree, TTree *clustersTree) const
111 {
112   // reconstruct clusters. If clustersTree is provided, write the tree
113   if (!digitsTree) return;
114   AliDebug(1,"ITSU Cluster finder (from digits tree) is initiated here \n");
115   //
116   // At the moment only pixel digits
117   TClonesArray *digArrPix = 0;
118   digitsTree->SetBranchAddress("ITSDigitsPix",&digArrPix);
119   //
120   // a new tree is created for each event: add each layer as separate branch
121   TBranch *lrBranch[fGeom->GetNLayers()];
122   //
123   for (int ilr=0;ilr<fGeom->GetNLayers();ilr++) {
124     lrBranch[ilr] = 0;
125     if (clustersTree) { // do we write clusters tree?
126       int tp = fGeom->GetLayerDetTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerDetType;
127       if (tp==AliITSUGeomTGeo::kDetTypePix) {
128         lrBranch[ilr] = clustersTree->Bronch(Form("ITSRecPoints%d",ilr),"TClonesArray",&fClusters[ilr]);
129       }
130       else {
131         AliFatal(Form("Detector type %d is not defined",tp));
132       }
133     }
134   }
135   //
136   AliITSUClusterizer* clFinder = 0;
137   AliMagF* field = dynamic_cast<AliMagF*>(TGeoGlobalMagField::Instance()->GetField());
138   double bz = 0;
139   if (field == 0) AliError("Cannot get magnetic field from TGeoGlobalMagField");
140   else bz = field->SolenoidField();
141   //
142   for (int ilr=0;ilr<fGeom->GetNLayers();ilr++) {
143     //
144     fClusters[ilr]->Clear();
145     clFinder = (AliITSUClusterizer*)fClusterFinders[ilr];
146     clFinder->SetSegmentation((AliITSUSegmentationPix*)fGeom->GetSegmentation(ilr));
147     clFinder->SetLayerID(ilr);
148     clFinder->SetClusters(fClusters[ilr]);
149     clFinder->SetRecoParam(GetRecoParam()); // RS: Do we need to set it for every event?
150     clFinder->PrepareLorentzAngleCorrection(bz);
151     //
152     int modF=fGeom->GetFirstModIndex(ilr);
153     int modL=fGeom->GetLastModIndex(ilr)+1;
154     for (int imod=modF;imod<modL;imod++) {
155       digitsTree->GetEntry(imod);   
156       int ndig  = digArrPix->GetEntries();
157       if (!ndig) continue;
158       clFinder->SetVolID(imod);
159       clFinder->SetDigits(digArrPix);
160       clFinder->Clusterize();
161     }
162     //
163     AliITSUClusterPix::SetSortMode( AliITSUClusterPix::SortModeIdTrkYZ());
164     fClusters[ilr]->Sort();
165     AliDebug(1,Form(" -> Lr%d : %d Cluster",ilr,fClusters[ilr]->GetEntries()));
166     if (clustersTree) lrBranch[ilr]->Fill();
167   }
168   if (clustersTree) clustersTree->SetEntries();
169   //
170 }
171
172 //_____________________________________________________________________________
173 AliTracker* AliITSUReconstructor::CreateTracker() const
174 {
175   // create a ITS tracker
176   AliITSUTrackerGlo* tracker = new AliITSUTrackerGlo((AliITSUReconstructor*)this);
177  
178   return tracker;
179  
180 }
181
182 //_____________________________________________________________________________
183 AliVertexer* AliITSUReconstructor::CreateVertexer() const
184 {
185   // create a ITS vertexer
186   // to be implemented for the upgrade
187
188   AliDebug(1,"ITSU vertexer should be initiated here\n");
189   return 0;
190
191 }
192
193 //_____________________________________________________________________________
194 AliTrackleter* AliITSUReconstructor::CreateMultFinder() const
195 {
196   // create the SPD trackeleter for mult. reconstruction
197   // to be implemented for the upgrade
198
199   AliDebug(1,"ITSU MultFinder  should be initiated here\n");
200   return 0;
201
202 }
203
204 //_____________________________________________________________________________
205 AliTracker* AliITSUReconstructor::CreateTrackleter() const
206 {
207   // create the SPD trackeleter (for SPD PlaneEfficiency evaluation)
208   // to be implemented for the upgrade
209
210   AliDebug(1,"ITSU Trackleter  should be initiated here\n");
211   return 0;
212
213 }
214
215 //_____________________________________________________________________________
216 Int_t AliITSUReconstructor::LoadClusters(TTree* treeRP) 
217 {
218   // read clusters from the tree, if it is provided
219   if (!treeRP) return 0;
220   for (int ilr=fGeom->GetNLayers();ilr--;) {
221     if (!fClusters[ilr]) AliFatal(Form("Clusters array for layer %d is not defined",ilr)); 
222     TBranch* br = treeRP->GetBranch(Form("ITSRecPoints%d",ilr));
223     if (!br) AliFatal(Form("Provided cluster tree does not contain branch for layer %d",ilr));
224     br->SetAddress(&fClusters[ilr]);
225   }
226   treeRP->GetEntry(0); // we are still in 1 ev/tree mode...
227   return 1;
228 }