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