]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDReconstructor.cxx
Fixing compilation problem in debug build.
[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 AliFMDRecPoint objects from of Digits
24 // This class reads either digits from a TClonesArray or raw data from 
25 // a DDL file (or similar), and stores the read ADC counts in an
26 // internal cache (fAdcs).   The rec-points are made via the naiive
27 // method. 
28 //
29 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
30 //  Latest changes by Christian Holm Christensen <cholm@nbi.dk>
31 //
32 //
33 //____________________________________________________________________
34
35 // #include <AliLog.h>                        // ALILOG_H
36 // #include <AliRun.h>                        // ALIRUN_H
37 #include "AliFMDDebug.h"
38 // to be removed as soon as we remove it from the base class
39 #include "AliRunLoader.h"
40 #include "AliFMDGeometry.h"                // ALIFMDGEOMETRY_H
41 #include "AliFMDParameters.h"              // ALIFMDPARAMETERS_H
42 #include "AliFMDDigit.h"                   // ALIFMDDIGIT_H
43 #include "AliFMDReconstructor.h"           // ALIFMDRECONSTRUCTOR_H
44 #include "AliFMDRawReader.h"               // ALIFMDRAWREADER_H
45 #include "AliFMDRecPoint.h"                // ALIFMDMULTNAIIVE_H
46 #include "AliESDEvent.h"                   // ALIESDEVENT_H
47 #include "AliESDVertex.h"                  // ALIESDVERTEX_H
48 #include <AliESDFMD.h>                     // ALIESDFMD_H
49 #include <TMath.h>
50 #include <TH1.h>
51 #include <TH2.h>
52 #include <TFile.h>
53 class AliRawReader;
54
55 //____________________________________________________________________
56 ClassImp(AliFMDReconstructor)
57 #if 0
58   ; // This is here to keep Emacs for indenting the next line
59 #endif
60
61 //____________________________________________________________________
62 AliFMDReconstructor::AliFMDReconstructor() 
63   : AliReconstructor(),
64     fMult(0x0), 
65     fNMult(0),
66     fTreeR(0x0),
67     fCurrentVertex(0),
68     fESDObj(0x0),
69     fNoiseFactor(0),
70     fAngleCorrect(kTRUE),
71     fVertexType(kNoVertex),
72     fESD(0x0),
73     fDiagnostics(kFALSE),
74     fDiagStep1(0), 
75     fDiagStep2(0),
76     fDiagStep3(0),
77     fDiagStep4(0),
78     fDiagAll(0)
79 {
80   // Make a new FMD reconstructor object - default CTOR.  
81   SetNoiseFactor();
82   SetAngleCorrect();
83 }
84   
85
86 //____________________________________________________________________
87 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other) 
88   : AliReconstructor(), 
89     fMult(other.fMult),
90     fNMult(other.fNMult),
91     fTreeR(other.fTreeR),
92     fCurrentVertex(other.fCurrentVertex),
93     fESDObj(other.fESDObj),
94     fNoiseFactor(other.fNoiseFactor),
95     fAngleCorrect(other.fAngleCorrect),
96     fVertexType(other.fVertexType),
97     fESD(other.fESD),
98     fDiagnostics(other.fDiagnostics),
99     fDiagStep1(other.fDiagStep1), 
100     fDiagStep2(other.fDiagStep2),
101     fDiagStep3(other.fDiagStep3),
102     fDiagStep4(other.fDiagStep4),
103     fDiagAll(other.fDiagAll) 
104 {
105   // Copy constructor 
106 }
107   
108
109 //____________________________________________________________________
110 AliFMDReconstructor&
111 AliFMDReconstructor::operator=(const AliFMDReconstructor& other) 
112 {
113   // Assignment operator
114   fMult          = other.fMult;
115   fNMult         = other.fNMult;
116   fTreeR         = other.fTreeR;
117   fCurrentVertex = other.fCurrentVertex;
118   fESDObj        = other.fESDObj;
119   fNoiseFactor   = other.fNoiseFactor;
120   fAngleCorrect  = other.fAngleCorrect;
121   fVertexType    = other.fVertexType;
122   fESD           = other.fESD;
123   fDiagnostics   = other.fDiagnostics;
124   fDiagStep1     = other.fDiagStep1;
125   fDiagStep2     = other.fDiagStep2;
126   fDiagStep3     = other.fDiagStep3;
127   fDiagStep4     = other.fDiagStep4;
128   fDiagAll       = other.fDiagAll;
129   return *this;
130 }
131
132 //____________________________________________________________________
133 AliFMDReconstructor::~AliFMDReconstructor() 
134 {
135   // Destructor 
136   if (fMult)   fMult->Delete();
137   if (fMult)   delete fMult;
138   if (fESDObj) delete fESDObj;
139 }
140
141 //____________________________________________________________________
142 void 
143 AliFMDReconstructor::Init(AliRunLoader* /*runLoader*/) 
144 {
145   // Initialize the reconstructor 
146
147   // Initialize the geometry 
148   AliFMDGeometry* geom = AliFMDGeometry::Instance();
149   geom->Init();
150   geom->InitTransformations();
151
152   // Initialize the parameters
153   AliFMDParameters* param = AliFMDParameters::Instance();
154   param->Init();
155   
156   // Current vertex position
157   fCurrentVertex = 0;
158   // Create array of reconstructed strip multiplicities 
159   fMult = new TClonesArray("AliFMDRecPoint", 51200);
160   // Create ESD output object 
161   fESDObj = new AliESDFMD;
162   
163   // Check if we need diagnostics histograms 
164   if (!fDiagnostics) return;
165   fDiagStep1   = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
166                         1024, -.5, 1023.5, 1024, -.5, 1023.5);
167   fDiagStep1->SetDirectory(0);
168   fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
169   fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)", 
170                                         fNoiseFactor));
171   fDiagStep2  = new TH2F("diagStep2",  "ADC vs Edep deduced",
172                         1024, -.5, 1023.5, 100, 0, 2);
173   fDiagStep2->SetDirectory(0);
174   fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
175   fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
176   fDiagStep3  = new TH2F("diagStep3",  "Edep vs Edep path corrected",
177                         100, 0., 2., 100, 0., 2.);
178   fDiagStep3->SetDirectory(0);
179   fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
180   fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
181   fDiagStep4  = new TH2F("diagStep4",  "Edep vs Multiplicity deduced", 
182                         100, 0., 2., 100, -.1, 19.9);
183   fDiagStep4->SetDirectory(0);
184   fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
185   fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
186   fDiagAll    = new TH2F("diagAll",    "Read ADC vs Multiplicity deduced", 
187                          1024, -.5, 1023.5, 100, -.1, 19.9);
188   fDiagAll->SetDirectory(0);
189   fDiagAll->GetXaxis()->SetTitle("ADC (read)");
190   fDiagAll->GetYaxis()->SetTitle("Multiplicity");
191 }
192
193 //____________________________________________________________________
194 void 
195 AliFMDReconstructor::ConvertDigits(AliRawReader* reader, 
196                                    TTree* digitsTree) const
197 {
198   // Convert Raw digits to AliFMDDigit's in a tree 
199   AliFMDDebug(2, ("Reading raw data into digits tree"));
200   AliFMDRawReader rawRead(reader, digitsTree);
201   // rawRead.SetSampleRate(fFMD->GetSampleRate());
202   rawRead.Exec();
203 }
204
205 //____________________________________________________________________
206 void 
207 AliFMDReconstructor::GetVertex() const
208 {
209   // Return the vertex to use. 
210   // This is obtained from the ESD object. 
211   // If not found, a warning is issued.
212   fVertexType    = kNoVertex;
213   fCurrentVertex = 0;
214   if (fESD) {
215     const AliESDVertex* vertex = fESD->GetVertex();
216     if (vertex) {
217       AliFMDDebug(2, ("Got vertex from ESD: %f", vertex->GetZv()));
218       fCurrentVertex = vertex->GetZv();
219       fVertexType    = kESDVertex;
220       return;
221     }
222   }
223   AliWarning("Didn't get any vertex from ESD or generator");
224 }
225   
226
227 //____________________________________________________________________
228 void 
229 AliFMDReconstructor::Reconstruct(TTree* digitsTree, 
230                                  TTree* clusterTree) const 
231 {
232   // Reconstruct event from digits in tree 
233   // Get the FMD branch holding the digits. 
234   // FIXME: The vertex may not be known yet, so we may have to move
235   // some of this to FillESD. 
236   AliFMDDebug(2, ("Reconstructing from digits in a tree"));
237   GetVertex();
238
239   TBranch *digitBranch = digitsTree->GetBranch("FMD");
240   if (!digitBranch) {
241     Error("Exec", "No digit branch for the FMD found");
242     return;
243   }
244   TClonesArray* digits = new TClonesArray("AliFMDDigit");
245   digitBranch->SetAddress(&digits);
246
247   if (fMult)   fMult->Clear();
248   if (fESDObj) fESDObj->Clear();
249   
250   fNMult = 0;
251   fTreeR = clusterTree;
252   fTreeR->Branch("FMD", &fMult);
253   
254   AliFMDDebug(5, ("Getting entry 0 from digit branch"));
255   digitBranch->GetEntry(0);
256   
257   AliFMDDebug(5, ("Processing digits"));
258   ProcessDigits(digits);
259
260   Int_t written = clusterTree->Fill();
261   AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
262   digits->Delete();
263   delete digits;
264 }
265  
266
267 //____________________________________________________________________
268 void
269 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
270 {
271   // For each digit, find the pseudo rapdity, azimuthal angle, and
272   // number of corrected ADC counts, and pass it on to the algorithms
273   // used. 
274   Int_t nDigits = digits->GetEntries();
275   AliFMDDebug(1, ("Got %d digits", nDigits));
276   fESDObj->SetNoiseFactor(fNoiseFactor);
277   fESDObj->SetAngleCorrected(fAngleCorrect);
278   for (Int_t i = 0; i < nDigits; i++) {
279     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
280     AliFMDParameters* param  = AliFMDParameters::Instance();
281     // Check that the strip is not marked as dead 
282     if (param->IsDead(digit->Detector(), digit->Ring(), 
283                       digit->Sector(), digit->Strip())) {
284       AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", digit->Detector(), 
285                         digit->Ring(), digit->Sector(), digit->Strip()));
286       continue;
287     }
288
289     // digit->Print();
290     // Get eta and phi 
291     Float_t eta, phi;
292     PhysicalCoordinates(digit, eta, phi);
293     
294     // Substract pedestal. 
295     UShort_t counts   = SubtractPedestal(digit);
296     
297     // Gain match digits. 
298     Double_t edep     = Adc2Energy(digit, eta, counts);
299     
300     // Make rough multiplicity 
301     Double_t mult     = Energy2Multiplicity(digit, edep);
302     
303     AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
304                       "ADC: %d, Counts: %d, Energy: %f, Mult: %f", 
305                       digit->Detector(), digit->Ring(), digit->Sector(),
306                       digit->Strip(), digit->Counts(), counts, edep, mult));
307     
308     // Create a `RecPoint' on the output branch. 
309     AliFMDRecPoint* m = 
310       new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(), 
311                                             digit->Ring(), 
312                                             digit->Sector(),
313                                             digit->Strip(),
314                                             eta, phi, 
315                                             edep, mult);
316     (void)m; // Suppress warnings about unused variables. 
317     fNMult++;
318
319     fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(), 
320                              digit->Sector(),  digit->Strip(), mult);
321     fESDObj->SetEta(digit->Detector(), digit->Ring(), 
322                     digit->Sector(),  digit->Strip(), eta);
323
324     if (fDiagAll) fDiagAll->Fill(digit->Counts(), mult);  
325   }
326 }
327
328 //____________________________________________________________________
329 UShort_t
330 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
331 {
332   // Member function to subtract the pedestal from a digit
333   // This implementation does nothing, but a derived class could over
334   // load this to subtract a pedestal that was given in a database or
335   // something like that. 
336
337   Int_t             counts = 0;
338   Int_t             adc    = 0;
339   AliFMDParameters* param  = AliFMDParameters::Instance();
340   Float_t           ped    = param->GetPedestal(digit->Detector(), 
341                                                 digit->Ring(), 
342                                                 digit->Sector(), 
343                                                 digit->Strip());
344   Float_t           noise  = param->GetPedestalWidth(digit->Detector(), 
345                                                      digit->Ring(), 
346                                                      digit->Sector(), 
347                                                      digit->Strip());
348   AliFMDDebug(15, ("Subtracting pedestal %f from signal %d", 
349                    ped, digit->Counts()));
350   if (digit->Count3() > 0)      adc = digit->Count3();
351   else if (digit->Count2() > 0) adc = digit->Count2();
352   else                          adc = digit->Count1();
353   counts = TMath::Max(Int_t(adc - ped), 0);
354   if (counts < noise * fNoiseFactor) counts = 0;
355   if (counts > 0) AliFMDDebug(15, ("Got a hit strip"));
356   if (fDiagStep1) fDiagStep1->Fill(adc, counts);
357   
358   return  UShort_t(counts);
359 }
360
361 //____________________________________________________________________
362 Float_t
363 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit, 
364                                 Float_t      eta, 
365                                 UShort_t     count) const
366 {
367   // Converts number of ADC counts to energy deposited. 
368   // Note, that this member function can be overloaded by derived
369   // classes to do strip-specific look-ups in databases or the like,
370   // to find the proper gain for a strip. 
371   // 
372   // In this simple version, we calculate the energy deposited as 
373   // 
374   //    EnergyDeposited = cos(theta) * gain * count
375   // 
376   // where 
377   // 
378   //           Pre_amp_MIP_Range
379   //    gain = ----------------- * Energy_deposited_per_MIP
380   //           ADC_channel_size    
381   // 
382   // is constant and the same for all strips.
383   if (count <= 0) return 0;
384   AliFMDParameters* param = AliFMDParameters::Instance();
385   Float_t           gain  = param->GetPulseGain(digit->Detector(), 
386                                                 digit->Ring(), 
387                                                 digit->Sector(), 
388                                                 digit->Strip());
389   AliFMDDebug(15, ("Converting counts %d to energy via factor %f", 
390                     count, gain));
391
392   Double_t edep  = count * gain;
393   if (fDiagStep2) fDiagStep2->Fill(count, edep);  
394   if (fAngleCorrect) {
395     Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
396     Double_t corr  = TMath::Abs(TMath::Cos(theta));
397     Double_t cedep = corr * edep;
398     AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
399                       edep, corr, cedep, eta, theta));
400     if (fDiagStep3) fDiagStep3->Fill(edep, cedep);  
401     edep = cedep;
402   }
403   return edep;
404 }
405
406 //____________________________________________________________________
407 Float_t
408 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */, 
409                                          Float_t      edep) const
410 {
411   // Converts an energy signal to number of particles. 
412   // Note, that this member function can be overloaded by derived
413   // classes to do strip-specific look-ups in databases or the like,
414   // to find the proper gain for a strip. 
415   // 
416   // In this simple version, we calculate the multiplicity as 
417   // 
418   //   multiplicity = Energy_deposited / Energy_deposited_per_MIP
419   // 
420   // where 
421   //
422   //   Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness 
423   // 
424   // is constant and the same for all strips 
425   AliFMDParameters* param   = AliFMDParameters::Instance();
426   Double_t          edepMIP = param->GetEdepMip();
427   Float_t           mult    = edep / edepMIP;
428   if (edep > 0) 
429     AliFMDDebug(15, ("Translating energy %f to multiplicity via "
430                      "divider %f->%f", edep, edepMIP, mult));
431   if (fDiagStep4) fDiagStep4->Fill(edep, mult);  
432   return mult;
433 }
434
435 //____________________________________________________________________
436 void
437 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit, 
438                                          Float_t& eta, 
439                                          Float_t& phi) const
440 {
441   // Get the eta and phi of a digit 
442   // 
443   // Get geometry. 
444   AliFMDGeometry* geom = AliFMDGeometry::Instance();
445   Double_t x, y, z, r, theta;
446   geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(), 
447                     digit->Strip(), x, y, z);
448   // Correct for vertex offset. 
449   z     += fCurrentVertex;
450   phi   =  TMath::ATan2(y, x);
451   r     =  TMath::Sqrt(y * y + x * x);
452   theta =  TMath::ATan2(r, z);
453   eta   = -TMath::Log(TMath::Tan(theta / 2));
454 }
455
456       
457
458 //____________________________________________________________________
459 void 
460 AliFMDReconstructor::FillESD(TTree*  /* digitsTree */, 
461                              TTree*  /* clusterTree */,
462                              AliESDEvent* esd) const
463 {
464   // nothing to be done
465   // FIXME: The vertex may not be known when Reconstruct is executed,
466   // so we may have to move some of that member function here. 
467   AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
468   // fESDObj->Print();
469
470   if (esd) { 
471     AliFMDDebug(2, ("Writing FMD data to ESD tree"));
472     esd->SetFMDData(fESDObj);
473   }
474
475   if (!fDiagnostics || !esd) return;
476   static bool first = true;
477   Int_t evno = esd->GetEventNumberInFile(); // This is most likely NOT the event number you'd like to use. It has nothing to do with the 'real' event number.
478   AliFMDDebug(1, ("Writing diagnostics histograms to FMD.Diag.root/%03d",
479                    evno));
480   TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
481   first = false;
482   f.cd(); 
483   TDirectory* d = f.mkdir(Form("%03d", evno), 
484                           Form("Diagnostics histograms for event # %d", evno));
485   d->cd();
486   if (fDiagStep1) fDiagStep1->Write();
487   if (fDiagStep2) fDiagStep2->Write();
488   if (fDiagStep3) fDiagStep3->Write();
489   if (fDiagStep4) fDiagStep4->Write();
490   if (fDiagAll)   fDiagAll->Write();
491   d->Write();
492   f.Write();
493   f.Close();
494
495   if (fDiagStep1) fDiagStep1->Reset();
496   if (fDiagStep2) fDiagStep2->Reset();
497   if (fDiagStep3) fDiagStep3->Reset();
498   if (fDiagStep4) fDiagStep4->Reset();
499   if (fDiagAll)   fDiagAll->Reset();
500 }
501
502
503 //____________________________________________________________________
504 void 
505 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const 
506 {
507   // Cannot be used.  See member function with same name but with 2
508   // TTree arguments.   Make sure you do local reconstrucion 
509   AliFMDDebug(2, ("Calling FillESD with loader and tree"));
510   AliError("MayNotUse");
511 }
512 //____________________________________________________________________
513 void 
514 AliFMDReconstructor::Reconstruct(AliRunLoader*) const 
515 {
516   // Cannot be used.  See member function with same name but with 2
517   // TTree arguments.   Make sure you do local reconstrucion 
518   AliFMDDebug(2, ("Calling FillESD with loader"));
519   AliError("MayNotUse");
520 }
521 //____________________________________________________________________
522 void 
523 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const 
524 {
525   // Cannot be used.  See member function with same name but with 2
526   // TTree arguments.   Make sure you do local reconstrucion 
527   AliFMDDebug(2, ("Calling FillESD with loader and raw reader"));
528   AliError("MayNotUse");
529 }
530 //____________________________________________________________________
531 void 
532 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESDEvent*) const 
533 {
534   // Cannot be used.  See member function with same name but with 2
535   // TTree arguments.   Make sure you do local reconstrucion 
536   AliFMDDebug(2, ("Calling FillESD with raw reader, tree, and ESD"));
537   AliError("MayNotUse");
538 }
539 //____________________________________________________________________
540 void 
541 AliFMDReconstructor::FillESD(AliRunLoader*,AliESDEvent*) const
542 {
543   // Cannot be used.  See member function with same name but with 2
544   // TTree arguments.   Make sure you do local reconstrucion 
545   AliFMDDebug(2, ("Calling FillESD with loader and ESD"));
546   AliError("MayNotUse");
547 }
548 //____________________________________________________________________
549 void 
550 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESDEvent*) const 
551 {
552   // Cannot be used.  See member function with same name but with 2
553   // TTree arguments.   Make sure you do local reconstrucion 
554   AliFMDDebug(2, ("Calling FillESD with loader, raw reader, and ESD"));
555   AliError("MayNotUse");
556 }
557
558 //____________________________________________________________________
559 //
560 // EOF
561 //