]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROReconstructor.cxx
Major update of the reconstruction. All the code from the ANALYSIS/AliTriggerAnalysis...
[u/mrichter/AliRoot.git] / VZERO / AliVZEROReconstructor.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
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 ///                                                                          //
20 /// class for VZERO reconstruction                                           //
21 ///                                                                          //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <TH1F.h>
25 #include <TF1.h>
26
27 #include "AliRunLoader.h"
28 #include "AliRawReader.h"
29 #include "AliGRPObject.h"
30 #include "AliCDBManager.h"
31 #include "AliCDBStorage.h"
32 #include "AliCDBEntry.h"
33 #include "AliVZEROReconstructor.h"
34 #include "AliVZERORawStream.h"
35 #include "AliVZEROConst.h"
36 #include "AliESDEvent.h"
37 #include "AliVZEROTriggerMask.h"
38 #include "AliESDfriend.h"
39 #include "AliESDVZEROfriend.h"
40 #include "AliVZEROdigit.h"
41 #include "AliVZEROCalibData.h"
42 #include "AliRunInfo.h"
43 #include "AliCTPTimeParams.h"
44
45 ClassImp(AliVZEROReconstructor)
46
47 //_____________________________________________________________________________
48 AliVZEROReconstructor:: AliVZEROReconstructor(): AliReconstructor(),
49                         fESDVZERO(0x0),
50                         fESD(0x0),
51                         fESDVZEROfriend(0x0),
52                         fCalibData(NULL),
53                         fTimeSlewing(NULL),
54                         fCollisionMode(0),
55                         fBeamEnergy(0.),
56                         fDigitsArray(0)
57 {
58   // Default constructor  
59   // Get calibration data
60   
61   fCalibData = GetCalibData();
62
63   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
64   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
65   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
66   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
67
68   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
69   if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
70   TH1F *delays = (TH1F*)entry2->GetObject();
71
72   AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeSlewing");
73   if (!entry3) AliFatal("VZERO time slewing function is not found in OCDB !");
74   fTimeSlewing = (TF1*)entry3->GetObject();
75
76   for(Int_t i = 0 ; i < 64; ++i) {
77     Int_t board = AliVZEROCalibData::GetBoardNumber(i);
78     fTimeOffset[i] = (((Float_t)fCalibData->GetTriggerCountOffset(board)-
79                         (Float_t)fCalibData->GetRollOver(board))*25.0+
80                        fCalibData->GetTimeOffset(i)+
81                        l1Delay+
82                        delays->GetBinContent(i+1)+
83                        kV0Offset);
84   }
85 }
86
87
88 //_____________________________________________________________________________
89 AliVZEROReconstructor& AliVZEROReconstructor::operator = 
90   (const AliVZEROReconstructor& /*reconstructor*/)
91 {
92 // assignment operator
93
94   Fatal("operator =", "assignment operator not implemented");
95   return *this;
96 }
97
98 //_____________________________________________________________________________
99 AliVZEROReconstructor::~AliVZEROReconstructor()
100 {
101 // destructor
102
103    delete fESDVZERO;
104    delete fESDVZEROfriend;
105    delete fDigitsArray;
106 }
107
108 //_____________________________________________________________________________
109 void AliVZEROReconstructor::Init()
110 {
111 // initializer
112
113   fESDVZERO  = new AliESDVZERO;
114   fESDVZEROfriend = new AliESDVZEROfriend;
115   
116   GetCollisionMode();  // fCollisionMode =1 for Pb-Pb simulated data
117 }
118
119 //______________________________________________________________________
120 void AliVZEROReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
121 {
122 // converts RAW to digits 
123
124   if (!digitsTree) {
125     AliError("No digits tree!");
126     return;
127   }
128
129   if (!fDigitsArray)
130     fDigitsArray = new TClonesArray("AliVZEROdigit", 64);
131   digitsTree->Branch("VZERODigit", &fDigitsArray);
132
133   fESDVZEROfriend->Reset();
134
135   rawReader->Reset();
136   AliVZERORawStream rawStream(rawReader);
137   if (rawStream.Next()) { 
138      Float_t adc[64]; 
139      Float_t time[64], width[64];  
140      Bool_t BBFlag[64], BGFlag[64], integrator[64];
141      for(Int_t i=0; i<64; i++) {
142        Int_t j   =  rawStream.GetOfflineChannel(i);
143        adc[j]    = 0.0;
144        time[j]   = 0.0;
145        width[j]  = 0.0;
146        BBFlag[j] = kFALSE;
147        BGFlag[j] = kFALSE;
148        integrator[j] = kFALSE;
149        // Search for the maximum charge in the train of 21 LHC clocks 
150        // regardless of the integrator which has been operated:
151        Float_t maxadc = 0;
152        Int_t imax = -1;
153        Float_t adcPedSub[21];
154        for(Int_t iClock=0; iClock<21; iClock++){
155          Bool_t iIntegrator = rawStream.GetIntegratorFlag(i,iClock);
156          Int_t k = j+64*iIntegrator;
157          adcPedSub[iClock] = rawStream.GetPedestal(i,iClock) - fCalibData->GetPedestal(k);
158          if(adcPedSub[iClock] <= GetRecoParam()->GetNSigmaPed()*fCalibData->GetSigma(k)) {
159            adcPedSub[iClock] = 0;
160            continue;
161          }
162          if(iClock < GetRecoParam()->GetStartClock() || iClock > GetRecoParam()->GetEndClock()) continue;
163          if(adcPedSub[iClock] > maxadc) {
164            maxadc = adcPedSub[iClock];
165            imax   = iClock;
166          }
167        }
168
169        AliDebug(2,Form("Channel %d (online), %d (offline)",i,j)); 
170          if (imax != -1) {
171            Int_t start = imax - GetRecoParam()->GetNPreClocks();
172            if (start < 0) start = 0;
173            Int_t end = imax + GetRecoParam()->GetNPostClocks();
174            if (end > 20) end = 20;
175            for(Int_t iClock = start; iClock <= end; iClock++) {
176              if (iClock >= imax) {
177                BBFlag[j] |= rawStream.GetBBFlag(i,iClock);
178                BGFlag[j] |= rawStream.GetBGFlag(i,iClock);
179              }
180              if (iClock == imax)
181                adc[j] += rawStream.GetPedestal(i,iClock);
182              else 
183                adc[j] += adcPedSub[iClock];
184
185              AliDebug(2,Form("clock = %d adc = %f",iClock,rawStream.GetPedestal(i,iClock))); 
186            }
187            // Convert i (FEE channel numbering) to j (aliroot channel numbering)
188
189            integrator[j] =  rawStream.GetIntegratorFlag(i,imax); 
190          }
191
192          Int_t board   = j / 8;
193          time[j]       =  rawStream.GetTime(i)/ (25./256.) * fCalibData->GetTimeResolution(board);
194          width[j]      =  rawStream.GetWidth(i) / 0.4 * fCalibData->GetWidthResolution(board);
195
196          // Filling the esd friend object
197          fESDVZEROfriend->SetBBScalers(j,rawStream.GetBBScalers(i));
198          fESDVZEROfriend->SetBGScalers(j,rawStream.GetBGScalers(i));
199          for (Int_t iBunch = 0; iBunch < AliESDVZEROfriend::kNBunches; iBunch++) {
200              fESDVZEROfriend->SetChargeMB(j,iBunch,rawStream.GetChargeMB(i,iBunch));
201              fESDVZEROfriend->SetIntMBFlag(j,iBunch,rawStream.GetIntMBFlag(i,iBunch));
202              fESDVZEROfriend->SetBBMBFlag(j,iBunch,rawStream.GetBBMBFlag(i,iBunch));
203              fESDVZEROfriend->SetBGMBFlag(j,iBunch,rawStream.GetBGMBFlag(i,iBunch));
204          }
205          for (Int_t iEv = 0; iEv < AliESDVZEROfriend::kNEvOfInt; iEv++) {
206              fESDVZEROfriend->SetPedestal(j,iEv,rawStream.GetPedestal(i,iEv));
207              fESDVZEROfriend->SetIntegratorFlag(j,iEv,rawStream.GetIntegratorFlag(i,iEv));
208              fESDVZEROfriend->SetBBFlag(j,iEv,rawStream.GetBBFlag(i,iEv));
209              fESDVZEROfriend->SetBGFlag(j,iEv,rawStream.GetBGFlag(i,iEv));
210          }
211          fESDVZEROfriend->SetTime(j,time[j]);
212          fESDVZEROfriend->SetWidth(j,width[j]);
213      }  
214
215      // Filling the esd friend object
216      fESDVZEROfriend->SetTriggerInputs(rawStream.GetTriggerInputs());
217      fESDVZEROfriend->SetTriggerInputsMask(rawStream.GetTriggerInputsMask());
218
219      for(Int_t iScaler = 0; iScaler < AliESDVZEROfriend::kNScalers; iScaler++)
220          fESDVZEROfriend->SetTriggerScalers(iScaler,rawStream.GetTriggerScalers(iScaler));
221
222      for(Int_t iBunch = 0; iBunch < AliESDVZEROfriend::kNBunches; iBunch++)
223          fESDVZEROfriend->SetBunchNumbersMB(iBunch,rawStream.GetBunchNumbersMB(iBunch));
224      
225
226      // Channels(aliroot numbering) will be ordered in the tree
227      for(Int_t iChannel = 0; iChannel < 64; iChannel++) {
228          if(fCalibData->IsChannelDead(iChannel)){
229             adc[iChannel]  = (Float_t) kInvalidADC; 
230             time[iChannel] = (Float_t) kInvalidTime;     
231          }
232          if (adc[iChannel] > 0)
233            new ((*fDigitsArray)[fDigitsArray->GetEntriesFast()])
234              AliVZEROdigit(iChannel, adc[iChannel], time[iChannel],
235                            width[iChannel], BBFlag[iChannel], BGFlag[iChannel],integrator[iChannel]);
236         
237      }          
238      digitsTree->Fill();
239   }
240
241   fDigitsArray->Clear();
242 }
243
244 //______________________________________________________________________
245 void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
246                                     AliESDEvent* esd) const
247 {
248 // fills multiplicities to the ESD - pedestal is now subtracted
249     
250   if (!digitsTree) {
251       AliError("No digits tree!");
252       return;
253   }
254
255   TBranch* digitBranch = digitsTree->GetBranch("VZERODigit");
256   digitBranch->SetAddress(&fDigitsArray);
257
258   Float_t   mult[64];  
259   Float_t    adc[64]; 
260   Float_t   time[64]; 
261   Float_t  width[64];
262   Bool_t  BBFlag[64];
263   Bool_t  BGFlag[64];
264    
265   for (Int_t i=0; i<64; i++){
266        adc[i]    = 0.0;
267        mult[i]   = 0.0;
268        time[i]   = kInvalidTime;
269        width[i]  = 0.0;
270        BBFlag[i] = kFALSE;
271        BGFlag[i] = kFALSE;
272   }
273      
274   // loop over VZERO entries to get multiplicity
275   Int_t nEntries = (Int_t)digitsTree->GetEntries();
276   for (Int_t e=0; e<nEntries; e++) {
277     digitsTree->GetEvent(e);
278
279     Int_t nDigits = fDigitsArray->GetEntriesFast();
280     
281     for (Int_t d=0; d<nDigits; d++) {    
282         AliVZEROdigit* digit = (AliVZEROdigit*) fDigitsArray->At(d);      
283         Int_t  pmNumber      = digit->PMNumber(); 
284         // Pedestal retrieval and suppression: 
285         Bool_t   integrator  = digit->Integrator();
286         Int_t k = pmNumber+64*integrator;
287         Float_t  pedestal    = fCalibData->GetPedestal(k);
288         adc[pmNumber]   =  digit->ADC() - pedestal; 
289         time[pmNumber]  =  CorrectLeadingTime(pmNumber,digit->Time(),adc[pmNumber]);
290         width[pmNumber] =  digit->Width();
291         BBFlag[pmNumber]=  digit->BBFlag();
292         BGFlag[pmNumber]=  digit->BGFlag();
293
294         AliDebug(2,Form("PM = %d ADC = %f TDC %f",pmNumber, digit->ADC(),digit->Time()));
295
296         if(adc[pmNumber] > (fCalibData->GetPedestal(k) + GetRecoParam()->GetNSigmaPed()*fCalibData->GetSigma(k))) {
297             mult[pmNumber] += adc[pmNumber]*fCalibData->GetMIPperADC(pmNumber);
298         }           
299     } // end of loop over digits
300   } // end of loop over events in digits tree
301          
302   fESDVZERO->SetBit(AliESDVZERO::kCorrectedLeadingTime,kTRUE);
303   fESDVZERO->SetMultiplicity(mult);
304   fESDVZERO->SetADC(adc);
305   fESDVZERO->SetTime(time);
306   fESDVZERO->SetWidth(width);
307   fESDVZERO->SetBBFlag(BBFlag);
308   fESDVZERO->SetBGFlag(BGFlag);
309
310   // now fill the V0 decision and channel flags
311   {
312     AliVZEROTriggerMask triggerMask;
313     triggerMask.FillMasks(fESDVZERO, fCalibData, fTimeSlewing);
314   }
315
316   if (esd) { 
317      AliDebug(1, Form("Writing VZERO data to ESD tree"));
318      esd->SetVZEROData(fESDVZERO);
319   }
320
321   if (esd) {
322      AliESDfriend *fr = (AliESDfriend*)esd->FindListObject("AliESDfriend");
323      if (fr) {
324         AliDebug(1, Form("Writing VZERO friend data to ESD tree"));
325         fr->SetVZEROfriend(fESDVZEROfriend);
326     }
327   }
328
329   fDigitsArray->Clear();
330 }
331
332 //_____________________________________________________________________________
333 AliCDBStorage* AliVZEROReconstructor::SetStorage(const char *uri) 
334 {
335 // Sets the storage  
336
337   Bool_t deleteManager = kFALSE;
338   
339   AliCDBManager *manager = AliCDBManager::Instance();
340   AliCDBStorage *defstorage = manager->GetDefaultStorage();
341   
342   if(!defstorage || !(defstorage->Contains("VZERO"))){ 
343      AliWarning("No default storage set or default storage doesn't contain VZERO!");
344      manager->SetDefaultStorage(uri);
345      deleteManager = kTRUE;
346   }
347  
348   AliCDBStorage *storage = manager->GetDefaultStorage();
349
350   if(deleteManager){
351      AliCDBManager::Instance()->UnsetDefaultStorage();
352      defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
353   }
354
355   return storage; 
356 }
357
358 //____________________________________________________________________________
359 void AliVZEROReconstructor::GetCollisionMode()
360 {
361   // Retrieval of collision mode 
362
363   TString beamType = GetRunInfo()->GetBeamType();
364   if(beamType==AliGRPObject::GetInvalidString()){
365      AliError("VZERO cannot retrieve beam type");
366      return;
367   }
368
369   if( (beamType.CompareTo("P-P") ==0)  || (beamType.CompareTo("p-p") ==0) ){
370     fCollisionMode=0;
371   }
372   else if( (beamType.CompareTo("Pb-Pb") ==0)  || (beamType.CompareTo("A-A") ==0) ){
373     fCollisionMode=1;
374   }
375     
376   fBeamEnergy = GetRunInfo()->GetBeamEnergy();
377   if(fBeamEnergy==AliGRPObject::GetInvalidFloat()) {
378      AliError("Missing value for the beam energy ! Using 0");
379      fBeamEnergy = 0.;
380   }
381   
382   AliDebug(1,Form("\n ++++++ Beam type and collision mode retrieved as %s %d @ %1.3f GeV ++++++\n\n",beamType.Data(), fCollisionMode, fBeamEnergy));
383
384 }
385
386 //_____________________________________________________________________________
387 AliVZEROCalibData* AliVZEROReconstructor::GetCalibData() const
388 {
389   // Gets calibration object for VZERO set
390
391   AliCDBManager *man = AliCDBManager::Instance();
392
393   AliCDBEntry *entry=0;
394
395   entry = man->Get("VZERO/Calib/Data");
396
397   AliVZEROCalibData *calibdata = 0;
398
399   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
400   if (!calibdata)  AliFatal("No calibration data from calibration database !");
401
402   return calibdata;
403 }
404
405 Float_t AliVZEROReconstructor::CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const
406 {
407   // Correct the leading time
408   // for slewing effect and
409   // misalignment of the channels
410   if (time < 1e-6) return kInvalidTime;
411
412   // Channel alignment and general offset subtraction
413   if (i < 32) time -= kV0CDelayCables;
414   time -= fTimeOffset[i];
415
416   // In case of pathological signals
417   if (adc < 1e-6) return time;
418
419   // Slewing correction
420   Float_t thr = fCalibData->GetDiscriThr(i);
421   time -= fTimeSlewing->Eval(adc/thr);
422
423   return time;
424 }