]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDReconstructor.cxx
Important changes to the reconstructor classes. Complete elimination of the run-loade...
[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 #include "AliFMDGeometry.h"                // ALIFMDGEOMETRY_H
39 #include "AliFMDParameters.h"              // ALIFMDPARAMETERS_H
40 #include "AliFMDDigit.h"                   // ALIFMDDIGIT_H
41 #include "AliFMDReconstructor.h"           // ALIFMDRECONSTRUCTOR_H
42 #include "AliFMDRawReader.h"               // ALIFMDRAWREADER_H
43 #include "AliFMDRecPoint.h"                // ALIFMDMULTNAIIVE_H
44 #include "AliESDEvent.h"                   // ALIESDEVENT_H
45 #include "AliESDVertex.h"                  // ALIESDVERTEX_H
46 #include <AliESDFMD.h>                     // ALIESDFMD_H
47 #include <TMath.h>
48 #include <TH1.h>
49 #include <TH2.h>
50 #include <TFile.h>
51 class AliRawReader;
52
53 //____________________________________________________________________
54 ClassImp(AliFMDReconstructor)
55 #if 0
56   ; // This is here to keep Emacs for indenting the next line
57 #endif
58
59 //____________________________________________________________________
60 AliFMDReconstructor::AliFMDReconstructor() 
61   : AliReconstructor(),
62     fMult(0x0), 
63     fNMult(0),
64     fTreeR(0x0),
65     fCurrentVertex(0),
66     fESDObj(0x0),
67     fNoiseFactor(0),
68     fAngleCorrect(kTRUE),
69     fVertexType(kNoVertex),
70     fESD(0x0),
71     fDiagnostics(kFALSE),
72     fDiagStep1(0), 
73     fDiagStep2(0),
74     fDiagStep3(0),
75     fDiagStep4(0),
76     fDiagAll(0)
77 {
78   // Make a new FMD reconstructor object - default CTOR.  
79   SetNoiseFactor();
80   SetAngleCorrect();
81 }
82   
83
84 //____________________________________________________________________
85 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other) 
86   : AliReconstructor(), 
87     fMult(other.fMult),
88     fNMult(other.fNMult),
89     fTreeR(other.fTreeR),
90     fCurrentVertex(other.fCurrentVertex),
91     fESDObj(other.fESDObj),
92     fNoiseFactor(other.fNoiseFactor),
93     fAngleCorrect(other.fAngleCorrect),
94     fVertexType(other.fVertexType),
95     fESD(other.fESD),
96     fDiagnostics(other.fDiagnostics),
97     fDiagStep1(other.fDiagStep1), 
98     fDiagStep2(other.fDiagStep2),
99     fDiagStep3(other.fDiagStep3),
100     fDiagStep4(other.fDiagStep4),
101     fDiagAll(other.fDiagAll) 
102 {
103   // Copy constructor 
104 }
105   
106
107 //____________________________________________________________________
108 AliFMDReconstructor&
109 AliFMDReconstructor::operator=(const AliFMDReconstructor& other) 
110 {
111   // Assignment operator
112   fMult          = other.fMult;
113   fNMult         = other.fNMult;
114   fTreeR         = other.fTreeR;
115   fCurrentVertex = other.fCurrentVertex;
116   fESDObj        = other.fESDObj;
117   fNoiseFactor   = other.fNoiseFactor;
118   fAngleCorrect  = other.fAngleCorrect;
119   fVertexType    = other.fVertexType;
120   fESD           = other.fESD;
121   fDiagnostics   = other.fDiagnostics;
122   fDiagStep1     = other.fDiagStep1;
123   fDiagStep2     = other.fDiagStep2;
124   fDiagStep3     = other.fDiagStep3;
125   fDiagStep4     = other.fDiagStep4;
126   fDiagAll       = other.fDiagAll;
127   return *this;
128 }
129
130 //____________________________________________________________________
131 AliFMDReconstructor::~AliFMDReconstructor() 
132 {
133   // Destructor 
134   if (fMult)   fMult->Delete();
135   if (fMult)   delete fMult;
136   if (fESDObj) delete fESDObj;
137 }
138
139 //____________________________________________________________________
140 void 
141 AliFMDReconstructor::Init() 
142 {
143   // Initialize the reconstructor 
144
145   // Initialize the geometry 
146   AliFMDGeometry* geom = AliFMDGeometry::Instance();
147   geom->Init();
148   geom->InitTransformations();
149
150   // Initialize the parameters
151   AliFMDParameters* param = AliFMDParameters::Instance();
152   param->Init();
153   
154   // Current vertex position
155   fCurrentVertex = 0;
156   // Create array of reconstructed strip multiplicities 
157   fMult = new TClonesArray("AliFMDRecPoint", 51200);
158   // Create ESD output object 
159   fESDObj = new AliESDFMD;
160   
161   // Check if we need diagnostics histograms 
162   if (!fDiagnostics) return;
163   fDiagStep1   = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
164                         1024, -.5, 1023.5, 1024, -.5, 1023.5);
165   fDiagStep1->SetDirectory(0);
166   fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
167   fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)", 
168                                         fNoiseFactor));
169   fDiagStep2  = new TH2F("diagStep2",  "ADC vs Edep deduced",
170                         1024, -.5, 1023.5, 100, 0, 2);
171   fDiagStep2->SetDirectory(0);
172   fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
173   fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
174   fDiagStep3  = new TH2F("diagStep3",  "Edep vs Edep path corrected",
175                         100, 0., 2., 100, 0., 2.);
176   fDiagStep3->SetDirectory(0);
177   fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
178   fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
179   fDiagStep4  = new TH2F("diagStep4",  "Edep vs Multiplicity deduced", 
180                         100, 0., 2., 100, -.1, 19.9);
181   fDiagStep4->SetDirectory(0);
182   fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
183   fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
184   fDiagAll    = new TH2F("diagAll",    "Read ADC vs Multiplicity deduced", 
185                          1024, -.5, 1023.5, 100, -.1, 19.9);
186   fDiagAll->SetDirectory(0);
187   fDiagAll->GetXaxis()->SetTitle("ADC (read)");
188   fDiagAll->GetYaxis()->SetTitle("Multiplicity");
189 }
190
191 //____________________________________________________________________
192 void 
193 AliFMDReconstructor::ConvertDigits(AliRawReader* reader, 
194                                    TTree* digitsTree) const
195 {
196   // Convert Raw digits to AliFMDDigit's in a tree 
197   AliFMDDebug(2, ("Reading raw data into digits tree"));
198   AliFMDRawReader rawRead(reader, digitsTree);
199   // rawRead.SetSampleRate(fFMD->GetSampleRate());
200   rawRead.Exec();
201 }
202
203 //____________________________________________________________________
204 void 
205 AliFMDReconstructor::GetVertex() const
206 {
207   // Return the vertex to use. 
208   // This is obtained from the ESD object. 
209   // If not found, a warning is issued.
210   fVertexType    = kNoVertex;
211   fCurrentVertex = 0;
212   if (fESD) {
213     const AliESDVertex* vertex = fESD->GetVertex();
214     if (vertex) {
215       AliFMDDebug(2, ("Got vertex from ESD: %f", vertex->GetZv()));
216       fCurrentVertex = vertex->GetZv();
217       fVertexType    = kESDVertex;
218       return;
219     }
220   }
221   AliWarning("Didn't get any vertex from ESD or generator");
222 }
223   
224
225 //____________________________________________________________________
226 void 
227 AliFMDReconstructor::Reconstruct(TTree* digitsTree, 
228                                  TTree* clusterTree) const 
229 {
230   // Reconstruct event from digits in tree 
231   // Get the FMD branch holding the digits. 
232   // FIXME: The vertex may not be known yet, so we may have to move
233   // some of this to FillESD. 
234   AliFMDDebug(2, ("Reconstructing from digits in a tree"));
235   GetVertex();
236
237   TBranch *digitBranch = digitsTree->GetBranch("FMD");
238   if (!digitBranch) {
239     Error("Exec", "No digit branch for the FMD found");
240     return;
241   }
242   TClonesArray* digits = new TClonesArray("AliFMDDigit");
243   digitBranch->SetAddress(&digits);
244
245   if (fMult)   fMult->Clear();
246   if (fESDObj) fESDObj->Clear();
247   
248   fNMult = 0;
249   fTreeR = clusterTree;
250   fTreeR->Branch("FMD", &fMult);
251   
252   AliFMDDebug(5, ("Getting entry 0 from digit branch"));
253   digitBranch->GetEntry(0);
254   
255   AliFMDDebug(5, ("Processing digits"));
256   ProcessDigits(digits);
257
258   Int_t written = clusterTree->Fill();
259   AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
260   digits->Delete();
261   delete digits;
262 }
263  
264
265 //____________________________________________________________________
266 void
267 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
268 {
269   // For each digit, find the pseudo rapdity, azimuthal angle, and
270   // number of corrected ADC counts, and pass it on to the algorithms
271   // used. 
272   Int_t nDigits = digits->GetEntries();
273   AliFMDDebug(1, ("Got %d digits", nDigits));
274   fESDObj->SetNoiseFactor(fNoiseFactor);
275   fESDObj->SetAngleCorrected(fAngleCorrect);
276   for (Int_t i = 0; i < nDigits; i++) {
277     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
278     AliFMDParameters* param  = AliFMDParameters::Instance();
279     // Check that the strip is not marked as dead 
280     if (param->IsDead(digit->Detector(), digit->Ring(), 
281                       digit->Sector(), digit->Strip())) {
282       AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", digit->Detector(), 
283                         digit->Ring(), digit->Sector(), digit->Strip()));
284       continue;
285     }
286
287     // digit->Print();
288     // Get eta and phi 
289     Float_t eta, phi;
290     PhysicalCoordinates(digit, eta, phi);
291     
292     // Substract pedestal. 
293     UShort_t counts   = SubtractPedestal(digit);
294     
295     // Gain match digits. 
296     Double_t edep     = Adc2Energy(digit, eta, counts);
297     
298     // Make rough multiplicity 
299     Double_t mult     = Energy2Multiplicity(digit, edep);
300     
301     AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
302                       "ADC: %d, Counts: %d, Energy: %f, Mult: %f", 
303                       digit->Detector(), digit->Ring(), digit->Sector(),
304                       digit->Strip(), digit->Counts(), counts, edep, mult));
305     
306     // Create a `RecPoint' on the output branch. 
307     if (fMult) {
308       AliFMDRecPoint* m = 
309         new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(), 
310                                               digit->Ring(), 
311                                               digit->Sector(),
312                                               digit->Strip(),
313                                               eta, phi, 
314                                               edep, mult);
315       (void)m; // Suppress warnings about unused variables. 
316       fNMult++;
317     }
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   // This is most likely NOT the event number you'd like to use. It
478   // has nothing to do with the 'real' event number. 
479   // - That's OK.  We just use it for the name of the directory -
480   // nothing else.  Christian
481   Int_t evno = esd->GetEventNumberInFile(); 
482   AliFMDDebug(1, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
483   TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
484   first = false;
485   f.cd(); 
486   TDirectory* d = f.mkdir(Form("%03d", evno), 
487                           Form("Diagnostics histograms for event # %d", evno));
488   d->cd();
489   if (fDiagStep1) fDiagStep1->Write();
490   if (fDiagStep2) fDiagStep2->Write();
491   if (fDiagStep3) fDiagStep3->Write();
492   if (fDiagStep4) fDiagStep4->Write();
493   if (fDiagAll)   fDiagAll->Write();
494   d->Write();
495   f.Write();
496   f.Close();
497
498   if (fDiagStep1) fDiagStep1->Reset();
499   if (fDiagStep2) fDiagStep2->Reset();
500   if (fDiagStep3) fDiagStep3->Reset();
501   if (fDiagStep4) fDiagStep4->Reset();
502   if (fDiagAll)   fDiagAll->Reset();
503 }
504
505 //____________________________________________________________________
506 //
507 // EOF
508 //