]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDclusterizer.cxx
Go back to previous version of simulation
[u/mrichter/AliRoot.git] / TRD / AliTRDclusterizer.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$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  TRD cluster finder base class                                            //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <TROOT.h>
25 #include <TTree.h>
26 #include <TFile.h>
27 #include <TObjArray.h>
28
29 #include "AliRunLoader.h"
30 #include "AliLoader.h"
31 #include "AliLog.h"
32
33 #include "AliTRDclusterizer.h"
34 #include "AliTRDcluster.h"
35 #include "AliTRDrecPoint.h"
36 #include "AliTRDgeometry.h"
37 #include "AliTRDcalibDB.h"
38 #include "AliTRDCommonParam.h"
39
40 ClassImp(AliTRDclusterizer)
41
42 //_____________________________________________________________________________
43 AliTRDclusterizer::AliTRDclusterizer()
44   :TNamed()
45   ,fRunLoader(NULL)
46   ,fClusterTree(NULL)
47   ,fRecPoints(NULL)
48 {
49   //
50   // AliTRDclusterizer default constructor
51   //
52
53 }
54
55 //_____________________________________________________________________________
56 AliTRDclusterizer::AliTRDclusterizer(const Text_t* name, const Text_t* title)
57   :TNamed(name,title)
58   ,fRunLoader(NULL)
59   ,fClusterTree(NULL)
60   ,fRecPoints(NULL)
61 {
62   //
63   // AliTRDclusterizer constructor
64   //
65
66 }
67
68 //_____________________________________________________________________________
69 AliTRDclusterizer::AliTRDclusterizer(const AliTRDclusterizer &c)
70   :TNamed(c)
71   ,fRunLoader(NULL)
72   ,fClusterTree(NULL)
73   ,fRecPoints(NULL)
74 {
75   //
76   // AliTRDclusterizer copy constructor
77   //
78
79 }
80
81 //_____________________________________________________________________________
82 AliTRDclusterizer::~AliTRDclusterizer()
83 {
84   //
85   // AliTRDclusterizer destructor
86   //
87
88   if (fRecPoints) {
89     fRecPoints->Delete();
90     delete fRecPoints;
91   }
92
93 }
94
95 //_____________________________________________________________________________
96 AliTRDclusterizer &AliTRDclusterizer::operator=(const AliTRDclusterizer &c)
97 {
98   //
99   // Assignment operator
100   //
101
102   if (this != &c) {
103     ((AliTRDclusterizer &) c).Copy(*this);
104   }
105   return *this;
106
107 }
108
109 //_____________________________________________________________________________
110 void AliTRDclusterizer::Copy(TObject &c) const
111 {
112   //
113   // Copy function
114   //
115
116   ((AliTRDclusterizer &) c).fClusterTree = NULL;
117   ((AliTRDclusterizer &) c).fRecPoints   = NULL;  
118
119 }
120
121 //_____________________________________________________________________________
122 Bool_t AliTRDclusterizer::Open(const Char_t *name, Int_t nEvent)
123 {
124   //
125   // Opens the AliROOT file. Output and input are in the same file
126   //
127
128   TString evfoldname = AliConfig::GetDefaultEventFolderName();
129   fRunLoader         = AliRunLoader::GetRunLoader(evfoldname);
130
131   if (!fRunLoader) {
132     fRunLoader = AliRunLoader::Open(name);
133   }
134
135   if (!fRunLoader) {
136     AliError(Form("Can not open session for file %s.",name));
137     return kFALSE;
138   }
139
140   OpenInput(nEvent);
141   OpenOutput();
142
143   return kTRUE;
144
145 }
146
147 //_____________________________________________________________________________
148 Bool_t AliTRDclusterizer::OpenOutput()
149 {
150   //
151   // Open the output file
152   //
153
154   TObjArray *ioArray = 0;
155
156   AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
157   loader->MakeTree("R");
158
159   fClusterTree = loader->TreeR();
160   fClusterTree->Branch("TRDcluster","TObjArray",&ioArray,32000,0);
161
162   return kTRUE;
163
164 }
165
166 //_____________________________________________________________________________
167 Bool_t AliTRDclusterizer::OpenOutput(TTree *clusterTree)
168 {
169   //
170   // Connect the output tree
171   //
172
173   TObjArray *ioArray = 0;
174
175   fClusterTree = clusterTree;
176   fClusterTree->Branch("TRDcluster","TObjArray",&ioArray,32000,0);
177
178   return kTRUE;
179
180 }
181
182 //_____________________________________________________________________________
183 Bool_t AliTRDclusterizer::OpenInput(Int_t nEvent)
184 {
185   //
186   // Opens a ROOT-file with TRD-hits and reads in the digits-tree
187   //
188
189   // Import the Trees for the event nEvent in the file
190   fRunLoader->GetEvent(nEvent);
191   
192   return kTRUE;
193
194 }
195
196 //_____________________________________________________________________________
197 Bool_t AliTRDclusterizer::WriteClusters(Int_t det)
198 {
199   //
200   // Fills TRDcluster branch in the tree with the clusters 
201   // found in detector = det. For det=-1 writes the tree. 
202   //
203
204   if ((det <                      -1) || 
205       (det >= AliTRDgeometry::Ndet())) {
206     AliError(Form("Unexpected detector index %d.\n",det));
207     return kFALSE;
208   }
209  
210   TBranch *branch = fClusterTree->GetBranch("TRDcluster");
211   if (!branch) {
212     TObjArray *ioArray = 0;
213     branch = fClusterTree->Branch("TRDcluster","TObjArray",&ioArray,32000,0);
214   }
215
216   if ((det >=                      0) && 
217       (det <  AliTRDgeometry::Ndet())) {
218
219     Int_t nRecPoints = RecPoints()->GetEntriesFast();
220     TObjArray *detRecPoints = new TObjArray(400);
221
222     for (Int_t i = 0; i < nRecPoints; i++) {
223       AliTRDcluster *c = (AliTRDcluster *) RecPoints()->UncheckedAt(i);
224       if (det == c->GetDetector()) {
225         detRecPoints->AddLast(c);
226       }
227       else {
228         AliError("Attempt to write a cluster with unexpected detector index\n");
229       }
230     }
231
232     branch->SetAddress(&detRecPoints);
233     fClusterTree->Fill();
234
235     delete detRecPoints;
236
237     return kTRUE;
238
239   }
240
241   if (det == -1) {
242
243     AliInfo(Form("Writing the cluster tree %s for event %d."
244                 ,fClusterTree->GetName(),fRunLoader->GetEventNumber()));
245
246     if (fRecPoints) {
247
248       branch->SetAddress(&fRecPoints);
249
250       AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
251       loader->WriteRecPoints("OVERWRITE");
252   
253     }
254     else {
255
256       AliError("Cluster tree does not exist. Cannot write clusters.\n");
257       return kFALSE;
258
259     }
260
261     return kTRUE;  
262
263   }
264
265   AliError(Form("Unexpected detector index %d.\n",det));
266  
267   return kFALSE;  
268   
269 }
270
271
272 //_____________________________________________________________________________
273 AliTRDcluster* AliTRDclusterizer::AddCluster(Double_t *pos, Int_t timebin
274                                            , Int_t det, Double_t amp
275                                            , Int_t *tracks, Double_t *sig
276                                            , Int_t iType, Float_t center)
277 {
278   //
279   // Add a cluster for the TRD
280   //
281
282   AliTRDcluster *c = new AliTRDcluster();
283
284   c->SetDetector(det);
285   c->SetQ(amp);
286   c->SetX(pos[2]);
287   c->SetY(pos[0]);
288   c->SetZ(pos[1]);
289   c->SetSigmaY2(sig[0]);   
290   c->SetSigmaZ2(sig[1]);
291   c->SetLocalTimeBin(timebin);
292   c->SetCenter(center);
293
294   if (tracks) {
295     c->AddTrackIndex(tracks);
296   }
297
298   switch (iType) {
299   case 0:
300     c->Set2pad();
301     break;
302   case 1:
303     c->Set3pad();
304     break;
305   case 2:
306     c->Set4pad();
307     break;
308   case 3:
309     c->Set5pad();
310     break;
311   case 4:
312     c->SetLarge();
313     break;
314   };
315
316   RecPoints()->Add(c);
317   return c;
318
319 }
320
321 //_____________________________________________________________________________
322 Double_t AliTRDclusterizer::CalcXposFromTimebin(Float_t timebin, Int_t idet
323                                               , Int_t col, Int_t row)
324 {
325   //
326   // Calculates the local x position in the detector from the timebin, 
327   // depends on the drift velocity and t0
328   //
329   
330   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
331   if (!calibration) {
332     AliError("Cannot find calibration object");
333     return -1;
334   }
335   AliTRDCommonParam *parCom      = AliTRDCommonParam::Instance();
336   if (!parCom) {
337     AliError("Could not get common parameters\n");
338     return kFALSE;
339   }
340
341   Float_t vdrift            = calibration->GetVdrift(idet,col,row);  
342   Float_t t0                = calibration->GetT0(idet,col,row);
343   Float_t samplingFrequency = parCom->GetSamplingFrequency();
344
345   timebin -= t0;
346
347   return timebin / samplingFrequency * vdrift;
348
349 }
350
351 //_____________________________________________________________________________
352 void AliTRDclusterizer::ResetRecPoints() 
353 {
354   //
355   // Resets the list of rec points
356   //
357
358   if (fRecPoints) {
359     fRecPoints->Delete();
360   }
361
362 }
363
364 //_____________________________________________________________________________
365 TObjArray* AliTRDclusterizer::RecPoints() 
366 {
367   //
368   // Returns the list of rec points
369   //
370
371   if (!fRecPoints) {
372     fRecPoints = new TObjArray(400);
373   }
374  
375   return fRecPoints;
376
377 }