]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDReconstructor.cxx
Correct raw data reconstruction in case of trigger
[u/mrichter/AliRoot.git] / FMD / AliFMDReconstructor.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 /* $Id$ */
16 /** @file    AliFMDReconstructor.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:47:09 2006
19     @brief   FMD reconstruction 
20 */
21 //____________________________________________________________________
22 //
23 // This is a class that constructs AliFMDMult (reconstructed
24 // multiplicity) from of Digits
25 //
26 // This class reads either digits from a TClonesArray or raw data from
27 // a DDL file (or similar), and stores the read ADC counts in an
28 // internal cache (fAdcs). 
29 //
30 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
31 //  Latest changes by Christian Holm Christensen <cholm@nbi.dk>
32 //
33 //
34 //____________________________________________________________________
35
36 #include <AliLog.h>                        // ALILOG_H
37 #include <AliRun.h>                        // ALIRUN_H
38 #include <AliRunLoader.h>                  // ALIRUNLOADER_H
39 #include <AliLoader.h>                     // ALILOADER_H
40 #include <AliHeader.h>                     // ALIHEADER_H
41 #include <AliRawReader.h>                  // ALIRAWREADER_H
42 #include <AliGenEventHeader.h>             // ALIGENEVENTHEADER_H
43 #include "AliFMD.h"                        // ALIFMD_H
44 #include "AliFMDGeometry.h"                // ALIFMDGEOMETRY_H
45 #include "AliFMDParameters.h"              // ALIFMDPARAMETERS_H
46 #include "AliFMDDetector.h"                // ALIFMDDETECTOR_H
47 #include "AliFMDRing.h"                    // ALIFMDRING_H
48 #include "AliFMDDigit.h"                   // ALIFMDDIGIT_H
49 #include "AliFMDReconstructor.h"           // ALIFMDRECONSTRUCTOR_H
50 #include "AliFMDRawStream.h"               // ALIFMDRAWSTREAM_H
51 #include "AliFMDRawReader.h"               // ALIFMDRAWREADER_H
52 #include "AliFMDRecPoint.h"                // ALIFMDMULTNAIIVE_H
53 #include "AliESD.h"                        // ALIESD_H
54 #include <AliESDFMD.h>                     // ALIESDFMD_H
55 #include <TFile.h>
56
57 //____________________________________________________________________
58 ClassImp(AliFMDReconstructor)
59 #if 0
60   ; // This is here to keep Emacs for indenting the next line
61 #endif
62
63 //____________________________________________________________________
64 AliFMDReconstructor::AliFMDReconstructor() 
65   : AliReconstructor(),
66     fMult(0x0), 
67     fNMult(0),
68     fTreeR(0x0),
69     fCurrentVertex(0),
70     fESDObj(0x0),
71     fESD(0x0)
72 {
73   // Make a new FMD reconstructor object - default CTOR.  
74 }
75   
76
77 //____________________________________________________________________
78 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other) 
79   : AliReconstructor(), 
80     fMult(other.fMult),
81     fNMult(other.fNMult),
82     fTreeR(other.fTreeR),
83     fCurrentVertex(other.fCurrentVertex),
84     fESDObj(other.fESDObj),
85     fESD(other.fESD)
86 {
87   // Copy constructor 
88 }
89   
90
91 //____________________________________________________________________
92 AliFMDReconstructor&
93 AliFMDReconstructor::operator=(const AliFMDReconstructor& other) 
94 {
95   // Assignment operator
96   fMult   = other.fMult;
97   fNMult = other.fNMult;
98   fTreeR = other.fTreeR;
99   fCurrentVertex = other.fCurrentVertex;
100   fESDObj = other.fESDObj;
101   fESD = other.fESD;
102   return *this;
103 }
104
105 //____________________________________________________________________
106 AliFMDReconstructor::~AliFMDReconstructor() 
107 {
108   // Destructor 
109   if (fMult)   fMult->Delete();
110   if (fMult)   delete fMult;
111   if (fESDObj) delete fESDObj;
112 }
113
114 //____________________________________________________________________
115 void 
116 AliFMDReconstructor::Init(AliRunLoader* runLoader) 
117 {
118   // Initialize the reconstructor 
119   AliDebug(1, Form("Init called with runloader 0x%x", runLoader));
120   // Initialize the geometry 
121   AliFMDGeometry* geom = AliFMDGeometry::Instance();
122   geom->Init();
123   geom->InitTransformations();
124
125   // Initialize the parameters
126   AliFMDParameters* param = AliFMDParameters::Instance();
127   param->Init();
128   
129   // Current vertex position
130   fCurrentVertex = 0;
131   // Create array of reconstructed strip multiplicities 
132   fMult = new TClonesArray("AliFMDRecPoint", 51200);
133   // Create ESD output object 
134   fESDObj = new AliESDFMD;
135   
136   // Check that we have a run loader
137   if (!runLoader) { 
138     Warning("Init", "No run loader");
139     return;
140   }
141
142   // Check if we can get the header tree 
143   AliHeader* header = runLoader->GetHeader();
144   if (!header) {
145     Warning("Init", "No header");
146     return;
147   }
148
149   // Check if we can get a simulation header 
150   AliGenEventHeader* eventHeader = header->GenEventHeader();
151   if (eventHeader) {
152     TArrayF vtx;
153     eventHeader->PrimaryVertex(vtx);
154     fCurrentVertex = vtx[2];
155     AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f", 
156                      header->GetRun(), header->GetEvent(), fCurrentVertex));
157     Warning("Init", "no generator event header");
158   }
159   else {
160     Warning("Init", "No generator event header - "
161             "perhaps we get the vertex from ESD?");
162   }
163 }
164
165 //____________________________________________________________________
166 void 
167 AliFMDReconstructor::ConvertDigits(AliRawReader* reader, 
168                                    TTree* digitsTree) const
169 {
170   // Convert Raw digits to AliFMDDigit's in a tree 
171   AliDebug(1, "Reading raw data into digits tree");
172   AliFMDRawReader rawRead(reader, digitsTree);
173   // rawRead.SetSampleRate(fFMD->GetSampleRate());
174   rawRead.Exec();
175 }
176
177 //____________________________________________________________________
178 void 
179 AliFMDReconstructor::Reconstruct(TTree* digitsTree, 
180                                  TTree* clusterTree) const 
181 {
182   // Reconstruct event from digits in tree 
183   // Get the FMD branch holding the digits. 
184   // FIXME: The vertex may not be known yet, so we may have to move
185   // some of this to FillESD. 
186   AliDebug(1, "Reconstructing from digits in a tree");
187 #if 1
188   if (fESD) {
189     const AliESDVertex* vertex = fESD->GetVertex();
190     if (vertex) {
191       AliDebug(1, Form("Got vertex from ESD: %f", vertex->GetZv()));
192       fCurrentVertex = vertex->GetZv();
193     }
194   }
195 #endif  
196   TBranch *digitBranch = digitsTree->GetBranch("FMD");
197   if (!digitBranch) {
198     Error("Exec", "No digit branch for the FMD found");
199     return;
200   }
201   TClonesArray* digits = new TClonesArray("AliFMDDigit");
202   digitBranch->SetAddress(&digits);
203
204   if (fMult)   fMult->Clear();
205   if (fESDObj) fESDObj->Clear();
206   
207   fNMult = 0;
208   fTreeR = clusterTree;
209   fTreeR->Branch("FMD", &fMult);
210   
211   AliDebug(5, "Getting entry 0 from digit branch");
212   digitBranch->GetEntry(0);
213   
214   AliDebug(5, "Processing digits");
215   ProcessDigits(digits);
216
217   Int_t written = clusterTree->Fill();
218   AliDebug(10, Form("Filled %d bytes into cluster tree", written));
219   digits->Delete();
220   delete digits;
221 }
222  
223
224 //____________________________________________________________________
225 void
226 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
227 {
228   // For each digit, find the pseudo rapdity, azimuthal angle, and
229   // number of corrected ADC counts, and pass it on to the algorithms
230   // used. 
231   Int_t nDigits = digits->GetEntries();
232   AliDebug(1, Form("Got %d digits", nDigits));
233   for (Int_t i = 0; i < nDigits; i++) {
234     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
235     AliFMDParameters* param  = AliFMDParameters::Instance();
236     // Check that the strip is not marked as dead 
237     if (param->IsDead(digit->Detector(), digit->Ring(), 
238                       digit->Sector(), digit->Strip())) continue;
239
240     // digit->Print();
241     // Get eta and phi 
242     Float_t eta, phi;
243     PhysicalCoordinates(digit, eta, phi);
244     
245     // Substract pedestal. 
246     UShort_t counts   = SubtractPedestal(digit);
247     
248     // Gain match digits. 
249     Double_t edep     = Adc2Energy(digit, eta, counts);
250     
251     // Make rough multiplicity 
252     Double_t mult     = Energy2Multiplicity(digit, edep);
253     
254     AliDebug(10, Form("FMD%d%c[%2d,%3d]: "
255                       "ADC: %d, Counts: %d, Energy: %f, Mult: %f", 
256                       digit->Detector(), digit->Ring(), digit->Sector(),
257                       digit->Strip(), digit->Counts(), counts, edep, mult));
258     
259     // Create a `RecPoint' on the output branch. 
260     AliFMDRecPoint* m = 
261       new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(), 
262                                             digit->Ring(), 
263                                             digit->Sector(),
264                                             digit->Strip(),
265                                             eta, phi, 
266                                             edep, mult);
267     (void)m; // Suppress warnings about unused variables. 
268     fNMult++;
269
270     fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(), 
271                              digit->Sector(),  digit->Strip(), mult);
272     fESDObj->SetEta(digit->Detector(), digit->Ring(), 
273                     digit->Sector(),  digit->Strip(), eta);
274   }
275 }
276
277 //____________________________________________________________________
278 UShort_t
279 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
280 {
281   // Member function to subtract the pedestal from a digit
282   // This implementation does nothing, but a derived class could over
283   // load this to subtract a pedestal that was given in a database or
284   // something like that. 
285
286   Int_t             counts = 0;
287   AliFMDParameters* param  = AliFMDParameters::Instance();
288   Float_t           pedM   = param->GetPedestal(digit->Detector(), 
289                                                 digit->Ring(), 
290                                                 digit->Sector(), 
291                                                 digit->Strip());
292   AliDebug(10, Form("Subtracting pedestal %f from signal %d", 
293                    pedM, digit->Counts()));
294   if (digit->Count3() > 0)      counts = digit->Count3();
295   else if (digit->Count2() > 0) counts = digit->Count2();
296   else                          counts = digit->Count1();
297   counts = TMath::Max(Int_t(counts - pedM), 0);
298   if (counts > 0) AliDebug(10, "Got a hit strip");
299   
300   return  UShort_t(counts);
301 }
302
303 //____________________________________________________________________
304 Float_t
305 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit, 
306                                 Float_t      /* eta */, 
307                                 UShort_t     count) const
308 {
309   // Converts number of ADC counts to energy deposited. 
310   // Note, that this member function can be overloaded by derived
311   // classes to do strip-specific look-ups in databases or the like,
312   // to find the proper gain for a strip. 
313   // 
314   // In this simple version, we calculate the energy deposited as 
315   // 
316   //    EnergyDeposited = cos(theta) * gain * count
317   // 
318   // where 
319   // 
320   //           Pre_amp_MIP_Range
321   //    gain = ----------------- * Energy_deposited_per_MIP
322   //           ADC_channel_size    
323   // 
324   // is constant and the same for all strips.
325
326   // Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
327   // Double_t edep  = TMath::Abs(TMath::Cos(theta)) * fGain * count;
328   AliFMDParameters* param = AliFMDParameters::Instance();
329   Float_t           gain  = param->GetPulseGain(digit->Detector(), 
330                                                 digit->Ring(), 
331                                                 digit->Sector(), 
332                                                 digit->Strip());
333   Double_t          edep  = count * gain;
334   AliDebug(10, Form("Converting counts %d to energy via factor %f", 
335                     count, gain));
336   return edep;
337 }
338
339 //____________________________________________________________________
340 Float_t
341 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */, 
342                                          Float_t      edep) const
343 {
344   // Converts an energy signal to number of particles. 
345   // Note, that this member function can be overloaded by derived
346   // classes to do strip-specific look-ups in databases or the like,
347   // to find the proper gain for a strip. 
348   // 
349   // In this simple version, we calculate the multiplicity as 
350   // 
351   //   multiplicity = Energy_deposited / Energy_deposited_per_MIP
352   // 
353   // where 
354   //
355   //   Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness 
356   // 
357   // is constant and the same for all strips 
358   AliFMDParameters* param   = AliFMDParameters::Instance();
359   Double_t          edepMIP = param->GetEdepMip();
360   Float_t           mult    = edep / edepMIP;
361   if (edep > 0) 
362     AliDebug(10, Form("Translating energy %f to multiplicity via "
363                      "divider %f->%f", edep, edepMIP, mult));
364   return mult;
365 }
366
367 //____________________________________________________________________
368 void
369 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit, 
370                                          Float_t& eta, 
371                                          Float_t& phi) const
372 {
373   // Get the eta and phi of a digit 
374   // 
375   // Get geometry. 
376   AliFMDGeometry* geom = AliFMDGeometry::Instance();
377   Double_t x, y, z, r, theta;
378   geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(), 
379                     digit->Strip(), x, y, z);
380   // Correct for vertex offset. 
381   z     += fCurrentVertex;
382   phi   =  TMath::ATan2(y, x);
383   r     =  TMath::Sqrt(y * y + x * x);
384   theta =  TMath::ATan2(r, z);
385   eta   = -TMath::Log(TMath::Tan(theta / 2));
386 }
387
388       
389
390 //____________________________________________________________________
391 void 
392 AliFMDReconstructor::FillESD(TTree*  /* digitsTree */, 
393                              TTree*  /* clusterTree */,
394                              AliESD* esd) const
395 {
396   // nothing to be done
397   // FIXME: The vertex may not be known when Reconstruct is executed,
398   // so we may have to move some of that member function here. 
399   AliDebug(1, Form("Calling FillESD with two trees and one ESD"));
400   // fESDObj->Print();
401
402   if (esd) { 
403     AliDebug(1, Form("Writing FMD data to ESD tree"));
404     esd->SetFMDData(fESDObj);
405   }
406 }
407
408
409 //____________________________________________________________________
410 void 
411 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const 
412 {
413   // Cannot be used.  See member function with same name but with 2
414   // TTree arguments.   Make sure you do local reconstrucion 
415   AliDebug(1, Form("Calling FillESD with loader and tree"));
416   AliError("MayNotUse");
417 }
418 //____________________________________________________________________
419 void 
420 AliFMDReconstructor::Reconstruct(AliRunLoader*) const 
421 {
422   // Cannot be used.  See member function with same name but with 2
423   // TTree arguments.   Make sure you do local reconstrucion 
424   AliDebug(1, Form("Calling FillESD with loader"));
425   AliError("MayNotUse");
426 }
427 //____________________________________________________________________
428 void 
429 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const 
430 {
431   // Cannot be used.  See member function with same name but with 2
432   // TTree arguments.   Make sure you do local reconstrucion 
433   AliDebug(1, Form("Calling FillESD with loader and raw reader"));
434   AliError("MayNotUse");
435 }
436 //____________________________________________________________________
437 void 
438 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESD*) const 
439 {
440   // Cannot be used.  See member function with same name but with 2
441   // TTree arguments.   Make sure you do local reconstrucion 
442   AliDebug(1, Form("Calling FillESD with raw reader, tree, and ESD"));
443   AliError("MayNotUse");
444 }
445 //____________________________________________________________________
446 void 
447 AliFMDReconstructor::FillESD(AliRunLoader*,AliESD*) const
448 {
449   // Cannot be used.  See member function with same name but with 2
450   // TTree arguments.   Make sure you do local reconstrucion 
451   AliDebug(1, Form("Calling FillESD with loader and ESD"));
452   AliError("MayNotUse");
453 }
454 //____________________________________________________________________
455 void 
456 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESD*) const 
457 {
458   // Cannot be used.  See member function with same name but with 2
459   // TTree arguments.   Make sure you do local reconstrucion 
460   AliDebug(1, Form("Calling FillESD with loader, raw reader, and ESD"));
461   AliError("MayNotUse");
462 }
463
464 //____________________________________________________________________
465 //
466 // EOF
467 //