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