]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDclusterizer.cxx
red.C is now Hdisp.C
[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::OpenInput(Int_t nEvent)
168 {
169   //
170   // Opens a ROOT-file with TRD-hits and reads in the digits-tree
171   //
172
173   // Import the Trees for the event nEvent in the file
174   fRunLoader->GetEvent(nEvent);
175   
176   return kTRUE;
177
178 }
179
180 //_____________________________________________________________________________
181 Bool_t AliTRDclusterizer::WriteClusters(Int_t det)
182 {
183   //
184   // Fills TRDcluster branch in the tree with the clusters 
185   // found in detector = det. For det=-1 writes the tree. 
186   //
187
188   if ((det <                      -1) || 
189       (det >= AliTRDgeometry::Ndet())) {
190     AliError(Form("Unexpected detector index %d.\n",det));
191     return kFALSE;
192   }
193  
194   TBranch *branch = fClusterTree->GetBranch("TRDcluster");
195   if (!branch) {
196     TObjArray *ioArray = 0;
197     branch = fClusterTree->Branch("TRDcluster","TObjArray",&ioArray,32000,0);
198   }
199
200   if ((det >=                      0) && 
201       (det <  AliTRDgeometry::Ndet())) {
202
203     Int_t nRecPoints = RecPoints()->GetEntriesFast();
204     TObjArray *detRecPoints = new TObjArray(400);
205
206     for (Int_t i = 0; i < nRecPoints; i++) {
207       AliTRDcluster *c = (AliTRDcluster *) RecPoints()->UncheckedAt(i);
208       if (det == c->GetDetector()) {
209         detRecPoints->AddLast(c);
210       }
211       else {
212         AliError("Attempt to write a cluster with unexpected detector index\n");
213       }
214     }
215
216     branch->SetAddress(&detRecPoints);
217     fClusterTree->Fill();
218
219     delete detRecPoints;
220
221     return kTRUE;
222
223   }
224
225   if (det == -1) {
226
227     AliInfo(Form("Writing the cluster tree %s for event %d."
228                 ,fClusterTree->GetName(),fRunLoader->GetEventNumber()));
229
230     if (fRecPoints) {
231
232       branch->SetAddress(&fRecPoints);
233
234       AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
235       loader->WriteRecPoints("OVERWRITE");
236   
237     }
238     else {
239
240       AliError("Cluster tree does not exist. Cannot write clusters.\n");
241       return kFALSE;
242
243     }
244
245     return kTRUE;  
246
247   }
248
249   AliError(Form("Unexpected detector index %d.\n",det));
250  
251   return kFALSE;  
252   
253 }
254
255
256 //_____________________________________________________________________________
257 AliTRDcluster* AliTRDclusterizer::AddCluster(Double_t *pos, Int_t timebin
258                                            , Int_t det, Double_t amp
259                                            , Int_t *tracks, Double_t *sig
260                                            , Int_t iType, Float_t center)
261 {
262   //
263   // Add a cluster for the TRD
264   //
265
266   AliTRDcluster *c = new AliTRDcluster();
267
268   c->SetDetector(det);
269   c->SetQ(amp);
270   c->SetX(pos[2]);
271   c->SetY(pos[0]);
272   c->SetZ(pos[1]);
273   c->SetSigmaY2(sig[0]);   
274   c->SetSigmaZ2(sig[1]);
275   c->SetLocalTimeBin(timebin);
276   c->SetCenter(center);
277
278   if (tracks) {
279     c->AddTrackIndex(tracks);
280   }
281
282   switch (iType) {
283   case 0:
284     c->Set2pad();
285     break;
286   case 1:
287     c->Set3pad();
288     break;
289   case 2:
290     c->Set4pad();
291     break;
292   case 3:
293     c->Set5pad();
294     break;
295   case 4:
296     c->SetLarge();
297     break;
298   };
299
300   RecPoints()->Add(c);
301   return c;
302
303 }
304
305 //_____________________________________________________________________________
306 Double_t AliTRDclusterizer::CalcXposFromTimebin(Float_t timebin, Int_t idet
307                                               , Int_t col, Int_t row)
308 {
309   //
310   // Calculates the local x position in the detector from the timebin, 
311   // depends on the drift velocity and t0
312   //
313   
314   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
315   if (!calibration) {
316     AliError("Cannot find calibration object");
317     return -1;
318   }
319   AliTRDCommonParam *parCom      = AliTRDCommonParam::Instance();
320   if (!parCom) {
321     AliError("Could not get common parameters\n");
322     return kFALSE;
323   }
324
325   Float_t vdrift            = calibration->GetVdrift(idet,col,row);  
326   Float_t t0                = calibration->GetT0(idet,col,row);
327   Float_t samplingFrequency = parCom->GetSamplingFrequency();
328
329   timebin -= t0;
330
331   return timebin / samplingFrequency * vdrift;
332
333 }
334
335 //_____________________________________________________________________________
336 void AliTRDclusterizer::ResetRecPoints() 
337 {
338   //
339   // Resets the list of rec points
340   //
341
342   if (fRecPoints) {
343     fRecPoints->Delete();
344   }
345
346 }
347
348 //_____________________________________________________________________________
349 TObjArray* AliTRDclusterizer::RecPoints() 
350 {
351   //
352   // Returns the list of rec points
353   //
354
355   if (!fRecPoints) {
356     fRecPoints = new TObjArray(400);
357   }
358  
359   return fRecPoints;
360
361 }