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