]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDReconstructor.cxx
Fixes for SDigit generation. First attempt at making SDigit->Digit
[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   if (AliDebugLevel() > 0) fDiagnostics = kTRUE;
82 }
83   
84
85 //____________________________________________________________________
86 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other) 
87   : AliReconstructor(), 
88     fMult(other.fMult),
89     fNMult(other.fNMult),
90     fTreeR(other.fTreeR),
91     fCurrentVertex(other.fCurrentVertex),
92     fESDObj(other.fESDObj),
93     fNoiseFactor(other.fNoiseFactor),
94     fAngleCorrect(other.fAngleCorrect),
95     fVertexType(other.fVertexType),
96     fESD(other.fESD),
97     fDiagnostics(other.fDiagnostics),
98     fDiagStep1(other.fDiagStep1), 
99     fDiagStep2(other.fDiagStep2),
100     fDiagStep3(other.fDiagStep3),
101     fDiagStep4(other.fDiagStep4),
102     fDiagAll(other.fDiagAll) 
103 {
104   // Copy constructor 
105 }
106   
107
108 //____________________________________________________________________
109 AliFMDReconstructor&
110 AliFMDReconstructor::operator=(const AliFMDReconstructor& other) 
111 {
112   // Assignment operator
113   fMult          = other.fMult;
114   fNMult         = other.fNMult;
115   fTreeR         = other.fTreeR;
116   fCurrentVertex = other.fCurrentVertex;
117   fESDObj        = other.fESDObj;
118   fNoiseFactor   = other.fNoiseFactor;
119   fAngleCorrect  = other.fAngleCorrect;
120   fVertexType    = other.fVertexType;
121   fESD           = other.fESD;
122   fDiagnostics   = other.fDiagnostics;
123   fDiagStep1     = other.fDiagStep1;
124   fDiagStep2     = other.fDiagStep2;
125   fDiagStep3     = other.fDiagStep3;
126   fDiagStep4     = other.fDiagStep4;
127   fDiagAll       = other.fDiagAll;
128   return *this;
129 }
130
131 //____________________________________________________________________
132 AliFMDReconstructor::~AliFMDReconstructor() 
133 {
134   // Destructor 
135   if (fMult)   fMult->Delete();
136   if (fMult)   delete fMult;
137   if (fESDObj) delete fESDObj;
138 }
139
140 //____________________________________________________________________
141 void 
142 AliFMDReconstructor::Init() 
143 {
144   // Initialize the reconstructor 
145
146   // Initialize the geometry 
147   AliFMDGeometry* geom = AliFMDGeometry::Instance();
148   geom->Init();
149   geom->InitTransformations();
150
151   // Initialize the parameters
152   AliFMDParameters* param = AliFMDParameters::Instance();
153   param->Init();
154   
155   // Current vertex position
156   fCurrentVertex = 0;
157   // Create array of reconstructed strip multiplicities 
158   fMult = new TClonesArray("AliFMDRecPoint", 51200);
159   // Create ESD output object 
160   fESDObj = new AliESDFMD;
161   
162   // Check if we need diagnostics histograms 
163   if (!fDiagnostics) return;
164   AliInfo("Making diagnostics histograms");
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(AliRawReader* /*reader*/, TTree*) const
230 {
231   // Reconstruct directly from raw data (no intermediate output on
232   // digit tree or rec point tree).  
233   // Parameters: 
234   //   reader   Raw event reader 
235   //   ctree    Not used. 
236   AliError("Method is not used");
237 #if 0
238   TClonesArray*   array = new TClonesArray("AliFMDDigit");
239   AliFMDRawReader rawRead(reader, 0);
240   rawRead.ReadAdcs(array);
241   ProcessDigits(array);
242   array->Delete();
243   delete array;
244 #endif
245 }
246
247 //____________________________________________________________________
248 void 
249 AliFMDReconstructor::Reconstruct(TTree* digitsTree, 
250                                  TTree* clusterTree) const 
251 {
252   // Reconstruct event from digits in tree 
253   // Get the FMD branch holding the digits. 
254   // FIXME: The vertex may not be known yet, so we may have to move
255   // some of this to FillESD. 
256   AliFMDDebug(2, ("Reconstructing from digits in a tree"));
257   GetVertex();
258
259   TBranch *digitBranch = digitsTree->GetBranch("FMD");
260   if (!digitBranch) {
261     Error("Exec", "No digit branch for the FMD found");
262     return;
263   }
264   TClonesArray* digits = new TClonesArray("AliFMDDigit");
265   digitBranch->SetAddress(&digits);
266
267   if (fMult)   fMult->Clear();
268   if (fESDObj) fESDObj->Clear();
269   
270   fNMult = 0;
271   fTreeR = clusterTree;
272   fTreeR->Branch("FMD", &fMult);
273   
274   AliFMDDebug(5, ("Getting entry 0 from digit branch"));
275   digitBranch->GetEntry(0);
276   
277   AliFMDDebug(5, ("Processing digits"));
278   ProcessDigits(digits);
279
280   Int_t written = clusterTree->Fill();
281   AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
282   digits->Delete();
283   delete digits;
284 }
285  
286
287 //____________________________________________________________________
288 void
289 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
290 {
291   // For each digit, find the pseudo rapdity, azimuthal angle, and
292   // number of corrected ADC counts, and pass it on to the algorithms
293   // used. 
294   Int_t nDigits = digits->GetEntries();
295   AliFMDDebug(1, ("Got %d digits", nDigits));
296   fESDObj->SetNoiseFactor(fNoiseFactor);
297   fESDObj->SetAngleCorrected(fAngleCorrect);
298   for (Int_t i = 0; i < nDigits; i++) {
299     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
300     AliFMDParameters* param  = AliFMDParameters::Instance();
301     // Check that the strip is not marked as dead 
302     if (param->IsDead(digit->Detector(), digit->Ring(), 
303                       digit->Sector(), digit->Strip())) {
304       AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", digit->Detector(), 
305                         digit->Ring(), digit->Sector(), digit->Strip()));
306       continue;
307     }
308
309     // digit->Print();
310     // Get eta and phi 
311     Float_t eta, phi;
312     PhysicalCoordinates(digit, eta, phi);
313     
314     // Substract pedestal. 
315     UShort_t counts   = SubtractPedestal(digit);
316     
317     // Gain match digits. 
318     Double_t edep     = Adc2Energy(digit, eta, counts);
319     
320     // Make rough multiplicity 
321     Double_t mult     = Energy2Multiplicity(digit, edep);
322     
323     AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
324                       "ADC: %d, Counts: %d, Energy: %f, Mult: %f", 
325                       digit->Detector(), digit->Ring(), digit->Sector(),
326                       digit->Strip(), digit->Counts(), counts, edep, mult));
327     
328     // Create a `RecPoint' on the output branch. 
329     if (fMult) {
330       AliFMDRecPoint* m = 
331         new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(), 
332                                               digit->Ring(), 
333                                               digit->Sector(),
334                                               digit->Strip(),
335                                               eta, phi, 
336                                               edep, mult);
337       (void)m; // Suppress warnings about unused variables. 
338       fNMult++;
339     }
340     
341     fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(), 
342                              digit->Sector(),  digit->Strip(), mult);
343     fESDObj->SetEta(digit->Detector(), digit->Ring(), 
344                     digit->Sector(),  digit->Strip(), eta);
345
346     if (fDiagAll) fDiagAll->Fill(digit->Counts(), mult);  
347   }
348 }
349
350 //____________________________________________________________________
351 UShort_t
352 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
353 {
354   // Member function to subtract the pedestal from a digit
355   // This implementation does nothing, but a derived class could over
356   // load this to subtract a pedestal that was given in a database or
357   // something like that. 
358
359   AliFMDParameters* param  = AliFMDParameters::Instance();
360   Float_t           ped    = param->GetPedestal(digit->Detector(), 
361                                                 digit->Ring(), 
362                                                 digit->Sector(), 
363                                                 digit->Strip());
364   Float_t           noise  = param->GetPedestalWidth(digit->Detector(), 
365                                                      digit->Ring(), 
366                                                      digit->Sector(), 
367                                                      digit->Strip());
368   AliFMDDebug(15, ("Subtracting pedestal %f from signal %d", 
369                    ped, digit->Counts()));
370   // if (digit->Count3() > 0)      adc = digit->Count3();
371   // else if (digit->Count2() > 0) adc = digit->Count2();
372   // else                          adc = digit->Count1();
373   Int_t adc    = digit->Counts();
374   Int_t counts = TMath::Max(Int_t(adc - ped), 0);
375   if (counts < noise * fNoiseFactor) counts = 0;
376   if (counts > 0) AliFMDDebug(15, ("Got a hit strip"));
377   if (fDiagStep1) fDiagStep1->Fill(adc, counts);
378   
379   return  UShort_t(counts);
380 }
381
382 //____________________________________________________________________
383 Float_t
384 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit, 
385                                 Float_t      eta, 
386                                 UShort_t     count) const
387 {
388   // Converts number of ADC counts to energy deposited. 
389   // Note, that this member function can be overloaded by derived
390   // classes to do strip-specific look-ups in databases or the like,
391   // to find the proper gain for a strip. 
392   // 
393   // In this simple version, we calculate the energy deposited as 
394   // 
395   //    EnergyDeposited = cos(theta) * gain * count
396   // 
397   // where 
398   // 
399   //           Pre_amp_MIP_Range
400   //    gain = ----------------- * Energy_deposited_per_MIP
401   //           ADC_channel_size    
402   // 
403   // is constant and the same for all strips.
404   if (count <= 0) return 0;
405   AliFMDParameters* param = AliFMDParameters::Instance();
406   Float_t           gain  = param->GetPulseGain(digit->Detector(), 
407                                                 digit->Ring(), 
408                                                 digit->Sector(), 
409                                                 digit->Strip());
410   AliFMDDebug(15, ("Converting counts %d to energy via factor %f", 
411                     count, gain));
412
413   Double_t edep  = count * gain;
414   if (fDiagStep2) fDiagStep2->Fill(count, edep);  
415   if (fAngleCorrect) {
416     Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
417     Double_t corr  = TMath::Abs(TMath::Cos(theta));
418     Double_t cedep = corr * edep;
419     AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
420                       edep, corr, cedep, eta, theta));
421     if (fDiagStep3) fDiagStep3->Fill(edep, cedep);  
422     edep = cedep;
423   }
424   return edep;
425 }
426
427 //____________________________________________________________________
428 Float_t
429 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */, 
430                                          Float_t      edep) const
431 {
432   // Converts an energy signal to number of particles. 
433   // Note, that this member function can be overloaded by derived
434   // classes to do strip-specific look-ups in databases or the like,
435   // to find the proper gain for a strip. 
436   // 
437   // In this simple version, we calculate the multiplicity as 
438   // 
439   //   multiplicity = Energy_deposited / Energy_deposited_per_MIP
440   // 
441   // where 
442   //
443   //   Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness 
444   // 
445   // is constant and the same for all strips 
446   AliFMDParameters* param   = AliFMDParameters::Instance();
447   Double_t          edepMIP = param->GetEdepMip();
448   Float_t           mult    = edep / edepMIP;
449   if (edep > 0) 
450     AliFMDDebug(15, ("Translating energy %f to multiplicity via "
451                      "divider %f->%f", edep, edepMIP, mult));
452   if (fDiagStep4) fDiagStep4->Fill(edep, mult);  
453   return mult;
454 }
455
456 //____________________________________________________________________
457 void
458 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit, 
459                                          Float_t& eta, 
460                                          Float_t& phi) const
461 {
462   // Get the eta and phi of a digit 
463   // 
464   // Get geometry. 
465   AliFMDGeometry* geom = AliFMDGeometry::Instance();
466   Double_t x, y, z, r, theta;
467   geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(), 
468                     digit->Strip(), x, y, z);
469   // Correct for vertex offset. 
470   z     += fCurrentVertex;
471   phi   =  TMath::ATan2(y, x);
472   r     =  TMath::Sqrt(y * y + x * x);
473   theta =  TMath::ATan2(r, z);
474   eta   = -TMath::Log(TMath::Tan(theta / 2));
475 }
476
477       
478
479 //____________________________________________________________________
480 void 
481 AliFMDReconstructor::FillESD(TTree*  /* digitsTree */, 
482                              TTree*  /* clusterTree */,
483                              AliESDEvent* esd) const
484 {
485   // nothing to be done
486   // FIXME: The vertex may not be known when Reconstruct is executed,
487   // so we may have to move some of that member function here. 
488   AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
489   // fESDObj->Print();
490
491   if (esd) { 
492     AliFMDDebug(2, ("Writing FMD data to ESD tree"));
493     esd->SetFMDData(fESDObj);
494   }
495
496   if (!fDiagnostics || !esd) return;
497   static bool first = true;
498   // This is most likely NOT the event number you'd like to use. It
499   // has nothing to do with the 'real' event number. 
500   // - That's OK.  We just use it for the name of the directory -
501   // nothing else.  Christian
502   Int_t evno = esd->GetEventNumberInFile(); 
503   AliFMDDebug(1, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
504   TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
505   first = false;
506   f.cd(); 
507   TDirectory* d = f.mkdir(Form("%03d", evno), 
508                           Form("Diagnostics histograms for event # %d", evno));
509   d->cd();
510   if (fDiagStep1) fDiagStep1->Write();
511   if (fDiagStep2) fDiagStep2->Write();
512   if (fDiagStep3) fDiagStep3->Write();
513   if (fDiagStep4) fDiagStep4->Write();
514   if (fDiagAll)   fDiagAll->Write();
515   d->Write();
516   f.Write();
517   f.Close();
518
519   if (fDiagStep1) fDiagStep1->Reset();
520   if (fDiagStep2) fDiagStep2->Reset();
521   if (fDiagStep3) fDiagStep3->Reset();
522   if (fDiagStep4) fDiagStep4->Reset();
523   if (fDiagAll)   fDiagAll->Reset();
524 }
525
526 //____________________________________________________________________
527 void
528 AliFMDReconstructor::FillESD(AliRawReader*, TTree* clusterTree, 
529                              AliESDEvent* esd) const
530 {
531   TTree* dummy = 0;
532   FillESD(dummy, clusterTree, esd);
533 }
534
535 //____________________________________________________________________
536 //
537 // EOF
538 //