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