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