]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROReconstructor.cxx
Introduction of the correction for signal-saturation in PbPb runs. The reco-param...
[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 #include <TParameter.h>
27
28 #include "AliRunLoader.h"
29 #include "AliRawReader.h"
30 #include "AliGRPObject.h"
31 #include "AliCDBManager.h"
32 #include "AliCDBStorage.h"
33 #include "AliCDBEntry.h"
34 #include "AliVZEROReconstructor.h"
35 #include "AliVZERORawStream.h"
36 #include "AliVZEROConst.h"
37 #include "AliESDEvent.h"
38 #include "AliVZEROTriggerMask.h"
39 #include "AliESDfriend.h"
40 #include "AliESDVZEROfriend.h"
41 #include "AliVZEROdigit.h"
42 #include "AliVZEROCalibData.h"
43 #include "AliRunInfo.h"
44 #include "AliCTPTimeParams.h"
45 #include "AliLHCClockPhase.h"
46
47 ClassImp(AliVZEROReconstructor)
48
49 //_____________________________________________________________________________
50 AliVZEROReconstructor:: AliVZEROReconstructor(): AliReconstructor(),
51                         fESDVZERO(0x0),
52                         fESD(0x0),
53                         fESDVZEROfriend(0x0),
54                         fCalibData(NULL),
55                         fTimeSlewing(NULL),
56                         fSaturationCorr(NULL),
57                         fCollisionMode(0),
58                         fBeamEnergy(0.),
59                         fDigitsArray(0)
60 {
61   // Default constructor  
62   // Get calibration data
63   
64   fCalibData = GetCalibData();
65
66   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
67   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
68   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
69   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
70
71   AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
72   if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
73   AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
74   l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
75
76   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
77   if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
78   TH1F *delays = (TH1F*)entry2->GetObject();
79
80   AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeSlewing");
81   if (!entry3) AliFatal("VZERO time slewing function is not found in OCDB !");
82   fTimeSlewing = (TF1*)entry3->GetObject();
83
84   AliCDBEntry *entry4 = AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
85   if (!entry4) AliFatal("LHC clock-phase shift is not found in OCDB !");
86   AliLHCClockPhase *phase = (AliLHCClockPhase*)entry4->GetObject();
87
88   for(Int_t i = 0 ; i < 64; ++i) {
89     Int_t board = AliVZEROCalibData::GetBoardNumber(i);
90     fTimeOffset[i] = (((Float_t)fCalibData->GetTriggerCountOffset(board)-
91                         (Float_t)fCalibData->GetRollOver(board))*25.0+
92                        fCalibData->GetTimeOffset(i)-
93                        l1Delay-
94                        phase->GetMeanPhase()+
95                        delays->GetBinContent(i+1)+
96                        kV0Offset);
97   }
98
99   AliCDBEntry *entry5 =  AliCDBManager::Instance()->Get("VZERO/Calib/Saturation");
100   if (!entry5) AliFatal("Saturation entry is not found in OCDB !");
101   fSaturationCorr = (TObjArray*)entry5->GetObject();
102 }
103
104
105 //_____________________________________________________________________________
106 AliVZEROReconstructor& AliVZEROReconstructor::operator = 
107   (const AliVZEROReconstructor& /*reconstructor*/)
108 {
109 // assignment operator
110
111   Fatal("operator =", "assignment operator not implemented");
112   return *this;
113 }
114
115 //_____________________________________________________________________________
116 AliVZEROReconstructor::~AliVZEROReconstructor()
117 {
118 // destructor
119
120   if(fESDVZERO)
121     delete fESDVZERO;
122   if(fESDVZEROfriend)
123    delete fESDVZEROfriend;
124   if(fDigitsArray)
125    delete fDigitsArray;
126 }
127
128 //_____________________________________________________________________________
129 void AliVZEROReconstructor::Init()
130 {
131 // initializer
132
133   fESDVZERO  = new AliESDVZERO;
134   fESDVZEROfriend = new AliESDVZEROfriend;
135   
136   GetCollisionMode();  // fCollisionMode =1 for Pb-Pb simulated data
137 }
138
139 //______________________________________________________________________
140 void AliVZEROReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
141 {
142 // converts RAW to digits 
143
144   if (!digitsTree) {
145     AliError("No digits tree!");
146     return;
147   }
148
149   if (!fDigitsArray)
150     fDigitsArray = new TClonesArray("AliVZEROdigit", 64);
151   digitsTree->Branch("VZERODigit", &fDigitsArray);
152
153   fESDVZEROfriend->Reset();
154
155   rawReader->Reset();
156   AliVZERORawStream rawStream(rawReader);
157   if (rawStream.Next()) { 
158
159     Int_t aBBflagsV0A = 0;
160     Int_t aBBflagsV0C = 0;
161     Int_t aBGflagsV0A = 0;
162     Int_t aBGflagsV0C = 0;
163
164     for(Int_t iChannel=0; iChannel < 64; ++iChannel) {
165       Int_t offlineCh = rawStream.GetOfflineChannel(iChannel);
166       // ADC charge samples
167       Short_t chargeADC[AliVZEROdigit::kNClocks];
168       for(Int_t iClock=0; iClock < AliVZEROdigit::kNClocks; ++iClock) {
169         chargeADC[iClock] = rawStream.GetPedestal(iChannel,iClock);
170       }
171       // Integrator flag
172       Bool_t integrator = rawStream.GetIntegratorFlag(iChannel,AliVZEROdigit::kNClocks/2);
173       // Beam-beam and beam-gas flags
174       if(offlineCh<32) {
175         if (rawStream.GetBBFlag(iChannel,AliVZEROdigit::kNClocks/2)) aBBflagsV0C |= (1 << offlineCh);
176         if (rawStream.GetBGFlag(iChannel,AliVZEROdigit::kNClocks/2)) aBGflagsV0C |= (1 << offlineCh);
177       } else {
178         if (rawStream.GetBBFlag(iChannel,AliVZEROdigit::kNClocks/2)) aBBflagsV0A |= (1 << (offlineCh-32));
179         if (rawStream.GetBGFlag(iChannel,AliVZEROdigit::kNClocks/2)) aBGflagsV0A |= (1 << (offlineCh-32));
180       }
181       // HPTDC data (leading time and width)
182       Int_t board = AliVZEROCalibData::GetBoardNumber(offlineCh);
183       Float_t time = rawStream.GetTime(iChannel)*fCalibData->GetTimeResolution(board);
184       Float_t width = rawStream.GetWidth(iChannel)*fCalibData->GetWidthResolution(board);
185       // Add a digit
186       if(!fCalibData->IsChannelDead(iChannel)){
187           new ((*fDigitsArray)[fDigitsArray->GetEntriesFast()])
188             AliVZEROdigit(offlineCh, time,
189                           width,integrator,
190                           chargeADC);
191       }
192
193       // Filling the part of esd friend object that is available only for raw data
194       fESDVZEROfriend->SetBBScalers(offlineCh,rawStream.GetBBScalers(iChannel));
195       fESDVZEROfriend->SetBGScalers(offlineCh,rawStream.GetBGScalers(iChannel));
196       for (Int_t iBunch = 0; iBunch < AliESDVZEROfriend::kNBunches; iBunch++) {
197         fESDVZEROfriend->SetChargeMB(offlineCh,iBunch,rawStream.GetChargeMB(iChannel,iBunch));
198         fESDVZEROfriend->SetIntMBFlag(offlineCh,iBunch,rawStream.GetIntMBFlag(iChannel,iBunch));
199         fESDVZEROfriend->SetBBMBFlag(offlineCh,iBunch,rawStream.GetBBMBFlag(iChannel,iBunch));
200         fESDVZEROfriend->SetBGMBFlag(offlineCh,iBunch,rawStream.GetBGMBFlag(iChannel,iBunch));
201       }
202       for (Int_t iEv = 0; iEv < AliESDVZEROfriend::kNEvOfInt; iEv++) {
203           fESDVZEROfriend->SetBBFlag(offlineCh,iEv,rawStream.GetBBFlag(iChannel,iEv));
204           fESDVZEROfriend->SetBGFlag(offlineCh,iEv,rawStream.GetBGFlag(iChannel,iEv));
205       }
206     }  
207
208     // Filling the global part of esd friend object that is available only for raw data
209     fESDVZEROfriend->SetTriggerInputs(rawStream.GetTriggerInputs());
210     fESDVZEROfriend->SetTriggerInputsMask(rawStream.GetTriggerInputsMask());
211
212     for(Int_t iScaler = 0; iScaler < AliESDVZEROfriend::kNScalers; iScaler++)
213       fESDVZEROfriend->SetTriggerScalers(iScaler,rawStream.GetTriggerScalers(iScaler));
214
215     for(Int_t iBunch = 0; iBunch < AliESDVZEROfriend::kNBunches; iBunch++)
216       fESDVZEROfriend->SetBunchNumbersMB(iBunch,rawStream.GetBunchNumbersMB(iBunch));
217      
218     // Store the BB and BG flags in the digits tree (user info)
219     digitsTree->GetUserInfo()->Add(new TParameter<int>("BBflagsV0A",aBBflagsV0A));
220     digitsTree->GetUserInfo()->Add(new TParameter<int>("BBflagsV0C",aBBflagsV0C));
221     digitsTree->GetUserInfo()->Add(new TParameter<int>("BGflagsV0A",aBGflagsV0A));
222     digitsTree->GetUserInfo()->Add(new TParameter<int>("BGflagsV0C",aBGflagsV0C));
223
224     digitsTree->Fill();
225   }
226
227   fDigitsArray->Clear();
228 }      
229
230 //______________________________________________________________________
231 void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
232                                     AliESDEvent* esd) const
233 {
234 // fills multiplicities to the ESD - pedestal is now subtracted
235     
236   if (!digitsTree) {
237       AliError("No digits tree!");
238       return;
239   }
240
241   TBranch* digitBranch = digitsTree->GetBranch("VZERODigit");
242   digitBranch->SetAddress(&fDigitsArray);
243
244   Float_t   mult[64];  
245   Float_t    adc[64]; 
246   Float_t   time[64]; 
247   Float_t  width[64];
248   Bool_t aBBflag[64];
249   Bool_t aBGflag[64];
250    
251   for (Int_t i=0; i<64; i++){
252        adc[i]    = 0.0;
253        mult[i]   = 0.0;
254        time[i]   = kInvalidTime;
255        width[i]  = 0.0;
256        aBBflag[i] = kFALSE;
257        aBGflag[i] = kFALSE;
258   }
259      
260   Int_t aBBflagsV0A = 0;
261   Int_t aBBflagsV0C = 0;
262   Int_t aBGflagsV0A = 0;
263   Int_t aBGflagsV0C = 0;
264
265   if (digitsTree->GetUserInfo()->FindObject("BBflagsV0A")) {
266     aBBflagsV0A = ((TParameter<int>*)digitsTree->GetUserInfo()->FindObject("BBflagsV0A"))->GetVal();
267   }
268   else
269     AliWarning("V0A beam-beam flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
270
271   if (digitsTree->GetUserInfo()->FindObject("BBflagsV0C")) {
272     aBBflagsV0C = ((TParameter<int>*)digitsTree->GetUserInfo()->FindObject("BBflagsV0C"))->GetVal();
273   }
274   else
275     AliWarning("V0C beam-beam flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
276
277   if (digitsTree->GetUserInfo()->FindObject("BGflagsV0A")) {
278     aBGflagsV0A = ((TParameter<int>*)digitsTree->GetUserInfo()->FindObject("BGflagsV0A"))->GetVal();
279   }
280   else
281     AliWarning("V0A beam-gas flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
282
283   if (digitsTree->GetUserInfo()->FindObject("BGflagsV0C")) {
284     aBGflagsV0C = ((TParameter<int>*)digitsTree->GetUserInfo()->FindObject("BGflagsV0C"))->GetVal();
285   }
286   else
287     AliWarning("V0C beam-gas flags not found in digits tree UserInfo! The flags will not be written to the raw-data stream!");
288   
289   // Beam-beam and beam-gas flags (hardware)
290   for (Int_t iChannel = 0; iChannel < 64; ++iChannel) {
291     if(iChannel < 32) {
292       aBBflag[iChannel] = (aBBflagsV0C >> iChannel) & 0x1;
293       aBGflag[iChannel] = (aBGflagsV0C >> iChannel) & 0x1;
294     }
295     else {
296       aBBflag[iChannel] = (aBBflagsV0A >> (iChannel-32)) & 0x1;
297       aBGflag[iChannel] = (aBGflagsV0A >> (iChannel-32)) & 0x1;
298     }
299   }
300
301   Int_t nEntries = (Int_t)digitsTree->GetEntries();
302   for (Int_t e=0; e<nEntries; e++) {
303     digitsTree->GetEvent(e);
304
305     Int_t nDigits = fDigitsArray->GetEntriesFast();
306     
307     for (Int_t d=0; d<nDigits; d++) {    
308         AliVZEROdigit* digit = (AliVZEROdigit*) fDigitsArray->At(d);      
309         Int_t  pmNumber = digit->PMNumber();
310
311         // Pedestal retrieval and suppression
312         Bool_t integrator = digit->Integrator();
313         Float_t maxadc = 0;
314         Int_t imax = -1;
315         Float_t adcPedSub[AliVZEROdigit::kNClocks];
316         for(Int_t iClock=0; iClock < AliVZEROdigit::kNClocks; ++iClock) {
317           Short_t charge = digit->ChargeADC(iClock);
318           Bool_t iIntegrator = (iClock%2 == 0) ? integrator : !integrator;
319           Int_t k = pmNumber + 64*iIntegrator;
320           adcPedSub[iClock] = (Float_t)charge - fCalibData->GetPedestal(k);
321           if(adcPedSub[iClock] <= GetRecoParam()->GetNSigmaPed()*fCalibData->GetSigma(k)) {
322             adcPedSub[iClock] = 0;
323             continue;
324           }
325           if(iClock < GetRecoParam()->GetStartClock() || iClock > GetRecoParam()->GetEndClock()) continue;
326           if(adcPedSub[iClock] > maxadc) {
327             maxadc = adcPedSub[iClock];
328             imax   = iClock;
329           }
330         }
331
332         if (imax != -1) {
333           Int_t start = imax - GetRecoParam()->GetNPreClocks();
334           if (start < 0) start = 0;
335           Int_t end = imax + GetRecoParam()->GetNPostClocks();
336           if (end > 20) end = 20;
337           for(Int_t iClock = start; iClock <= end; iClock++) {
338             adc[pmNumber] += adcPedSub[iClock];
339           }
340         }
341
342         // HPTDC leading time and width
343         // Correction for slewing and various time delays
344         time[pmNumber]  =  CorrectLeadingTime(pmNumber,digit->Time(),adc[pmNumber]);
345         width[pmNumber] =  digit->Width();
346
347         if (adc[pmNumber] > 0) {
348           AliDebug(1,Form("PM = %d ADC = %f TDC %f (%f)   Int %d (%d %d %d %d %d)    %f %f   %f %f    %d %d",pmNumber, adc[pmNumber],
349                           digit->Time(),time[pmNumber],
350                           integrator,
351                           digit->ChargeADC(8),digit->ChargeADC(9),digit->ChargeADC(10),
352                           digit->ChargeADC(11),digit->ChargeADC(12),
353                           fCalibData->GetPedestal(pmNumber),fCalibData->GetSigma(pmNumber),
354                           fCalibData->GetPedestal(pmNumber+64),fCalibData->GetSigma(pmNumber+64),
355                           aBBflag[pmNumber],aBGflag[pmNumber]));
356             };
357
358         TF1 *saturationFunc = (TF1*)fSaturationCorr->UncheckedAt(pmNumber);
359         if (!saturationFunc) AliFatal(Form("Saturation correction for channel %d is not found!",pmNumber));
360         AliDebug(1,Form("Saturation PM=%d   %f %f",pmNumber,adc[pmNumber],saturationFunc->Eval(adc[pmNumber])));
361         mult[pmNumber] = saturationFunc->Eval(adc[pmNumber])*fCalibData->GetMIPperADC(pmNumber);
362
363         // Fill ESD friend object
364         for (Int_t iEv = 0; iEv < AliESDVZEROfriend::kNEvOfInt; iEv++) {
365           fESDVZEROfriend->SetPedestal(pmNumber,iEv,(Float_t)digit->ChargeADC(iEv));
366           fESDVZEROfriend->SetIntegratorFlag(pmNumber,iEv,(iEv%2 == 0) ? integrator : !integrator);
367         }
368         fESDVZEROfriend->SetTime(pmNumber,digit->Time());
369         fESDVZEROfriend->SetWidth(pmNumber,digit->Width());
370
371     } // end of loop over digits
372   } // end of loop over events in digits tree
373          
374   fESDVZERO->SetBit(AliESDVZERO::kCorrectedLeadingTime,kTRUE);
375   fESDVZERO->SetMultiplicity(mult);
376   fESDVZERO->SetADC(adc);
377   fESDVZERO->SetTime(time);
378   fESDVZERO->SetWidth(width);
379   fESDVZERO->SetBit(AliESDVZERO::kOnlineBitsFilled,kTRUE);
380   fESDVZERO->SetBBFlag(aBBflag);
381   fESDVZERO->SetBGFlag(aBGflag);
382
383   // now fill the V0 decision and channel flags
384   {
385     AliVZEROTriggerMask triggerMask;
386     triggerMask.FillMasks(fESDVZERO, fCalibData, fTimeSlewing);
387   }
388
389   if (esd) { 
390      AliDebug(1, Form("Writing VZERO data to ESD tree"));
391      esd->SetVZEROData(fESDVZERO);
392   }
393
394   if (esd) {
395      AliESDfriend *fr = (AliESDfriend*)esd->FindListObject("AliESDfriend");
396      if (fr) {
397         AliDebug(1, Form("Writing VZERO friend data to ESD tree"));
398         fr->SetVZEROfriend(fESDVZEROfriend);
399     }
400   }
401
402   fDigitsArray->Clear();
403 }
404
405 //_____________________________________________________________________________
406 AliCDBStorage* AliVZEROReconstructor::SetStorage(const char *uri) 
407 {
408 // Sets the storage  
409
410   Bool_t deleteManager = kFALSE;
411   
412   AliCDBManager *manager = AliCDBManager::Instance();
413   AliCDBStorage *defstorage = manager->GetDefaultStorage();
414   
415   if(!defstorage || !(defstorage->Contains("VZERO"))){ 
416      AliWarning("No default storage set or default storage doesn't contain VZERO!");
417      manager->SetDefaultStorage(uri);
418      deleteManager = kTRUE;
419   }
420  
421   AliCDBStorage *storage = manager->GetDefaultStorage();
422
423   if(deleteManager){
424      AliCDBManager::Instance()->UnsetDefaultStorage();
425      defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
426   }
427
428   return storage; 
429 }
430
431 //____________________________________________________________________________
432 void AliVZEROReconstructor::GetCollisionMode()
433 {
434   // Retrieval of collision mode 
435
436   TString beamType = GetRunInfo()->GetBeamType();
437   if(beamType==AliGRPObject::GetInvalidString()){
438      AliError("VZERO cannot retrieve beam type");
439      return;
440   }
441
442   if( (beamType.CompareTo("P-P") ==0)  || (beamType.CompareTo("p-p") ==0) ){
443     fCollisionMode=0;
444   }
445   else if( (beamType.CompareTo("Pb-Pb") ==0)  || (beamType.CompareTo("A-A") ==0) ){
446     fCollisionMode=1;
447   }
448     
449   fBeamEnergy = GetRunInfo()->GetBeamEnergy();
450   if(fBeamEnergy==AliGRPObject::GetInvalidFloat()) {
451      AliError("Missing value for the beam energy ! Using 0");
452      fBeamEnergy = 0.;
453   }
454   
455   AliDebug(1,Form("\n ++++++ Beam type and collision mode retrieved as %s %d @ %1.3f GeV ++++++\n\n",beamType.Data(), fCollisionMode, fBeamEnergy));
456
457 }
458
459 //_____________________________________________________________________________
460 AliVZEROCalibData* AliVZEROReconstructor::GetCalibData() const
461 {
462   // Gets calibration object for VZERO set
463
464   AliCDBManager *man = AliCDBManager::Instance();
465
466   AliCDBEntry *entry=0;
467
468   entry = man->Get("VZERO/Calib/Data");
469
470   AliVZEROCalibData *calibdata = 0;
471
472   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
473   if (!calibdata)  AliFatal("No calibration data from calibration database !");
474
475   return calibdata;
476 }
477
478 Float_t AliVZEROReconstructor::CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const
479 {
480   // Correct the leading time
481   // for slewing effect and
482   // misalignment of the channels
483   if (time < 1e-6) return kInvalidTime;
484
485   // Channel alignment and general offset subtraction
486   if (i < 32) time -= kV0CDelayCables;
487   time -= fTimeOffset[i];
488
489   // In case of pathological signals
490   if (adc < 1e-6) return time;
491
492   // Slewing correction
493   Float_t thr = fCalibData->GetCalibDiscriThr(i,kTRUE);
494   time -= fTimeSlewing->Eval(adc/thr);
495
496   return time;
497 }