]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROReconstructor.cxx
Pedestal suppression implemented
[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 "AliVZEROReconstructor.h"
27 #include "AliVZERORawStream.h"
28 #include "AliESDEvent.h"
29
30 ClassImp(AliVZEROReconstructor)
31
32 //_____________________________________________________________________________
33 AliVZEROReconstructor:: AliVZEROReconstructor(): AliReconstructor(),
34    fESDVZERO(0x0),
35    fESD(0x0),
36    fCalibData(GetCalibData())
37 {
38   // Default constructor  
39   // Get calibration data
40   
41   // fCalibData = GetCalibData(); 
42 }
43
44
45 //_____________________________________________________________________________
46 AliVZEROReconstructor& AliVZEROReconstructor::operator = 
47   (const AliVZEROReconstructor& /*reconstructor*/)
48 {
49 // assignment operator
50
51   Fatal("operator =", "assignment operator not implemented");
52   return *this;
53 }
54
55 //_____________________________________________________________________________
56 AliVZEROReconstructor::~AliVZEROReconstructor()
57 {
58 // destructor
59    delete fESDVZERO; 
60    
61 }
62
63 //_____________________________________________________________________________
64 void AliVZEROReconstructor::Init()
65 {
66 // initializer
67
68   fESDVZERO  = new AliESDVZERO;
69 }
70
71 //______________________________________________________________________
72 void AliVZEROReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
73 {
74 // converts to digits
75
76 //  retrieval of calibration information 
77 //  Float_t     adc_pedestal[128],adc_sigma[128];
78 //  for(Int_t  i=0; i<128; i++){ adc_pedestal[i] = fCalibData->GetPedestal(i); 
79 //                               adc_sigma[i]    = fCalibData->GetSigma(i);                            
80 //                               printf(" i = %d pedestal = %f sigma = %f \n\n", 
81 //                                        i, adc_pedestal[i], adc_sigma[i] );} 
82               
83   if (!digitsTree) {
84     AliError("No digits tree!");
85     return;
86   }
87
88   TClonesArray* digitsArray = new TClonesArray("AliVZEROdigit");
89   digitsTree->Branch("VZERODigit", &digitsArray);
90
91   rawReader->Reset();
92   AliVZERORawStream rawStream(rawReader);
93   if (rawStream.Next()) {
94     for(Int_t iChannel = 0; iChannel < 64; iChannel++) {
95         Int_t adc;
96         if(!rawStream.GetIntegratorFlag(iChannel,10))
97            {adc = rawStream.GetADC(iChannel) - Int_t(fCalibData->GetPedestal(iChannel));} //even integrator
98         else 
99            {adc = rawStream.GetADC(iChannel) - Int_t(fCalibData->GetPedestal(iChannel+64));} // odd
100         Int_t time = rawStream.GetTime(iChannel);
101         if (adc >=0) new ((*digitsArray)[digitsArray->GetEntriesFast()])
102                      AliVZEROdigit(iChannel,adc,time);
103     }
104   }
105
106   digitsTree->Fill();
107 }
108
109 //______________________________________________________________________
110 void AliVZEROReconstructor::FillESD(TTree* digitsTree, TTree* /*clustersTree*/,
111                                     AliESDEvent* esd) const
112 {
113 // fills multiplicities to the ESD
114
115   if (!digitsTree) {
116     AliError("No digits tree!");
117     return;
118   }
119   
120   TClonesArray* digitsArray = NULL;
121   TBranch* digitBranch = digitsTree->GetBranch("VZERODigit");
122   digitBranch->SetAddress(&digitsArray);
123
124   Int_t   nbPMV0A = 0;
125   Int_t   nbPMV0C = 0;
126   Int_t   mTotV0A = 0;
127   Int_t   mTotV0C = 0;
128   Float_t adcV0A  = 0.0;
129   Float_t adcV0C  = 0.0;
130   Float_t multV0A[4];
131   Float_t multV0C[4];
132   Int_t   mRingV0A[4];
133   Int_t   mRingV0C[4];
134   
135   Int_t   adc[64]; 
136   Float_t mip[64];
137   for (Int_t i=0; i<64; i++){
138        adc[i] = 0;
139        mip[i] = 110.0;}
140   for (Int_t j=0; j<4; j++){
141        multV0A[j]  = 0.0;
142        multV0C[j]  = 0.0;
143        mRingV0A[j] = 0;
144        mRingV0C[j] = 0;}
145      
146   // loop over VZERO entries
147   Int_t nEntries = (Int_t)digitsTree->GetEntries();
148   for (Int_t e=0; e<nEntries; e++) {
149     digitsTree->GetEvent(e);
150
151     Int_t nDigits = digitsArray->GetEntriesFast();
152     
153     for (Int_t d=0; d<nDigits; d++) {    
154       AliVZEROdigit* digit = (AliVZEROdigit*)digitsArray->At(d);      
155       Int_t  pmNumber      = digit->PMNumber();  
156       adc[pmNumber] = digit->ADC(); 
157       // cut of ADC at MIP/2
158       if  (adc[pmNumber] > (mip[pmNumber]/2)) { 
159         if (pmNumber<=31) {
160           if (pmNumber<=7) multV0C[0]=multV0C[0]+ float(adc[pmNumber])/mip[pmNumber];
161           if (pmNumber>=8  && pmNumber<=15) multV0C[1]=multV0C[1]+ float(adc[pmNumber])/mip[pmNumber];
162           if (pmNumber>=16 && pmNumber<=23) multV0C[2]=multV0C[2]+ float(adc[pmNumber])/mip[pmNumber];
163           if (pmNumber>=24 && pmNumber<=31) multV0C[3]=multV0C[3]+ float(adc[pmNumber])/mip[pmNumber];
164           adcV0C = adcV0C + float(adc[pmNumber])/mip[pmNumber];
165           nbPMV0C++;
166         }       
167         if (pmNumber>=32 ) {
168           if (pmNumber>=32 && pmNumber<=39) multV0A[0]=multV0A[0]+ float(adc[pmNumber])/mip[pmNumber];
169           if (pmNumber>=40 && pmNumber<=47) multV0A[1]=multV0A[1]+ float(adc[pmNumber])/mip[pmNumber];
170           if (pmNumber>=48 && pmNumber<=55) multV0A[2]=multV0A[2]+ float(adc[pmNumber])/mip[pmNumber];
171           if (pmNumber>=56 && pmNumber<=63) multV0A[3]=multV0A[3]+ float(adc[pmNumber])/mip[pmNumber];
172           adcV0A = adcV0A + float(adc[pmNumber])/mip[pmNumber];
173           nbPMV0A++;
174         }
175       }
176     } // end of loop over digits
177     
178   } // end of loop over events in digits tree
179   
180   mTotV0A = int(adcV0A + 0.5);
181   mTotV0C = int(adcV0C + 0.5);
182   for (Int_t j=0; j<4; j++){       
183        mRingV0A[j] = int(multV0A[j] + 0.5);
184        mRingV0C[j] = int(multV0C[j] + 0.5);}
185      
186   AliDebug(1,Form("VZERO multiplicities : %d (V0A) %d (V0C)", mTotV0A, mTotV0C));
187   AliDebug(1,Form("Number of PMs fired  : %d (V0A) %d (V0C)", nbPMV0A, nbPMV0C));
188
189   fESDVZERO->SetNbPMV0A(nbPMV0A);
190   fESDVZERO->SetNbPMV0C(nbPMV0C);
191   fESDVZERO->SetMTotV0A(mTotV0A);
192   fESDVZERO->SetMTotV0C(mTotV0C);
193   fESDVZERO->SetMRingV0A(mRingV0A);
194   fESDVZERO->SetMRingV0C(mRingV0C);
195   
196   if (esd) { 
197       AliDebug(1, Form("Writing VZERO data to ESD tree"));
198       esd->SetVZEROData(fESDVZERO);
199   }
200 }
201
202 //_____________________________________________________________________________
203 AliCDBStorage* AliVZEROReconstructor::SetStorage(const char *uri) 
204 {
205 // Sets the storage  
206
207   Bool_t deleteManager = kFALSE;
208   
209   AliCDBManager *manager = AliCDBManager::Instance();
210   AliCDBStorage *defstorage = manager->GetDefaultStorage();
211   
212   if(!defstorage || !(defstorage->Contains("VZERO"))){ 
213      AliWarning("No default storage set or default storage doesn't contain VZERO!");
214      manager->SetDefaultStorage(uri);
215      deleteManager = kTRUE;
216   }
217  
218   AliCDBStorage *storage = manager->GetDefaultStorage();
219
220   if(deleteManager){
221     AliCDBManager::Instance()->UnsetDefaultStorage();
222     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
223   }
224
225   return storage; 
226 }
227
228 //_____________________________________________________________________________
229 AliVZEROCalibData* AliVZEROReconstructor::GetCalibData() const
230 {
231   // Gets calibration object for VZERO set
232
233   AliCDBManager *man = AliCDBManager::Instance();
234
235   AliCDBEntry *entry=0;
236
237   entry = man->Get("VZERO/Calib/Data");
238
239 //   if(!entry){
240 //     AliWarning("Load of calibration data from default storage failed!");
241 //     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
242 //     Int_t runNumber = man->GetRun();
243 //     entry = man->GetStorage("local://$ALICE_ROOT")
244 //       ->Get("VZERO/Calib/Data",runNumber);
245 //      
246 //   }
247
248   // Retrieval of data in directory VZERO/Calib/Data:
249
250   AliVZEROCalibData *calibdata = 0;
251
252   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
253   if (!calibdata)  AliFatal("No calibration data from calibration database !");
254
255   return calibdata;
256 }