]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCReconstructor.cxx
Access to the headers from RAW
[u/mrichter/AliRoot.git] / ZDC / AliZDCReconstructor.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 ZDC reconstruction                                              //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24
25 #include <TF1.h>
26
27 #include "AliRunLoader.h"
28 #include "AliRawReader.h"
29 #include "AliESD.h"
30 #include "AliZDCDigit.h"
31 #include "AliZDCRawStream.h"
32 #include "AliZDCReco.h"
33 #include "AliZDCReconstructor.h"
34 #include "AliZDCCalibData.h"
35
36
37 ClassImp(AliZDCReconstructor)
38
39
40 //_____________________________________________________________________________
41 AliZDCReconstructor:: AliZDCReconstructor()
42 {
43   // **** Default constructor
44   fStorage =  0;
45
46   //  ---      Number of generated spectator nucleons and impact parameter
47   // --------------------------------------------------------------------------------------------------
48   // [1] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
49   // Fit results for neutrons (Nspectator n true vs. EZN)
50   fZNCen = new TF1("fZNCen",
51       "(-2.287920+sqrt(2.287920*2.287920-4*(-0.007629)*(11.921710-x)))/(2*(-0.007629))",0.,164.);
52   fZNPer = new TF1("fZNPer",
53       "(-37.812280-sqrt(37.812280*37.812280-4*(-0.190932)*(-1709.249672-x)))/(2*(-0.190932))",0.,164.);
54   // Fit results for protons (Nspectator p true vs. EZP)
55   fZPCen = new TF1("fZPCen",
56        "(-1.321353+sqrt(1.321353*1.321353-4*(-0.007283)*(3.550697-x)))/(2*(-0.007283))",0.,60.);
57   fZPPer = new TF1("fZPPer",
58       "(-42.643308-sqrt(42.643308*42.643308-4*(-0.310786)*(-1402.945615-x)))/(2*(-0.310786))",0.,60.);
59   // Fit results for total number of spectators (Nspectators true vs. EZDC)
60   fZDCCen = new TF1("fZDCCen",
61       "(-1.934991+sqrt(1.934991*1.934991-4*(-0.004080)*(15.111124-x)))/(2*(-0.004080))",0.,225.);
62   fZDCPer = new TF1("fZDCPer",
63       "(-34.380639-sqrt(34.380639*34.380639-4*(-0.104251)*(-2612.189017-x)))/(2*(-0.104251))",0.,225.);
64   // --------------------------------------------------------------------------------------------------
65   // Fit results for b (b vs. EZDC)
66   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
67   fbCen = new TF1("fbCen","-0.056923+0.079703*x-0.0004301*x*x+0.000001366*x*x*x",0.,220.);
68   fbPer = new TF1("fbPer","17.943998-0.046846*x+0.000074*x*x",0.,220.);
69   // --------------------------------------------------------------------------------------------------
70   // Evaluating Nspectators and b from ZEM energy
71   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
72   fZEMn  = new TF1("fZEMn","126.2-0.05399*x+0.000005679*x*x",0.,4000.);
73   fZEMp  = new TF1("fZEMp","82.49-0.03611*x+0.00000385*x*x",0.,4000.);
74   fZEMsp = new TF1("fZEMsp","208.7-0.09006*x+0.000009526*x*x",0.,4000.);
75   fZEMb  = new TF1("fZEMb","16.06-0.01633*x+1.44e-5*x*x-6.778e-9*x*x*x+1.438e-12*x*x*x*x-1.112e-16*x*x*x*x*x",0.,4000.);
76   
77   // Setting storage
78   fStorage =  SetStorage("local://$ALICE_ROOT");
79
80   // Get calibration data
81   int runNumber = 0;
82   fCalibData = GetCalibData(runNumber); 
83 }
84
85 //_____________________________________________________________________________
86 AliZDCReconstructor::AliZDCReconstructor(const AliZDCReconstructor& 
87                                          reconstructor):
88   AliReconstructor(reconstructor)
89 {
90 // copy constructor
91
92   Fatal("AliZDCReconstructor", "copy constructor not implemented");
93 }
94
95 //_____________________________________________________________________________
96 AliZDCReconstructor& AliZDCReconstructor::operator = 
97   (const AliZDCReconstructor& /*reconstructor*/)
98 {
99 // assignment operator
100
101   Fatal("operator =", "assignment operator not implemented");
102   return *this;
103 }
104
105 //_____________________________________________________________________________
106 AliZDCReconstructor::~AliZDCReconstructor()
107 {
108 // destructor
109
110   delete fZNCen;
111   delete fZNPer;
112   delete fZPCen;
113   delete fZPPer;
114   delete fZDCCen;
115   delete fZDCPer;
116   delete fbCen;
117   delete fbPer;
118   delete fZEMn;
119   delete fZEMp;
120   delete fZEMsp;
121   delete fZEMb;
122 }
123
124
125 //_____________________________________________________________________________
126 void AliZDCReconstructor::Reconstruct(AliRunLoader* runLoader) const
127 {
128   // *** Local ZDC reconstruction for digits
129     
130   Float_t meanPed[47];
131   for(Int_t jj=0; jj<47; jj++) meanPed[jj] = fCalibData->GetMeanPed(jj);
132
133   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
134   if (!loader) return;
135   loader->LoadDigits("read");
136   loader->LoadRecPoints("recreate");
137   AliZDCDigit digit;
138   AliZDCDigit* pdigit = &digit;
139
140   // Event loop
141   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
142     runLoader->GetEvent(iEvent);
143
144     // load digits
145     loader->LoadDigits();
146     TTree* treeD = loader->TreeD();
147     if (!treeD) continue;
148     treeD->SetBranchAddress("ZDC", &pdigit);
149
150     // loop over digits
151     Float_t zncorr=0, zpcorr=0, zemcorr=0;
152     for (Int_t iDigit = 0; iDigit < treeD->GetEntries(); iDigit++) {
153       treeD->GetEntry(iDigit);
154       if (!pdigit) continue;
155
156       if(digit.GetSector(0) == 1)
157          zncorr  += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)]);     // ped 4 high gain ZN ADCs
158       else if(digit.GetSector(0) == 2)
159          zpcorr  += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+10]);  // ped 4 high gain ZP ADCs
160       else if(digit.GetSector(0) == 3){
161          if(digit.GetSector(1)==1)      
162            zemcorr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+20]); // ped 4 high gain ZEM1 ADCs
163          else if(digit.GetSector(1)==2) 
164            zemcorr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+22]); // ped 4 high gain ZEM2 ADCs
165       }
166     }
167     if(zncorr<0)  zncorr=0;
168     if(zpcorr<0)  zpcorr=0;
169     if(zemcorr<0) zemcorr=0;
170
171     // reconstruct the event
172     //printf("\n \t ZDCReco from digits-> Ev.#%d ZN = %.0f, ZP = %.0f, ZEM = %.0f\n",iEvent,zncorr,zpcorr,zemcorr);
173     ReconstructEvent(loader, zncorr, zpcorr, zemcorr);
174   }
175
176   loader->UnloadDigits();
177   loader->UnloadRecPoints();
178 }
179
180 //_____________________________________________________________________________
181 void AliZDCReconstructor::Reconstruct(AliRunLoader* runLoader, 
182                                       AliRawReader* rawReader) const
183 {
184   // *** Local ZDC reconstruction for raw data
185   
186   Float_t meanPed[47];
187   for(Int_t jj=0; jj<47; jj++) meanPed[jj] = fCalibData->GetMeanPed(jj);
188
189   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
190   if (!loader) return;
191   loader->LoadRecPoints("recreate");
192   // Event loop
193   Int_t iEvent = 0;
194   while (rawReader->NextEvent()) {
195     runLoader->GetEvent(iEvent++);
196
197     // loop over raw data digits
198     Float_t zncorr=0, zpcorr=0, zemcorr=0;
199     AliZDCRawStream digit(rawReader);
200     while (digit.Next()) {
201       if(digit.IsADCDataWord()){
202         if(digit.GetADCGain() == 0){
203           if(digit.GetSector(0) == 1)      zncorr  += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)]); // pedestals for high gain ZN ADCs;
204           else if(digit.GetSector(0) == 2) zpcorr  += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+10]); // pedestals for high gain ZP ADCs;
205           else if(digit.GetSector(0) == 3) zemcorr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+20]); // pedestals for high gain ZEM ADCs;
206         }
207       }
208     }
209     if(zncorr<0)  zncorr=0;
210     if(zpcorr<0)  zpcorr=0;
211     if(zemcorr<0) zemcorr=0;
212     
213     // reconstruct the event
214     //printf("\n\t ZDCReco from raw-> Ev.#%d ZN = %.0f, ZP = %.0f, ZEM = %.0f\n",iEvent,zncorr,zpcorr,zemcorr);
215     ReconstructEvent(loader, zncorr, zpcorr, zemcorr);
216   }
217
218   loader->UnloadRecPoints();
219 }
220
221 //_____________________________________________________________________________
222 void AliZDCReconstructor::ReconstructEvent(AliLoader* loader, Float_t zncorr,
223         Float_t zpcorr, Float_t zemcorr) const
224 {
225   // ***** Reconstruct one event
226   
227   //  ---      ADCchannel -> photoelectrons
228   // NB-> PM gain = 10^(5), ADC resolution = 6.4*10^(-7)
229   // Move to V965 (E.S.,15/09/04) NB-> PM gain = 10^(5), ADC resolution = 8*10^(-7)
230   Float_t znphe, zpphe, zemphe, convFactor = 0.08;
231   znphe  = zncorr/convFactor;
232   zpphe  = zpcorr/convFactor;
233   zemphe = zemcorr/convFactor;
234   //if AliDebug(1,Form("\n    znphe = %f, zpphe = %f, zemphe = %f\n",znphe, zpphe, zemphe);
235   
236   //  ---      Energy calibration
237   // Conversion factors for hadronic ZDCs goes from phe yield to TRUE 
238   // incident energy (conversion from GeV to TeV is included); while for EM 
239   // calos conversion is from light yield to detected energy calculated by
240   // GEANT NB -> ZN and ZP conversion factors are constant since incident
241   // spectators have all the same energy, ZEM energy is obtained through a
242   // fit over the whole range of incident particle energies 
243   // (obtained with full HIJING simulations) 
244   Float_t znenergy, zpenergy, zemenergy, zdcenergy;
245   Float_t znphexTeV=329., zpphexTeV=369.;
246   znenergy  = znphe/znphexTeV;
247   zpenergy  = zpphe/zpphexTeV;
248   zdcenergy = znenergy+zpenergy;
249   zemenergy = -4.81+0.3238*zemphe;
250   if(zemenergy<0) zemenergy=0;
251   //  if AliDebug(1,Form("    znenergy = %f TeV, zpenergy = %f TeV, zdcenergy = %f GeV, "
252   //                       "\n          zemenergy = %f TeV\n", znenergy, zpenergy, 
253   //                       zdcenergy, zemenergy);
254   //  if(zdcenergy==0)
255   //    if AliDebug(1,Form("\n\n        ###     ATTENZIONE!!! -> ev# %d: znenergy = %f TeV, zpenergy = %f TeV, zdcenergy = %f GeV, "
256   //                         " zemenergy = %f TeV\n\n", fMerger->EvNum(), znenergy, zpenergy, zdcenergy, zemenergy); 
257
258   //  ---      Number of incident spectator nucleons
259   Int_t nDetSpecN, nDetSpecP;
260   nDetSpecN = (Int_t) (znenergy/2.760);
261   nDetSpecP = (Int_t) (zpenergy/2.760);
262 //  if AliDebug(1,Form("\n    nDetSpecN = %d, nDetSpecP = %d\n",nDetSpecN, nDetSpecP);
263
264   Int_t nGenSpecN=0, nGenSpecP=0, nGenSpec=0;
265   Double_t impPar=0;
266   // Cut value for Ezem (GeV)
267   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
268   Float_t eZEMCut = 420.;
269   Float_t deltaEZEMSup = 690.; 
270   Float_t deltaEZEMInf = 270.; 
271   if(zemenergy > (eZEMCut+deltaEZEMSup)){
272     nGenSpecN = (Int_t) (fZNCen->Eval(znenergy));
273     nGenSpecP = (Int_t) (fZPCen->Eval(zpenergy));
274     nGenSpec  = (Int_t) (fZDCCen->Eval(zdcenergy));
275     impPar    = fbCen->Eval(zdcenergy);
276   }
277   else if(zemenergy < (eZEMCut-deltaEZEMInf)){
278     nGenSpecN = (Int_t) (fZNPer->Eval(znenergy)); 
279     nGenSpecP = (Int_t) (fZPPer->Eval(zpenergy));
280     nGenSpec  = (Int_t) (fZDCPer->Eval(zdcenergy));
281     impPar    = fbPer->Eval(zdcenergy);
282   }
283   else if(zemenergy >= (eZEMCut-deltaEZEMInf) && zemenergy <= (eZEMCut+deltaEZEMSup)){
284     nGenSpecN = (Int_t) (fZEMn->Eval(zemenergy));
285     nGenSpecP = (Int_t) (fZEMp->Eval(zemenergy));
286     nGenSpec  = (Int_t)(fZEMsp->Eval(zemenergy));
287     impPar    =  fZEMb->Eval(zemenergy);
288   }
289   // [2] ### Results from a new production  -> 0<b<18 fm (Apr 2002)
290   if(znenergy>162.)  nGenSpecN = (Int_t) (fZEMn->Eval(zemenergy));
291   if(zpenergy>59.75)  nGenSpecP = (Int_t) (fZEMp->Eval(zemenergy));
292   if(zdcenergy>221.5) nGenSpec  = (Int_t)(fZEMsp->Eval(zemenergy));
293   if(zdcenergy>220.)  impPar    =  fZEMb->Eval(zemenergy);
294   
295   if(nGenSpecN>125)    nGenSpecN=125;
296   else if(nGenSpecN<0) nGenSpecN=0;
297   if(nGenSpecP>82)     nGenSpecP=82;
298   else if(nGenSpecP<0) nGenSpecP=0;
299   if(nGenSpec>207)     nGenSpec=207;
300   else if(nGenSpec<0)  nGenSpec=0;
301   
302   //  ---      Number of participants
303   Int_t nPart, nPartTot;
304   nPart = 207-nGenSpecN-nGenSpecP;
305   nPartTot = 207-nGenSpec;
306   //printf("\t  ZDCeventReco-> ZNEn = %.0f GeV, ZPEn = %.0f GeV, ZEMEn = %.0f GeV\n",
307   //    znenergy, zpenergy, zemenergy);
308
309   // create the output tree
310   loader->MakeTree("R");
311   TTree* treeR = loader->TreeR();
312   AliZDCReco reco(znenergy, zpenergy, zdcenergy, zemenergy,
313                   nDetSpecN, nDetSpecP, nGenSpecN, nGenSpecP, nGenSpec,
314                   nPartTot, impPar);
315   AliZDCReco* preco = &reco;
316   const Int_t kBufferSize = 4000;
317   treeR->Branch("ZDC", "AliZDCReco", &preco, kBufferSize);
318
319   // write the output tree
320   treeR->Fill();
321   loader->WriteRecPoints("OVERWRITE");
322 }
323
324 //_____________________________________________________________________________
325 void AliZDCReconstructor::FillESD(AliRunLoader* runLoader, 
326                                   AliESD* esd) const
327 {
328 // fill energies and number of participants to the ESD
329
330   AliLoader* loader = runLoader->GetLoader("ZDCLoader");
331   if (!loader) return;
332   loader->LoadRecPoints();
333
334   TTree* treeR = loader->TreeR();
335   if (!treeR) return;
336   AliZDCReco reco;
337   AliZDCReco* preco = &reco;
338   treeR->SetBranchAddress("ZDC", &preco);
339
340   treeR->GetEntry(0);
341   esd->SetZDC(reco.GetZNenergy(), reco.GetZPenergy(), reco.GetZEMenergy(),
342               reco.GetNPart());
343
344   loader->UnloadRecPoints();
345 }
346
347 //_____________________________________________________________________________
348 AliCDBStorage* AliZDCReconstructor::SetStorage(const char *uri) 
349 {
350   //printf("\n\t AliZDCReconstructor::SetStorage \n");
351
352   Bool_t deleteManager = kFALSE;
353   
354   AliCDBManager *manager = AliCDBManager::Instance();
355   AliCDBStorage *defstorage = manager->GetDefaultStorage();
356   
357   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
358      AliWarning("No default storage set or default storage doesn't contain ZDC!");
359      manager->SetDefaultStorage(uri);
360      deleteManager = kTRUE;
361   }
362  
363   AliCDBStorage *storage = manager->GetDefaultStorage();
364
365   if(deleteManager){
366     AliCDBManager::Instance()->UnsetDefaultStorage();
367     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
368   }
369
370   return storage; 
371 }
372
373 //_____________________________________________________________________________
374 AliZDCCalibData* AliZDCReconstructor::GetCalibData(int runNumber) const
375 {
376
377   //printf("\n\t AliZDCReconstructor::GetCalibData \n");
378       
379   AliCDBEntry  *entry = fStorage->Get("ZDC/Calib/Data",runNumber);  
380   AliZDCCalibData *calibdata = (AliZDCCalibData*) entry->GetObject();
381     
382   if (!calibdata)  AliWarning("No calibration data from calibration database !");
383
384   return calibdata;
385 }