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