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