]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUReconstructor.cxx
minor fix
[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 "AliITSUVertexer.h"
40 #include "AliMagF.h"
41
42 ClassImp(AliITSUReconstructor)
43
44 //___________________________________________________________________________
45 AliITSUReconstructor::AliITSUReconstructor() 
46 :  AliReconstructor()
47   ,fGeom(0)
48   ,fITS(0)
49   ,fClusterFinders(0)
50   ,fClusters(0)
51 {
52   // Default constructor
53
54 }
55
56 //___________________________________________________________________________
57 AliITSUReconstructor::~AliITSUReconstructor()
58 {
59   // destructor
60   //
61   if (!fGeom) return; // was not initialized
62   //
63   // same cluster finders and recpoint arrays might be attached to different layers
64   for (int i=fGeom->GetNLayers();i--;) {
65     TObject* clFinder = fClusterFinders.At(i);
66     if (clFinder) {
67       while (fClusterFinders.Remove(clFinder)) {}
68       delete clFinder;
69     }
70     //
71     delete fClusters[i];
72   }
73   delete[] fClusters;
74   //
75   delete fITS;
76   delete fGeom;
77
78
79 //______________________________________________________________________
80 void AliITSUReconstructor::Init() 
81 {
82   // Initalize this reconstructor 
83   // Note: fITS cannot be initialized here since it requires RecoParams (not available ar
84   // the moment of reconstructors initialization)
85   //
86   AliInfo("Initializing");
87   if (fGeom) AliFatal("was already done, something is wrong...");
88   //
89   fGeom = new AliITSUGeomTGeo(kTRUE,kTRUE);
90   AliITSUClusterPix::SetGeom(fGeom);
91   //  
92   AliITSUClusterizer* clusPIX = 0;
93   fClusters = new TClonesArray*[fGeom->GetNLayers()];
94   //
95   for (int ilr=fGeom->GetNLayers();ilr--;) {
96     fClusters[ilr] = 0;
97     int tpDet = fGeom->GetLayerChipTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerChipType;
98     if (tpDet == AliITSUGeomTGeo::kChipTypePix) {
99       if (!clusPIX)    clusPIX    = new AliITSUClusterizer();
100       fClusterFinders.AddAtAndExpand(clusPIX, ilr);
101       fClusters[ilr] = new TClonesArray(AliITSUClusterPix::Class());
102       //
103       // to expand the buffers to max.size
104       clusPIX->SetSegmentation((AliITSUSegmentationPix*)fGeom->GetSegmentation(ilr)); 
105       continue;
106     }
107     else {
108       AliFatal(Form("ClusterFinder for detector type %d is not defined",tpDet));
109     }
110   }
111   //
112   return;
113 }
114
115 //_____________________________________________________________________________
116 void AliITSUReconstructor::Reconstruct(TTree *digitsTree, TTree *clustersTree) const
117 {
118   // reconstruct clusters. If clustersTree is provided, write the tree
119   if (!digitsTree) return;
120   AliDebug(1,"ITSU Cluster finder (from digits tree) is initiated here \n");
121   //
122   // At the moment only pixel digits
123   TClonesArray *digArrPix = 0;
124   digitsTree->SetBranchAddress("ITSDigitsPix",&digArrPix);
125   //
126   // a new tree is created for each event: add each layer as separate branch
127   TBranch *lrBranch[fGeom->GetNLayers()];
128   //
129   for (int ilr=0;ilr<fGeom->GetNLayers();ilr++) {
130     lrBranch[ilr] = 0;
131     if (clustersTree) { // do we write clusters tree?
132       int tp = fGeom->GetLayerChipTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerChipType;
133       if (tp==AliITSUGeomTGeo::kChipTypePix) {
134         lrBranch[ilr] = clustersTree->Bronch(Form("ITSRecPoints%d",ilr),"TClonesArray",&fClusters[ilr]);
135       }
136       else {
137         AliFatal(Form("Detector type %d is not defined",tp));
138       }
139     }
140   }
141   //
142   AliITSUClusterizer* clFinder = 0;
143   AliMagF* field = dynamic_cast<AliMagF*>(TGeoGlobalMagField::Instance()->GetField());
144   double bz = 0;
145   if (field == 0) AliError("Cannot get magnetic field from TGeoGlobalMagField");
146   else bz = field->SolenoidField();
147   const AliITSURecoParam* recPar = GetRecoParam();
148   //
149   for (int ilr=0;ilr<fGeom->GetNLayers();ilr++) {
150     //
151     fClusters[ilr]->Clear();
152     clFinder = (AliITSUClusterizer*)fClusterFinders[ilr];
153     clFinder->SetSegmentation((AliITSUSegmentationPix*)fGeom->GetSegmentation(ilr));
154     clFinder->SetLayerID(ilr);
155     clFinder->SetClusters(fClusters[ilr]);
156     clFinder->SetRecoParam(recPar); // RS: Do we need to set it for every event?
157     clFinder->SetAllowDiagonalClusterization(recPar->GetAllowDiagonalClusterization(ilr));
158     clFinder->PrepareLorentzAngleCorrection(bz);
159     //
160     int modF=fGeom->GetFirstChipIndex(ilr);
161     int modL=fGeom->GetLastChipIndex(ilr)+1;
162     for (int imod=modF;imod<modL;imod++) {
163       digitsTree->GetEntry(imod);   
164       int ndig  = digArrPix->GetEntries();
165       if (!ndig) continue;
166       clFinder->SetVolID(imod);
167       clFinder->SetDigits(digArrPix);
168       clFinder->Clusterize();
169     }
170     //
171     AliITSUClusterPix::SetSortMode( AliITSUClusterPix::SortModeIdTrkYZ());
172     fClusters[ilr]->Sort();
173     AliDebug(1,Form(" -> Lr%d : %d Cluster",ilr,fClusters[ilr]->GetEntries()));
174     if (clustersTree) lrBranch[ilr]->Fill();
175   }
176   if (clustersTree) clustersTree->SetEntries();
177   //
178 }
179
180 //_____________________________________________________________________________
181 AliTracker* AliITSUReconstructor::CreateTracker() const
182 {
183   // create a ITS tracker
184   AliITSUTrackerGlo* tracker = new AliITSUTrackerGlo((AliITSUReconstructor*)this);
185  
186   return tracker;
187  
188 }
189
190 //_____________________________________________________________________________
191 AliVertexer* AliITSUReconstructor::CreateVertexer() const
192 {
193   // create a ITS vertexer
194   // to be implemented for the upgrade
195   AliInfo("Creating dummy vertexer");
196   //  AliDebug(1,"ITSU vertexer should be initiated here\n");
197   return new AliITSUVertexer();
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 }