]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCReconstructor.cxx
9ede80ceb3fd0826386c890996a11f0888c196c5
[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 //                  Author: Chiara.Oppedisano@to.infn.it                     //
22 //                                                                           //
23 // NOTATIONS ADOPTED TO IDENTIFY DETECTORS (used in different ages!):        //
24 //   (ZN1,ZP1) or (ZNC, ZPC) or RIGHT refers to side C (RB26)                //
25 //   (ZN2,ZP2) or (ZNA, ZPA) or LEFT refers to side A (RB24)                 //
26 //                                                                           //
27 ///////////////////////////////////////////////////////////////////////////////
28
29
30 #include <TH2F.h>
31 #include <TH1D.h>
32 #include <TAxis.h>
33 #include <TMap.h>
34
35 #include "AliRawReader.h"
36 #include "AliGRPObject.h"
37 #include "AliESDEvent.h"
38 #include "AliESDZDC.h"
39 #include "AliZDCDigit.h"
40 #include "AliZDCRawStream.h"
41 #include "AliZDCReco.h"
42 #include "AliZDCReconstructor.h"
43 #include "AliZDCPedestals.h"
44 #include "AliZDCEnCalib.h"
45 #include "AliZDCTowerCalib.h"
46 #include "AliZDCRecoParam.h"
47 #include "AliZDCRecoParampp.h"
48 #include "AliZDCRecoParamPbPb.h"
49 #include "AliRunInfo.h"
50
51
52 ClassImp(AliZDCReconstructor)
53 AliZDCRecoParam *AliZDCReconstructor::fRecoParam=0;  //reconstruction parameters
54
55 //_____________________________________________________________________________
56 AliZDCReconstructor:: AliZDCReconstructor() :
57   fPedData(GetPedestalData()),
58   fEnCalibData(GetEnergyCalibData()),
59   fTowCalibData(GetTowerCalibData()),
60   fRecoMode(0),
61   fBeamEnergy(0.),
62   fNRun(0),
63   fIsCalibrationMB(kFALSE),
64   fPedSubMode(0),
65   fSignalThreshold(7)
66 {
67   // **** Default constructor
68
69 }
70
71
72 //_____________________________________________________________________________
73 AliZDCReconstructor::~AliZDCReconstructor()
74 {
75 // destructor
76 //   if(fRecoParam)    delete fRecoParam;
77    if(fPedData)      delete fPedData;    
78    if(fEnCalibData)  delete fEnCalibData;
79    if(fTowCalibData) delete fTowCalibData;
80 }
81
82 //____________________________________________________________________________
83 void AliZDCReconstructor::Init()
84 {
85   // Setting reconstruction mode
86   // Getting beam type and beam energy from GRP calibration object
87   
88   TString runType = GetRunInfo()->GetRunType();
89
90   if((runType.CompareTo("CALIBRATION_MB")) == 0){
91     fIsCalibrationMB = kTRUE;
92   }
93     
94   TString beamType = GetRunInfo()->GetBeamType();
95   // This is a temporary solution to allow reconstruction in tests without beam
96   if(((beamType.CompareTo("UNKNOWN"))==0) && ((runType.CompareTo("PHYSICS")) == 0)){
97     fRecoMode=1;
98   }
99   else if(beamType==AliGRPObject::GetInvalidString()){
100     AliWarning("GRP/GRP/Data entry:  missing value for the beam type !");
101     AliError("\t ZDC does not reconstruct event 4 UNKNOWN beam type\n");
102     return;
103   }
104
105   if(((beamType.CompareTo("pp"))==0) || ((beamType.CompareTo("p-p"))==0)
106      ||((beamType.CompareTo("PP"))==0) || ((beamType.CompareTo("P-P"))==0)) fRecoMode=1;
107   else if((beamType.CompareTo("A-A")) == 0){
108     fRecoMode=2;
109     if(fIsCalibrationMB == kTRUE){ 
110       fRecoParam = new AliZDCRecoParamPbPb();
111       //
112       TH2F* hZDCvsZEM = new TH2F("hZDCvsZEM","hZDCvsZEM",100,0.,10.,100,0.,1000.);
113       hZDCvsZEM->SetXTitle("E_{ZEM} (TeV)"); hZDCvsZEM->SetYTitle("E_{ZDC} (TeV)");
114       fRecoParam->SetZDCvsZEM(hZDCvsZEM);
115       //
116       TH2F* hZDCCvsZEM = new TH2F("hZDCCvsZEM","hZDCCvsZEM",100,0.,10.,100,0.,500.);
117       hZDCCvsZEM->SetXTitle("E_{ZEM} (TeV)"); hZDCCvsZEM->SetYTitle("E_{ZDCC} (TeV)");
118       fRecoParam->SetZDCCvsZEM(hZDCCvsZEM);
119       //
120       TH2F* hZDCAvsZEM = new TH2F("hZDCAvsZEM","hZDCAvsZEM",100,0.,10.,100,0.,500.);
121       hZDCAvsZEM->SetXTitle("E_{ZEM} (TeV)"); hZDCAvsZEM->SetYTitle("E_{ZDCA} (TeV)"); 
122       fRecoParam->SetZDCAvsZEM(hZDCAvsZEM);
123       //
124       AliInfo("\n ***** CALIBRATION_MB data -> building AliZDCRecoParamPbPb object *****");
125     }
126   }
127     
128   fBeamEnergy = GetRunInfo()->GetBeamEnergy();
129   if(fBeamEnergy==AliGRPObject::GetInvalidFloat()){
130     AliWarning("GRP/GRP/Data entry:  missing value for the beam energy ! Using 0.");
131     fBeamEnergy = 0.;
132   }
133     
134   if(fIsCalibrationMB==kFALSE)  
135     printf("\n\n ***** ZDC reconstruction initialized for %s @ %1.0f GeV *****\n\n",beamType.Data(), fBeamEnergy);
136   
137 }
138
139 //_____________________________________________________________________________
140 void AliZDCReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
141 {
142   // *** Local ZDC reconstruction for digits
143   // Works on the current event
144     
145   // Retrieving calibration data  
146   // Parameters for mean value pedestal subtraction
147   int const kNch = 24;
148   Float_t meanPed[2*kNch];    
149   for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedData->GetMeanPed(jj);
150   // Parameters pedestal subtraction through correlation with out-of-time signals
151   Float_t corrCoeff0[2*kNch], corrCoeff1[2*kNch];
152   for(Int_t jj=0; jj<2*kNch; jj++){
153      corrCoeff0[jj] =  fPedData->GetPedCorrCoeff0(jj);
154      corrCoeff1[jj] =  fPedData->GetPedCorrCoeff1(jj);
155   }
156
157   // get digits
158   AliZDCDigit digit;
159   AliZDCDigit* pdigit = &digit;
160   digitsTree->SetBranchAddress("ZDC", &pdigit);
161   //printf("\n\t # of digits in tree: %d\n",(Int_t) digitsTree->GetEntries());
162
163   // loop over digits
164   Float_t tZN1Corr[10], tZP1Corr[10], tZN2Corr[10], tZP2Corr[10]; 
165   Float_t dZEM1Corr[2], dZEM2Corr[2], sPMRef1[2], sPMRef2[2]; 
166   for(Int_t i=0; i<10; i++){
167      tZN1Corr[i] = tZP1Corr[i] = tZN2Corr[i] = tZP2Corr[i] = 0.;
168      if(i<2) dZEM1Corr[i] = dZEM2Corr[i] = sPMRef1[i] = sPMRef2[i] = 0.;
169   }  
170   
171   Int_t digNentries = digitsTree->GetEntries();
172   Float_t ootDigi[kNch];
173   // -- Reading out-of-time signals (last kNch entries) for current event
174   if(fPedSubMode==1){
175     for(Int_t iDigit=kNch; iDigit<digNentries; iDigit++){
176        ootDigi[iDigit] = digitsTree->GetEntry(iDigit);
177     }
178   }
179   
180   for(Int_t iDigit=0; iDigit<(digNentries/2); iDigit++) {
181    digitsTree->GetEntry(iDigit);
182    if (!pdigit) continue;
183    //  
184    Int_t det = digit.GetSector(0);
185    Int_t quad = digit.GetSector(1);
186    Int_t pedindex = -1;
187    Float_t ped2SubHg=0., ped2SubLg=0.;
188    if(quad!=5){
189      if(det==1)      pedindex = quad;
190      else if(det==2) pedindex = quad+5;
191      else if(det==3) pedindex = quad+9;
192      else if(det==4) pedindex = quad+12;
193      else if(det==5) pedindex = quad+17;
194    }
195    else pedindex = (det-1)/3+22;
196    //
197    if(fPedSubMode==0){
198      ped2SubHg = meanPed[pedindex];
199      ped2SubLg = meanPed[pedindex+kNch];
200    }
201    else if(fPedSubMode==1){
202      ped2SubHg = corrCoeff1[pedindex]*ootDigi[pedindex]+corrCoeff0[pedindex];
203      ped2SubLg = corrCoeff1[pedindex+kNch]*ootDigi[pedindex+kNch]+corrCoeff0[pedindex+kNch];
204    }
205       
206    if(quad != 5){ // ZDC (not reference PTMs!)
207     if(det == 1){ // *** ZNC
208        tZN1Corr[quad] = (Float_t) (digit.GetADCValue(0)-ped2SubHg);
209        tZN1Corr[quad+5] = (Float_t) (digit.GetADCValue(1)-ped2SubLg);
210        if(tZN1Corr[quad]<0.) tZN1Corr[quad] = 0.;
211        if(tZN1Corr[quad+5]<0.) tZN1Corr[quad+5] = 0.;
212        // Ch. debug
213        //printf("\t pedindex %d tZN1Corr[%d] = %1.0f tZN1Corr[%d] = %1.0f", 
214        //       pedindex, quad, tZN1Corr[quad], quad+5, tZN1Corr[quad+5]);
215     }
216     else if(det == 2){ // *** ZP1
217        tZP1Corr[quad] = (Float_t) (digit.GetADCValue(0)-ped2SubHg);
218        tZP1Corr[quad+5] = (Float_t) (digit.GetADCValue(1)-ped2SubLg);
219        if(tZP1Corr[quad]<0.) tZP1Corr[quad] = 0.;
220        if(tZP1Corr[quad+5]<0.) tZP1Corr[quad+5] = 0.;
221        // Ch. debug
222        //printf("\t pedindex %d tZP1Corr[%d] = %1.0f tZP1Corr[%d] = %1.0f", 
223        //       pedindex, quad, tZP1Corr[quad], quad+5, tZP1Corr[quad+5]);
224     }
225     else if(det == 3){
226        if(quad == 1){       // *** ZEM1  
227          dZEM1Corr[0] += (Float_t) (digit.GetADCValue(0)-ped2SubHg); 
228          dZEM1Corr[1] += (Float_t) (digit.GetADCValue(1)-ped2SubLg); 
229          if(dZEM1Corr[0]<0.) dZEM1Corr[0] = 0.;
230          if(dZEM1Corr[1]<0.) dZEM1Corr[1] = 0.;
231          // Ch. debug
232          //printf("\t pedindex %d tZEM1Corr[%d] = %1.0f tZEM1Corr[%d] = %1.0f", 
233          //     pedindex, quad, tZEM1Corr[quad], quad+1, tZEM1Corr[quad+1]);
234        }
235        else if(quad == 2){  // *** ZEM2
236          dZEM2Corr[0] += (Float_t) (digit.GetADCValue(0)-ped2SubHg); 
237          dZEM2Corr[1] += (Float_t) (digit.GetADCValue(1)-ped2SubLg); 
238          if(dZEM2Corr[0]<0.) dZEM2Corr[0] = 0.;
239          if(dZEM2Corr[1]<0.) dZEM2Corr[1] = 0.;
240          // Ch. debug
241          //printf("\t pedindex %d tZEM2Corr[%d] = %1.0f tZEM2Corr[%d] = %1.0f", 
242          //     pedindex, quad, tZEM2Corr[quad], quad+1, tZEM2Corr[quad+1]);
243        }
244     }
245     else if(det == 4){  // *** ZN2
246        tZN2Corr[quad] = (Float_t) (digit.GetADCValue(0)-ped2SubHg);
247        tZN2Corr[quad+5] = (Float_t) (digit.GetADCValue(1)-ped2SubLg);
248        if(tZN2Corr[quad]<0.) tZN2Corr[quad] = 0.;
249        if(tZN2Corr[quad+5]<0.) tZN2Corr[quad+5] = 0.;
250        // Ch. debug
251        //printf("\t pedindex %d tZN2Corr[%d] = %1.0f tZN2Corr[%d] = %1.0f\n", 
252        //       pedindex, quad, tZN2Corr[quad], quad+5, tZN2Corr[quad+5]);
253     }
254     else if(det == 5){  // *** ZP2 
255        tZP2Corr[quad] = (Float_t) (digit.GetADCValue(0)-ped2SubHg);
256        tZP2Corr[quad+5] = (Float_t) (digit.GetADCValue(1)-ped2SubLg);
257        if(tZP2Corr[quad]<0.) tZP2Corr[quad] = 0.;
258        if(tZP2Corr[quad+5]<0.) tZP2Corr[quad+5] = 0.;
259        // Ch. debug
260        //printf("\t pedindex %d tZP2Corr[%d] = %1.0f tZP2Corr[%d] = %1.0f\n", 
261        //       pedindex, quad, tZP2Corr[quad], quad+5, tZP2Corr[quad+5]);
262     }
263    }
264    else{ // Reference PMs
265      if(det == 1){
266        sPMRef1[0] = (Float_t) (digit.GetADCValue(0)-ped2SubHg);
267        sPMRef1[1] = (Float_t) (digit.GetADCValue(1)-ped2SubLg);
268        // Ch. debug
269        if(sPMRef1[0]<0.) sPMRef1[0] = 0.;
270        if(sPMRef2[1]<0.) sPMRef1[1] = 0.;
271      }
272      else if(det == 4){
273        sPMRef2[0] = (Float_t) (digit.GetADCValue(0)-ped2SubHg);
274        sPMRef2[1] = (Float_t) (digit.GetADCValue(1)-ped2SubLg);
275        // Ch. debug
276        if(sPMRef2[0]<0.) sPMRef2[0] = 0.;
277        if(sPMRef2[1]<0.) sPMRef2[1] = 0.;
278      }
279    }
280
281    // Ch. debug
282    /*printf(" - AliZDCReconstructor -> digit #%d det %d quad %d pedHG %1.0f pedLG %1.0f\n",
283          iDigit, det, quad, ped2SubHg, ped2SubLg);
284    printf("   HGChain -> RawDig %d DigCorr %1.2f\n", digit.GetADCValue(0), digit.GetADCValue(0)-ped2SubHg); 
285    printf("   LGChain -> RawDig %d DigCorr %1.2f\n", digit.GetADCValue(1), digit.GetADCValue(1)-ped2SubLg); 
286    */
287   }//digits loop
288  
289   
290   // If CALIBRATION_MB run build the RecoParam object 
291   if(fIsCalibrationMB){
292     Float_t ZDCC=0., ZDCA=0., ZEM=0;
293     ZEM += dZEM1Corr[0] + dZEM2Corr[0];
294     for(Int_t jkl=0; jkl<5; jkl++){
295        ZDCC += tZN1Corr[jkl] + tZP1Corr[jkl];
296        ZDCA += tZN2Corr[jkl] + tZP2Corr[jkl];
297     }
298     // Using energies in TeV in fRecoParam object
299     BuildRecoParam(ZDCC/1000., ZDCA/1000., ZEM/1000.);
300   }
301   Bool_t recFlag1 = kFALSE, recFlag2 = kFALSE, recFlag3 = kFALSE;
302   // reconstruct the event
303   if(fRecoMode==1)
304     ReconstructEventpp(clustersTree, tZN1Corr, tZP1Corr, tZN2Corr, tZP2Corr, 
305       dZEM1Corr, dZEM2Corr, sPMRef1, sPMRef2, recFlag1, recFlag2, recFlag3);
306   else if(fRecoMode==2)
307     ReconstructEventPbPb(clustersTree, tZN1Corr, tZP1Corr, tZN2Corr, tZP2Corr, 
308       dZEM1Corr, dZEM2Corr, sPMRef1, sPMRef2, recFlag1, recFlag2, recFlag3);
309 }
310
311 //_____________________________________________________________________________
312 void AliZDCReconstructor::Reconstruct(AliRawReader* rawReader, TTree* clustersTree) const
313 {
314   // *** ZDC raw data reconstruction
315   // Works on the current event
316   
317   // Retrieving calibration data  
318   // Parameters for pedestal subtraction
319   int const kNch = 24;
320   Float_t meanPed[2*kNch];    
321   for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedData->GetMeanPed(jj);
322   // Parameters pedestal subtraction through correlation with out-of-time signals
323   Float_t corrCoeff0[2*kNch], corrCoeff1[2*kNch];
324   for(Int_t jj=0; jj<2*kNch; jj++){
325      corrCoeff0[jj] =  fPedData->GetPedCorrCoeff0(jj);
326      corrCoeff1[jj] =  fPedData->GetPedCorrCoeff1(jj);
327   }
328
329   Int_t adcZN1[5], adcZN1oot[5], adcZN1lg[5], adcZN1ootlg[5];
330   Int_t adcZP1[5], adcZP1oot[5], adcZP1lg[5], adcZP1ootlg[5];
331   Int_t adcZN2[5], adcZN2oot[5], adcZN2lg[5], adcZN2ootlg[5];
332   Int_t adcZP2[5], adcZP2oot[5], adcZP2lg[5], adcZP2ootlg[5];
333   Int_t adcZEM[2], adcZEMoot[2], adcZEMlg[2], adcZEMootlg[2];
334   Int_t pmRef[2], pmRefoot[2], pmReflg[2], pmRefootlg[2];
335   for(Int_t ich=0; ich<5; ich++){
336     adcZN1[ich] = adcZN1oot[ich] = adcZN1lg[ich] = adcZN1ootlg[ich] = 0;
337     adcZP1[ich] = adcZP1oot[ich] = adcZP1lg[ich] = adcZP1ootlg[ich] = 0;
338     adcZN2[ich] = adcZN2oot[ich] = adcZN2lg[ich] = adcZN2ootlg[ich] = 0;
339     adcZP2[ich] = adcZP2oot[ich] = adcZP2lg[ich] = adcZP2ootlg[ich] = 0;
340     if(ich<2){
341       adcZEM[ich] = adcZEMoot[ich] = adcZEMlg[ich] = adcZEMootlg[ich] = 0;
342       pmRef[ich] = pmRefoot[ich] = pmReflg[ich] = pmRefootlg[ich] = 0;
343     }
344   }
345   
346   // loop over raw data
347   Float_t tZN1Corr[10], tZP1Corr[10], tZN2Corr[10], tZP2Corr[10]; 
348   Float_t dZEM1Corr[2], dZEM2Corr[2], sPMRef1[2], sPMRef2[2]; 
349   for(Int_t i=0; i<10; i++){
350      tZN1Corr[i] = tZP1Corr[i] = tZN2Corr[i] = tZP2Corr[i] = 0.;
351      if(i<2) dZEM1Corr[i] = dZEM2Corr[i] = sPMRef1[i] = sPMRef2[i] = 0.;
352   }  
353
354   //fNRun = (Int_t) rawReader->GetRunNumber();
355   Bool_t chOff=kFALSE, isUndflw=kFALSE, isOvflw=kFALSE, isADCEvGood=kFALSE;
356   Int_t kFirstADCGeo=0, kLastADCGeo=3, kScalerGeo=8, kPUGeo=29, kTrigScales=30, kTrigHistory=31;
357   //
358   rawReader->Reset();
359   AliZDCRawStream rawData(rawReader);
360   while(rawData.Next()){
361    //if(rawData.IsCalibration() == kFALSE){ // Reading ADCs
362    if((rawData.GetADCModule()>=kFirstADCGeo) && (rawData.GetADCModule()<=kLastADCGeo)){    
363     //printf(" **** Reading ADC raw data from module %d **** \n",rawData.GetADCModule());
364     //
365     if((rawData.IsADCDataWord()) && (rawData.GetNChannelsOn()<48))  chOff = kTRUE;
366     if((rawData.IsADCDataWord()) && (rawData.IsOverflow() == kTRUE)) isUndflw = kTRUE;
367     if((rawData.IsADCDataWord()) && (rawData.IsUnderflow() == kTRUE)) isOvflw = kTRUE;
368     if((rawData.IsADCDataWord()) && (rawData.IsADCEventGood() == kTRUE)) isADCEvGood = kTRUE;
369     
370     if((rawData.IsADCDataWord()) && (isUndflw==kFALSE) 
371         && (isOvflw==kFALSE) && (isADCEvGood==kTRUE)){
372      
373       Int_t adcMod = rawData.GetADCModule();
374       Int_t det = rawData.GetSector(0);
375       Int_t quad = rawData.GetSector(1);
376       Int_t gain = rawData.GetADCGain();
377       Int_t pedindex=0;
378       //
379       // Mean pedestal value subtraction -------------------------------------------------------
380       if(fPedSubMode == 0){
381        // Not interested in o.o.t. signals (ADC modules 2, 3)
382        if(adcMod == 2 || adcMod == 3) continue;
383        //
384        if(quad != 5){ // ZDCs (not reference PTMs)
385         if(det == 1){    
386           pedindex = quad;
387           if(gain == 0) tZN1Corr[quad]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]); 
388           else tZN1Corr[quad+5]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex+kNch]); 
389         }
390         else if(det == 2){ 
391           pedindex = quad+5;
392           if(gain == 0) tZP1Corr[quad]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]); 
393           else tZP1Corr[quad+5]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex+kNch]); 
394         }
395         else if(det == 3){ 
396           pedindex = quad+9;
397           if(quad==1){     
398             if(gain == 0) dZEM1Corr[0] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]); 
399             else dZEM1Corr[1] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex+kNch]); 
400           }
401           else if(quad==2){ 
402             if(gain == 0) dZEM2Corr[0] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]); 
403             else dZEM2Corr[1] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex+kNch]); 
404           }
405         }
406         else if(det == 4){       
407           pedindex = quad+12;
408           if(gain == 0) tZN2Corr[quad]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]); 
409           else tZN2Corr[quad+5]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex+kNch]); 
410         }
411         else if(det == 5){
412           pedindex = quad+17;
413           if(gain == 0) tZP2Corr[quad]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]); 
414           else tZP2Corr[quad+5]  += (Float_t) (rawData.GetADCValue()-meanPed[pedindex+kNch]); 
415         }
416        }
417        else{ // reference PM
418          pedindex = (det-1)/3 + 22;
419          if(det == 1){
420            if(gain==0) sPMRef1[0] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]);
421          else sPMRef1[1] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]);
422          }
423          else if(det == 4){
424            if(gain==0) sPMRef2[0] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]);
425            else sPMRef2[1] += (Float_t) (rawData.GetADCValue()-meanPed[pedindex]);
426          }
427        }
428        // Ch. debug
429        /*printf(" -> AliZDCReconstructor: det %d quad %d res %d -> Pedestal[%d] %1.0f\n", 
430          det,quad,gain, pedindex, meanPed[pedindex]);
431        printf(" -> AliZDCReconstructor: RawADC %1.0f ADCCorr %1.0f\n", 
432          rawData.GetADCValue(), rawData.GetADCValue()-meanPed[pedindex]);*/
433          
434       }// mean pedestal subtraction
435       // Pedestal subtraction from correlation ------------------------------------------------
436       else if(fPedSubMode == 1){
437        // In time signals
438        if(adcMod==0 || adcMod==1){
439          if(quad != 5){ // signals from ZDCs
440            if(det == 1){
441              if(gain==0) adcZN1[quad] = rawData.GetADCValue();
442              else adcZN1lg[quad] = rawData.GetADCValue();
443            }
444            else if(det == 2){
445              if(gain==0) adcZP1[quad] = rawData.GetADCValue();
446              else adcZP1lg[quad] = rawData.GetADCValue();
447            }
448            else if(det == 3){
449              if(gain==0) adcZEM[quad-1] = rawData.GetADCValue();
450              else adcZEMlg[quad-1] = rawData.GetADCValue();
451            }
452            else if(det == 4){
453              if(gain==0) adcZN2[quad] = rawData.GetADCValue();
454              else adcZN2lg[quad] = rawData.GetADCValue();
455            }
456            else if(det == 5){
457              if(gain==0) adcZP2[quad] = rawData.GetADCValue();
458              else adcZP2lg[quad] = rawData.GetADCValue();
459            }
460          }
461          else{ // signals from reference PM
462             if(gain==0) pmRef[quad-1] = rawData.GetADCValue();
463             else pmReflg[quad-1] = rawData.GetADCValue();
464          }
465        }
466        // Out-of-time pedestals
467        else if(adcMod==2 || adcMod==3){
468          if(quad != 5){ // signals from ZDCs
469            if(det == 1){
470              if(gain==0) adcZN1oot[quad] = rawData.GetADCValue();
471              else adcZN1ootlg[quad] = rawData.GetADCValue();
472            }
473            else if(det == 2){
474              if(gain==0) adcZP1oot[quad] = rawData.GetADCValue();
475              else adcZP1ootlg[quad] = rawData.GetADCValue();
476            }
477            else if(det == 3){
478              if(gain==0) adcZEMoot[quad-1] = rawData.GetADCValue();
479              else adcZEMootlg[quad-1] = rawData.GetADCValue();
480            }
481            else if(det == 4){
482              if(gain==0) adcZN2oot[quad] = rawData.GetADCValue();
483              else adcZN2ootlg[quad] = rawData.GetADCValue();
484            }
485            else if(det == 5){
486              if(gain==0) adcZP2oot[quad] = rawData.GetADCValue();
487              else adcZP2ootlg[quad] = rawData.GetADCValue();
488            }
489          }
490          else{ // signals from reference PM
491             if(gain==0) pmRefoot[quad-1] = rawData.GetADCValue();
492             else pmRefootlg[quad-1] = rawData.GetADCValue();
493          }
494        }
495       } // pedestal subtraction from correlation
496       // Ch. debug
497       //printf("\t AliZDCReconstructor - det %d quad %d res %d -> Ped[%d] = %1.0f\n", 
498       //  det,quad,gain, pedindex, meanPed[pedindex]);
499     }//IsADCDataWord
500    }// Not raw data from calibration run!
501    else{
502    }
503   }//loop on raw data
504   
505   if(fPedSubMode==1){
506     for(Int_t t=0; t<5; t++){
507        tZN1Corr[t] = adcZN1[t] - (corrCoeff1[t]*adcZN1oot[t]+corrCoeff0[t]);
508        tZN1Corr[t+5] = adcZN1lg[t] - (corrCoeff1[t+kNch]*adcZN1ootlg[t]+corrCoeff0[t+kNch]);
509        //
510        tZP1Corr[t] = adcZP1[t] - (corrCoeff1[t+5]*adcZP1oot[t]+corrCoeff0[t+5]);
511        tZP1Corr[t+5] = adcZP1lg[t] - (corrCoeff1[t+5+kNch]*adcZP1ootlg[t]+corrCoeff0[t+5+kNch]);
512        //
513        tZN2Corr[t] = adcZN2[t] - (corrCoeff1[t+12]*adcZN2oot[t]+corrCoeff0[t+12]);
514        tZN2Corr[t+5] = adcZN2lg[t] - (corrCoeff1[t+12+kNch]*adcZN2ootlg[t]+corrCoeff0[t+12+kNch]);
515        //
516        tZP2Corr[t] = adcZP2[t] - (corrCoeff1[t+17]*adcZP2oot[t]+corrCoeff0[t+17]);
517        tZP2Corr[t+5] = adcZP2lg[t] - (corrCoeff1[t+17+kNch]*adcZP2ootlg[t]+corrCoeff0[t+17+kNch]);
518        // 0---------0 Ch. debug 0---------0
519 /*       printf("\n\n ---------- Debug of pedestal subtraction from correlation ----------\n");
520        printf("\tCorrCoeff0\tCorrCoeff1\n");
521        printf(" ZN1 %d\t%1.0f\t%1.0f\n",t,corrCoeff0[t],corrCoeff1[t]);
522        printf(" ZN1lg %d\t%1.0f\t%1.0f\n",t+kNch,corrCoeff0[t+kNch],corrCoeff1[t+kNch]);
523        printf(" ZP1 %d\t%1.0f\t%1.0f\n",t+5,corrCoeff0[t+5],corrCoeff1[t+5]);
524        printf(" ZP1lg %d\t%1.0f\t%1.0f\n",t+5+kNch,corrCoeff0[t+5+kNch],corrCoeff1[t+5+kNch]);
525        printf(" ZN2 %d\t%1.0f\t%1.0f\n",t+12,corrCoeff0[t+12],corrCoeff1[t+12]);
526        printf(" ZN2lg %d\t%1.0f\t%1.0f\n",t+12+kNch,corrCoeff0[t+12+kNch],corrCoeff1[t+12+kNch]);
527        printf(" ZP2 %d\t%1.0f\t%1.0f\n",t+17,corrCoeff0[t+17],corrCoeff1[t+17]);
528        printf(" ZP2lg %d\t%1.0f\t%1.0f\n",t+17+kNch,corrCoeff0[t+17+kNch],corrCoeff1[t+17+kNch]);
529        
530        printf("ZN1 ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
531                 adcZN1[t],(corrCoeff1[t]*adcZN1oot[t]+corrCoeff0[t]),tZN1Corr[t]);
532        printf(" lg ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
533                 adcZN1lg[t],(corrCoeff1[t+kNch]*adcZN1ootlg[t]+corrCoeff0[t+kNch]),tZN1Corr[t+5]);
534        //
535        printf("ZP1 ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
536                 adcZP1[t],(corrCoeff1[t+5]*adcZP1oot[t]+corrCoeff0[t+5]),tZP1Corr[t]);
537        printf(" lg ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
538                 adcZP1lg[t],(corrCoeff1[t+5+kNch]*adcZP1ootlg[t]+corrCoeff0[t+5+kNch]),tZP1Corr[t+5]);
539        //
540        printf("ZN2 ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
541                 adcZN2[t],(corrCoeff1[t+12]*adcZN2oot[t]+corrCoeff0[t+12]),tZN2Corr[t]);
542        printf(" lg ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
543                 adcZN2lg[t],(corrCoeff1[t+12+kNch]*adcZN2ootlg[t]+corrCoeff0[t+12+kNch]),tZN2Corr[t+5]);
544        //
545        printf("ZP2 ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
546                 adcZP2[t],(corrCoeff1[t+17]*adcZP2oot[t]+corrCoeff0[t+17]),tZP2Corr[t]);
547        printf(" lg ->  rawADC %d\tpedestal%1.2f\tcorrADC%1.2f\n",
548                 adcZP2lg[t],(corrCoeff1[t+17+kNch]*adcZP2ootlg[t]+corrCoeff0[t+17+kNch]),tZP2Corr[t+5]);
549 */
550     }
551     dZEM1Corr[0] = adcZEM[0]   - (corrCoeff1[9]*adcZEMoot[0]+corrCoeff0[9]);
552     dZEM1Corr[1] = adcZEMlg[0] - (corrCoeff1[9+kNch]*adcZEMootlg[0]+corrCoeff0[9+kNch]);
553     dZEM2Corr[0] = adcZEM[1]   - (corrCoeff1[10]*adcZEMoot[1]+corrCoeff0[10]);
554     dZEM2Corr[1] = adcZEMlg[1] - (corrCoeff1[10+kNch]*adcZEMootlg[1]+corrCoeff0[10+kNch]);
555     //
556     sPMRef1[0] = pmRef[0]   - (corrCoeff1[22]*pmRefoot[0]+corrCoeff0[22]);
557     sPMRef1[1] = pmReflg[0] - (corrCoeff1[22+kNch]*pmRefootlg[0]+corrCoeff0[22+kNch]);
558     sPMRef2[0] = pmRef[0]   - (corrCoeff1[23]*pmRefoot[1]+corrCoeff0[23]);
559     sPMRef2[1] = pmReflg[0] - (corrCoeff1[23+kNch]*pmRefootlg[1]+corrCoeff0[23+kNch]);
560   }
561     
562   // If CALIBRATION_MB run build the RecoParam object 
563   if(fIsCalibrationMB){
564     Float_t ZDCC=0., ZDCA=0., ZEM=0;
565     ZEM += dZEM1Corr[0] + dZEM2Corr[0];
566     for(Int_t jkl=0; jkl<5; jkl++){
567        ZDCC += tZN1Corr[jkl] + tZP1Corr[jkl];
568        ZDCA += tZN2Corr[jkl] + tZP2Corr[jkl];
569     }
570     BuildRecoParam(ZDCC/100., ZDCA/100., ZEM/100.);
571   }
572   // reconstruct the event
573   else{
574     if(fRecoMode==1) // p-p data
575       ReconstructEventpp(clustersTree, tZN1Corr, tZP1Corr, tZN2Corr, tZP2Corr, 
576         dZEM1Corr, dZEM2Corr, sPMRef1, sPMRef2, chOff, isUndflw, isOvflw);
577     else if(fRecoMode==2) // Pb-Pb data
578       ReconstructEventPbPb(clustersTree, tZN1Corr, tZP1Corr, tZN2Corr, tZP2Corr, 
579         dZEM1Corr, dZEM2Corr, sPMRef1, sPMRef2, chOff, isUndflw, isOvflw);
580   }
581 }
582
583 //_____________________________________________________________________________
584 void AliZDCReconstructor::ReconstructEventpp(TTree *clustersTree, Float_t* corrADCZN1, 
585         Float_t* corrADCZP1, Float_t* corrADCZN2, Float_t* corrADCZP2,
586         Float_t* corrADCZEM1, Float_t* corrADCZEM2, Float_t* sPMRef1, Float_t* sPMRef2,
587         Bool_t channelsOff, Bool_t chUnderflow, Bool_t chOverflow) const
588 {
589   // ****************** Reconstruct one event ******************
590     
591   // ---- Setting reco flags
592   UInt_t recoFlag=0;
593   UInt_t rFlags[32];
594   for(Int_t ifl=0; ifl<32; ifl++) rFlags[ifl]=0;
595
596   if(corrADCZN2[0]>fSignalThreshold)  rFlags[0] = 0x1;
597   if(corrADCZP2[0]>fSignalThreshold)  rFlags[1] = 0x1;
598   if(corrADCZEM1[0]>fSignalThreshold) rFlags[2] = 0x1;
599   if(corrADCZEM2[0]>fSignalThreshold) rFlags[3] = 0x1;
600   if(corrADCZN1[0]>fSignalThreshold)  rFlags[4] = 0x1;
601   if(corrADCZP1[0]>fSignalThreshold)  rFlags[5] = 0x1;
602   //
603   if(channelsOff==kTRUE)   rFlags[8] = 0x1;
604   if(chUnderflow == kTRUE) rFlags[9] = 0x1;
605   if(chOverflow==kTRUE)    rFlags[10] = 0x1;
606   recoFlag = rFlags[10] << 10 | rFlags[9] << 9 | rFlags[8] << 8 |
607              rFlags[5] << 5 | rFlags[4] << 4 | rFlags[3] << 3 |
608              rFlags[2] << 2 | rFlags[1] << 1 | rFlags[0];
609
610   // ******     Retrieving calibration data 
611   // --- Equalization coefficients ---------------------------------------------
612   Float_t equalCoeffZN1[5], equalCoeffZP1[5], equalCoeffZN2[5], equalCoeffZP2[5];
613   for(Int_t ji=0; ji<5; ji++){
614      equalCoeffZN1[ji] = fTowCalibData->GetZN1EqualCoeff(ji);
615      equalCoeffZP1[ji] = fTowCalibData->GetZP1EqualCoeff(ji); 
616      equalCoeffZN2[ji] = fTowCalibData->GetZN2EqualCoeff(ji); 
617      equalCoeffZP2[ji] = fTowCalibData->GetZP2EqualCoeff(ji); 
618   }
619   // --- Energy calibration factors ------------------------------------
620   Float_t calibEne[4];
621   // **** Energy calibration coefficient set to 1 
622   // **** (no trivial way to calibrate in p-p runs)
623   for(Int_t ij=0; ij<4; ij++) calibEne[ij] = fEnCalibData->GetEnCalib(ij);
624   
625   // ******     Equalization of detector responses
626   Float_t equalTowZN1[10], equalTowZN2[10], equalTowZP1[10], equalTowZP2[10];
627   for(Int_t gi=0; gi<10; gi++){
628      if(gi<5){
629        equalTowZN1[gi] = corrADCZN1[gi]*equalCoeffZN1[gi];
630        equalTowZP1[gi] = corrADCZP1[gi]*equalCoeffZP1[gi];
631        equalTowZN2[gi] = corrADCZN2[gi]*equalCoeffZN2[gi];
632        equalTowZP2[gi] = corrADCZP2[gi]*equalCoeffZP2[gi];
633      }
634      else{
635        equalTowZN1[gi] = corrADCZN1[gi]*equalCoeffZN1[gi-5];
636        equalTowZP1[gi] = corrADCZP1[gi]*equalCoeffZP1[gi-5];
637        equalTowZN2[gi] = corrADCZN2[gi]*equalCoeffZN2[gi-5];
638        equalTowZP2[gi] = corrADCZP2[gi]*equalCoeffZP2[gi-5];
639      }
640   }
641   
642   // ******     Summed response for hadronic calorimeter (SUMMED and then CALIBRATED!)
643   Float_t calibSumZN1[]={0,0}, calibSumZN2[]={0,0}, calibSumZP1[]={0,0}, calibSumZP2[]={0,0};
644   for(Int_t gi=0; gi<5; gi++){
645        calibSumZN1[0] += equalTowZN1[gi];
646        calibSumZP1[0] += equalTowZP1[gi];
647        calibSumZN2[0] += equalTowZN2[gi];
648        calibSumZP2[0] += equalTowZP2[gi];
649        //
650        calibSumZN1[1] += equalTowZN1[gi+5];
651        calibSumZP1[1] += equalTowZP1[gi+5];
652        calibSumZN2[1] += equalTowZN2[gi+5];
653        calibSumZP2[1] += equalTowZP2[gi+5];
654   }
655   // High gain chain
656   calibSumZN1[0] = calibSumZN1[0]*calibEne[0]/8.;
657   calibSumZP1[0] = calibSumZP1[0]*calibEne[1]/8.;
658   calibSumZN2[0] = calibSumZN2[0]*calibEne[2]/8.;
659   calibSumZP2[0] = calibSumZP2[0]*calibEne[3]/8.;
660   // Low gain chain
661   calibSumZN1[1] = calibSumZN1[1]*calibEne[0];
662   calibSumZP1[1] = calibSumZP1[1]*calibEne[1];
663   calibSumZN2[1] = calibSumZN2[1]*calibEne[2];
664   calibSumZP2[1] = calibSumZP2[1]*calibEne[3];
665   
666   // ******     Energy calibration of detector responses
667   Float_t calibTowZN1[10], calibTowZN2[10], calibTowZP1[10], calibTowZP2[10];
668   for(Int_t gi=0; gi<5; gi++){
669      // High gain chain
670      calibTowZN1[gi] = equalTowZN1[gi]*calibEne[0]/8.;
671      calibTowZP1[gi] = equalTowZP1[gi]*calibEne[1]/8.;
672      calibTowZN2[gi] = equalTowZN2[gi]*calibEne[2]/8.;
673      calibTowZP2[gi] = equalTowZP2[gi]*calibEne[3]/8.;
674      // Low gain chain
675      calibTowZN1[gi+5] = equalTowZN1[gi+5]*calibEne[0];
676      calibTowZP1[gi+5] = equalTowZP1[gi+5]*calibEne[1];
677      calibTowZN2[gi+5] = equalTowZN2[gi+5]*calibEne[2];
678      calibTowZP2[gi+5] = equalTowZP2[gi+5]*calibEne[3];
679   }
680   //
681   Float_t sumZEM[]={0,0}, calibZEM1[]={0,0}, calibZEM2[]={0,0};
682   calibZEM1[0] = corrADCZEM1[0]*calibEne[5]/8.;
683   calibZEM1[1] = corrADCZEM1[1]*calibEne[5];
684   calibZEM2[0] = corrADCZEM2[0]*calibEne[5]/8.;
685   calibZEM2[1] = corrADCZEM2[1]*calibEne[5];
686   for(Int_t k=0; k<2; k++) sumZEM[k] = calibZEM1[k] + calibZEM2[k];
687   
688   //  ******    No. of spectator and participants nucleons
689   //  Variables calculated to comply with ESD structure
690   //  *** N.B. -> They have a meaning only in Pb-Pb!!!!!!!!!!!!
691   Int_t nDetSpecNLeft=0, nDetSpecPLeft=0, nDetSpecNRight=0, nDetSpecPRight=0;
692   Int_t nGenSpec=0, nGenSpecLeft=0, nGenSpecRight=0;
693   Int_t nPart=0, nPartTotLeft=0, nPartTotRight=0;
694   Double_t impPar=0., impPar1=0., impPar2=0.;
695   
696   // create the output tree
697   AliZDCReco* reco = new AliZDCReco(calibSumZN1, calibSumZP1, calibSumZN2, calibSumZP2, 
698                    calibTowZN1, calibTowZP1, calibTowZN2, calibTowZP2, 
699                    calibZEM1, calibZEM2, sPMRef1, sPMRef2,
700                    nDetSpecNLeft, nDetSpecPLeft, nDetSpecNRight, nDetSpecPRight, 
701                    nGenSpec, nGenSpecLeft, nGenSpecRight, 
702                    nPart, nPartTotLeft, nPartTotRight, 
703                    impPar, impPar1, impPar2,
704                    recoFlag);
705                   
706   const Int_t kBufferSize = 4000;
707   clustersTree->Branch("ZDC", "AliZDCReco", &reco, kBufferSize);
708   // write the output tree
709   clustersTree->Fill();
710 }
711
712 //_____________________________________________________________________________
713 void AliZDCReconstructor::ReconstructEventPbPb(TTree *clustersTree, 
714         Float_t* corrADCZN1, Float_t* corrADCZP1, Float_t* corrADCZN2, Float_t* corrADCZP2,
715         Float_t* corrADCZEM1, Float_t* corrADCZEM2, Float_t* sPMRef1, Float_t* sPMRef2,
716         Bool_t channelsOff, Bool_t chUnderflow, Bool_t chOverflow) const 
717 {
718   // ****************** Reconstruct one event ******************
719   
720   // ---- Setting reco flags
721   UInt_t recoFlag=0;
722   UInt_t rFlags[32];
723   for(Int_t ifl=0; ifl<32; ifl++) rFlags[ifl]=0;
724
725   if(corrADCZN2[0]>fSignalThreshold)  rFlags[0] = 0x1;
726   if(corrADCZP2[0]>fSignalThreshold)  rFlags[1] = 0x1;
727   if(corrADCZEM1[0]>fSignalThreshold) rFlags[2] = 0x1;
728   if(corrADCZEM2[0]>fSignalThreshold) rFlags[3] = 0x1;
729   if(corrADCZN1[0]>fSignalThreshold)  rFlags[4] = 0x1;
730   if(corrADCZP1[0]>fSignalThreshold)  rFlags[5] = 0x1;
731   //
732   if(channelsOff==kTRUE)   rFlags[8] = 0x1;
733   if(chUnderflow == kTRUE) rFlags[9] = 0x1;
734   if(chOverflow==kTRUE)    rFlags[10] = 0x1;
735   recoFlag = rFlags[10] << 10 | rFlags[9] << 9 | rFlags[8] << 8 |
736              rFlags[5] << 5 | rFlags[4] << 4 | rFlags[3] << 3 |
737              rFlags[2] << 2 | rFlags[1] << 1 | rFlags[0];
738
739   // ******     Retrieving calibration data 
740   // --- Equalization coefficients ---------------------------------------------
741   Float_t equalCoeffZN1[5], equalCoeffZP1[5], equalCoeffZN2[5], equalCoeffZP2[5];
742   for(Int_t ji=0; ji<5; ji++){
743      equalCoeffZN1[ji] = fTowCalibData->GetZN1EqualCoeff(ji);
744      equalCoeffZP1[ji] = fTowCalibData->GetZP1EqualCoeff(ji); 
745      equalCoeffZN2[ji] = fTowCalibData->GetZN2EqualCoeff(ji); 
746      equalCoeffZP2[ji] = fTowCalibData->GetZP2EqualCoeff(ji); 
747   }
748   // --- Energy calibration factors ------------------------------------
749   Float_t valFromOCDB[6], calibEne[6];
750   for(Int_t ij=0; ij<6; ij++){
751     valFromOCDB[ij] = fEnCalibData->GetEnCalib(ij);
752     if(ij<4){
753       if(valFromOCDB[ij]!=0) calibEne[ij] = fBeamEnergy/valFromOCDB[ij];
754       else AliWarning(" Value from OCDB for E calibration = 0 !!!\n");
755     } 
756     else calibEne[ij] = valFromOCDB[ij];
757   }
758   
759   // ******     Equalization of detector responses
760   Float_t equalTowZN1[10], equalTowZN2[10], equalTowZP1[10], equalTowZP2[10];
761   for(Int_t gi=0; gi<10; gi++){
762      if(gi<5){
763        equalTowZN1[gi] = corrADCZN1[gi]*equalCoeffZN1[gi];
764        equalTowZP1[gi] = corrADCZP1[gi]*equalCoeffZP1[gi];
765        equalTowZN2[gi] = corrADCZN2[gi]*equalCoeffZN2[gi];
766        equalTowZP2[gi] = corrADCZP2[gi]*equalCoeffZP2[gi];
767      }
768      else{
769        equalTowZN1[gi] = corrADCZN1[gi]*equalCoeffZN1[gi-5];
770        equalTowZP1[gi] = corrADCZP1[gi]*equalCoeffZP1[gi-5];
771        equalTowZN2[gi] = corrADCZN2[gi]*equalCoeffZN2[gi-5];
772        equalTowZP2[gi] = corrADCZP2[gi]*equalCoeffZP2[gi-5];
773      }
774   }
775   
776   // ******     Summed response for hadronic calorimeter (SUMMED and then CALIBRATED!)
777   Float_t calibSumZN1[]={0,0}, calibSumZN2[]={0,0}, calibSumZP1[]={0,0}, calibSumZP2[]={0,0};
778   for(Int_t gi=0; gi<5; gi++){
779        calibSumZN1[0] += equalTowZN1[gi];
780        calibSumZP1[0] += equalTowZP1[gi];
781        calibSumZN2[0] += equalTowZN2[gi];
782        calibSumZP2[0] += equalTowZP2[gi];
783        //
784        calibSumZN1[1] += equalTowZN1[gi+5];
785        calibSumZP1[1] += equalTowZP1[gi+5];
786        calibSumZN2[1] += equalTowZN2[gi+5];
787        calibSumZP2[1] += equalTowZP2[gi+5];
788   }
789   // High gain chain
790   calibSumZN1[0] = calibSumZN1[0]*calibEne[0]/8.;
791   calibSumZP1[0] = calibSumZP1[0]*calibEne[1]/8.;
792   calibSumZN2[0] = calibSumZN2[0]*calibEne[2]/8.;
793   calibSumZP2[0] = calibSumZP2[0]*calibEne[3]/8.;
794   // Low gain chain
795   calibSumZN1[1] = calibSumZN1[1]*calibEne[0];
796   calibSumZP1[1] = calibSumZP1[1]*calibEne[1];
797   calibSumZN2[1] = calibSumZN2[1]*calibEne[2];
798   calibSumZP2[1] = calibSumZP2[1]*calibEne[3];
799   //
800   Float_t sumZEM[]={0,0}, calibZEM1[]={0,0}, calibZEM2[]={0,0};
801   calibZEM1[0] = corrADCZEM1[0]*calibEne[5]/8.;
802   calibZEM1[1] = corrADCZEM1[1]*calibEne[5];
803   calibZEM2[0] = corrADCZEM2[0]*calibEne[5]/8.;
804   calibZEM2[1] = corrADCZEM2[1]*calibEne[5];
805   for(Int_t k=0; k<2; k++) sumZEM[k] = calibZEM1[k] + calibZEM2[k];
806     
807   // ******     Energy calibration of detector responses
808   Float_t calibTowZN1[10], calibTowZN2[10], calibTowZP1[10], calibTowZP2[10];
809   for(Int_t gi=0; gi<5; gi++){
810      // High gain chain
811      calibTowZN1[gi] = equalTowZN1[gi]*calibEne[0]/8.;
812      calibTowZP1[gi] = equalTowZP1[gi]*calibEne[1]/8.;
813      calibTowZN2[gi] = equalTowZN2[gi]*calibEne[2]/8.;
814      calibTowZP2[gi] = equalTowZP2[gi]*calibEne[3]/8.;
815      // Low gain chain
816      calibTowZN1[gi+5] = equalTowZN1[gi+5]*calibEne[0];
817      calibTowZP1[gi+5] = equalTowZP1[gi+5]*calibEne[1];
818      calibTowZN2[gi+5] = equalTowZN2[gi+5]*calibEne[2];
819      calibTowZP2[gi+5] = equalTowZP2[gi+5]*calibEne[3];
820   }
821   
822   //  ******    Number of detected spectator nucleons
823   Int_t nDetSpecNLeft=0, nDetSpecPLeft=0, nDetSpecNRight=0, nDetSpecPRight=0;
824   if(fBeamEnergy!=0){
825     nDetSpecNLeft = (Int_t) (calibSumZN1[0]/fBeamEnergy);
826     nDetSpecPLeft = (Int_t) (calibSumZP1[0]/fBeamEnergy);
827     nDetSpecNRight = (Int_t) (calibSumZN2[0]/fBeamEnergy);
828     nDetSpecPRight = (Int_t) (calibSumZP2[0]/fBeamEnergy);
829   }
830   else AliWarning(" ATTENTION!!! fBeamEnergy=0 -> N_spec will be ZERO!!! \n");
831   /*printf("\n\t AliZDCReconstructor -> nDetSpecNLeft %d, nDetSpecPLeft %d,"
832     " nDetSpecNRight %d, nDetSpecPRight %d\n",nDetSpecNLeft, nDetSpecPLeft, 
833     nDetSpecNRight, nDetSpecPRight);*/
834   
835   Int_t nGenSpec=0, nGenSpecA=0, nGenSpecC=0;
836   Int_t nPart=0, nPartA=0, nPartC=0;
837   Double_t b=0., bA=0., bC=0.;
838   
839   if(fIsCalibrationMB == kFALSE){
840     // ******   Reconstruction parameters ------------------ 
841     // Ch. debug
842     //fRecoParam->Print("");
843     //
844  
845     if (!fRecoParam) fRecoParam = const_cast<AliZDCRecoParam*>(GetRecoParam()); 
846  
847     TH2F *hZDCvsZEM  = fRecoParam->GethZDCvsZEM();
848     TH2F *hZDCCvsZEM = fRecoParam->GethZDCCvsZEM();
849     TH2F *hZDCAvsZEM = fRecoParam->GethZDCAvsZEM();
850     TH1D *hNpartDist = fRecoParam->GethNpartDist();
851     TH1D *hbDist = fRecoParam->GethbDist();    
852     Float_t ClkCenter = fRecoParam->GetClkCenter();
853     //
854     Double_t xHighEdge = hZDCvsZEM->GetXaxis()->GetXmax();
855     Double_t origin = xHighEdge*ClkCenter;
856     // Ch. debug
857     //printf("\n\n  xHighEdge %1.2f, origin %1.4f \n", xHighEdge, origin);
858     //
859     // ====> Summed ZDC info (sideA+side C)
860     TF1 *line = new TF1("line","[0]*x+[1]",0.,xHighEdge);
861     Float_t y = (calibSumZN1[0]+calibSumZP1[0]+calibSumZN2[0]+calibSumZP2[0])/1000.;
862     Float_t x = (calibZEM1[0]+calibZEM2[0])/1000.;
863     line->SetParameter(0, y/(x-origin));
864     line->SetParameter(1, -origin*y/(x-origin));
865     // Ch. debug
866     //printf("  ***************** Summed ZDC info (sideA+side C) \n");
867     //printf("  E_{ZEM} %1.4f, E_{ZDC} %1.2f, TF1: %1.2f*x + %1.2f   ", x, y,y/(x-origin),-origin*y/(x-origin));
868     //
869     Double_t countPerc=0;
870     Double_t xBinCenter=0, yBinCenter=0;
871     for(Int_t nbinx=1; nbinx<=hZDCvsZEM->GetNbinsX(); nbinx++){
872       for(Int_t nbiny=1; nbiny<=hZDCvsZEM->GetNbinsY(); nbiny++){
873          xBinCenter = hZDCvsZEM->GetXaxis()->GetBinCenter(nbinx);
874          yBinCenter = hZDCvsZEM->GetYaxis()->GetBinCenter(nbiny);
875          //
876          if(line->GetParameter(0)>0){
877            if(yBinCenter < (line->GetParameter(0)*xBinCenter + line->GetParameter(1))){
878              countPerc += hZDCvsZEM->GetBinContent(nbinx,nbiny);
879              // Ch. debug
880              /*printf(" xBinCenter  %1.3f, yBinCenter %1.0f,  countPerc %1.0f\n", 
881                 xBinCenter, yBinCenter, countPerc);*/
882            }
883          }
884          else{
885            if(yBinCenter > (line->GetParameter(0)*xBinCenter + line->GetParameter(1))){
886              countPerc += hZDCvsZEM->GetBinContent(nbinx,nbiny);
887              // Ch. debug
888              /*printf(" xBinCenter  %1.3f, yBinCenter %1.0f,  countPerc %1.0f\n", 
889                 xBinCenter, yBinCenter, countPerc);*/
890            }
891          }
892       }
893     }
894     //
895     Double_t xSecPerc = 0.;
896     if(hZDCvsZEM->GetEntries()!=0){ 
897       xSecPerc = countPerc/hZDCvsZEM->GetEntries();
898     }
899     else{
900       AliWarning("  Histogram hZDCvsZEM from OCDB has no entries!!!");
901     }
902     // Ch. debug
903     //printf("  xSecPerc %1.4f  \n", xSecPerc);
904
905     // ====> side C
906     TF1 *lineC = new TF1("lineC","[0]*x+[1]",0.,xHighEdge);
907     Float_t yC = (calibSumZN1[0]+calibSumZP1[0])/1000.;
908     lineC->SetParameter(0, yC/(x-origin));
909     lineC->SetParameter(1, -origin*yC/(x-origin));
910     // Ch. debug
911     //printf("  ***************** Side C \n");
912     //printf("  E_{ZEM} %1.4f, E_{ZDCC} %1.2f, TF1: %1.2f*x + %1.2f   ", x, yC,yC/(x-origin),-origin*yC/(x-origin));
913     //
914     Double_t countPercC=0;
915     Double_t xBinCenterC=0, yBinCenterC=0;
916     for(Int_t nbinx=1; nbinx<=hZDCCvsZEM->GetNbinsX(); nbinx++){
917       for(Int_t nbiny=1; nbiny<=hZDCCvsZEM->GetNbinsY(); nbiny++){
918          xBinCenterC = hZDCCvsZEM->GetXaxis()->GetBinCenter(nbinx);
919          yBinCenterC = hZDCCvsZEM->GetYaxis()->GetBinCenter(nbiny);
920          if(lineC->GetParameter(0)>0){
921            if(yBinCenterC < (lineC->GetParameter(0)*xBinCenterC + lineC->GetParameter(1))){
922              countPercC += hZDCCvsZEM->GetBinContent(nbinx,nbiny);
923            }
924          }
925          else{
926            if(yBinCenterC > (lineC->GetParameter(0)*xBinCenterC + lineC->GetParameter(1))){
927              countPercC += hZDCCvsZEM->GetBinContent(nbinx,nbiny);
928            }
929          }
930       }
931     }
932     //
933     Double_t xSecPercC = 0.;
934     if(hZDCCvsZEM->GetEntries()!=0){ 
935       xSecPercC = countPercC/hZDCCvsZEM->GetEntries();
936     }
937     else{
938       AliWarning("  Histogram hZDCCvsZEM from OCDB has no entries!!!");
939     }
940     // Ch. debug
941     //printf("  xSecPercC %1.4f  \n", xSecPercC);
942     
943     // ====> side A
944     TF1 *lineA = new TF1("lineA","[0]*x+[1]",0.,xHighEdge);
945     Float_t yA = (calibSumZN2[0]+calibSumZP2[0])/1000.;
946     lineA->SetParameter(0, yA/(x-origin));
947     lineA->SetParameter(1, -origin*yA/(x-origin));
948     //
949     // Ch. debug
950     //printf("  ***************** Side A \n");
951     //printf("  E_{ZEM} %1.4f, E_{ZDCA} %1.2f, TF1: %1.2f*x + %1.2f   ", x, yA,yA/(x-origin),-origin*yA/(x-origin));
952     //
953     Double_t countPercA=0;
954     Double_t xBinCenterA=0, yBinCenterA=0;
955     for(Int_t nbinx=1; nbinx<=hZDCAvsZEM->GetNbinsX(); nbinx++){
956       for(Int_t nbiny=1; nbiny<=hZDCAvsZEM->GetNbinsY(); nbiny++){
957          xBinCenterA = hZDCAvsZEM->GetXaxis()->GetBinCenter(nbinx);
958          yBinCenterA = hZDCAvsZEM->GetYaxis()->GetBinCenter(nbiny);
959          if(lineA->GetParameter(0)>0){
960            if(yBinCenterA < (lineA->GetParameter(0)*xBinCenterA + lineA->GetParameter(1))){
961              countPercA += hZDCAvsZEM->GetBinContent(nbinx,nbiny);
962            }
963          }
964          else{
965            if(yBinCenterA > (lineA->GetParameter(0)*xBinCenterA + lineA->GetParameter(1))){
966              countPercA += hZDCAvsZEM->GetBinContent(nbinx,nbiny);
967            }
968          }
969       }
970     }
971     //
972     Double_t xSecPercA = 0.;
973     if(hZDCAvsZEM->GetEntries()!=0){ 
974       xSecPercA = countPercA/hZDCAvsZEM->GetEntries();
975     }
976     else{
977       AliWarning("  Histogram hZDCAvsZEM from OCDB has no entries!!!");
978     }
979     // Ch. debug
980     //printf("  xSecPercA %1.4f  \n", xSecPercA);
981     
982     //  ******    Number of participants (from E_ZDC vs. E_ZEM correlation)
983     Int_t nPart=0, nPartC=0, nPartA=0;
984     Double_t nPartFrac=0., nPartFracC=0., nPartFracA=0.;
985     for(Int_t npbin=1; npbin<hNpartDist->GetNbinsX(); npbin++){
986       nPartFrac += (hNpartDist->GetBinContent(npbin))/(hNpartDist->GetEntries());
987       if((1.-nPartFrac) < xSecPerc){
988         nPart = (Int_t) hNpartDist->GetBinLowEdge(npbin);
989         // Ch. debug
990         //printf("  ***************** Summed ZDC info (sideA+side C) \n");
991         //printf("  nPartFrac %1.4f, nPart %d\n", nPartFrac, nPart);
992         break;
993       }
994     }
995     if(nPart<0) nPart=0;
996     //
997     for(Int_t npbin=1; npbin<hNpartDist->GetNbinsX(); npbin++){
998       nPartFracC += (hNpartDist->GetBinContent(npbin))/(hNpartDist->GetEntries());
999       if((1.-nPartFracC) < xSecPercC){
1000         nPartC = (Int_t) hNpartDist->GetBinLowEdge(npbin);
1001         // Ch. debug
1002         //printf("  ***************** Side C \n");
1003         //printf("  nPartFracC %1.4f, nPartC %d\n", nPartFracC, nPartC);
1004         break;
1005     }
1006     }
1007     if(nPartC<0) nPartC=0;
1008     //
1009     for(Int_t npbin=1; npbin<hNpartDist->GetNbinsX(); npbin++){
1010       nPartFracA += (hNpartDist->GetBinContent(npbin))/(hNpartDist->GetEntries());
1011       if((1.-nPartFracA) < xSecPercA){
1012         nPartA = (Int_t) hNpartDist->GetBinLowEdge(npbin);
1013         // Ch. debug
1014         //printf("  ***************** Side A \n");
1015         //printf("  nPartFracA %1.4f, nPartA %d\n\n", nPartFracA, nPartA);
1016         break;
1017       }
1018     }
1019     if(nPartA<0) nPartA=0;
1020     
1021     //  ******    Impact parameter (from E_ZDC vs. E_ZEM correlation)
1022     Float_t b=0, bC=0, bA=0;
1023     Double_t bFrac=0., bFracC=0., bFracA=0.;
1024     for(Int_t ibbin=1; ibbin<hbDist->GetNbinsX(); ibbin++){
1025       bFrac += (hbDist->GetBinContent(ibbin))/(hbDist->GetEntries());
1026       if((1.-bFrac) < xSecPerc){
1027         b = hbDist->GetBinLowEdge(ibbin);
1028         break;
1029       }
1030     }
1031     //
1032     for(Int_t ibbin=1; ibbin<hbDist->GetNbinsX(); ibbin++){
1033       bFracC += (hbDist->GetBinContent(ibbin))/(hbDist->GetEntries());
1034       if((1.-bFracC) < xSecPercC){
1035         bC = hbDist->GetBinLowEdge(ibbin);
1036         break;
1037       }
1038     }
1039     //
1040     for(Int_t ibbin=1; ibbin<hbDist->GetNbinsX(); ibbin++){
1041       bFracA += (hbDist->GetBinContent(ibbin))/(hNpartDist->GetEntries());
1042       if((1.-bFracA) < xSecPercA){
1043         bA = hbDist->GetBinLowEdge(ibbin);
1044         break;
1045       }
1046     }
1047
1048     //  ******  Number of spectator nucleons 
1049     Int_t nGenSpec=0, nGenSpecC=0, nGenSpecA=0;
1050     //
1051     nGenSpec = 416 - nPart;
1052     nGenSpecC = 416 - nPartC;
1053     nGenSpecA = 416 - nPartA;
1054     if(nGenSpec>416) nGenSpec=416; if(nGenSpec<0) nGenSpec=0;
1055     if(nGenSpecC>416) nGenSpecC=416; if(nGenSpecC<0) nGenSpecC=0;
1056     if(nGenSpecA>416) nGenSpecA=416; if(nGenSpecA<0) nGenSpecA=0;
1057   
1058     //  Ch. debug
1059     /*printf("\n\t AliZDCReconstructor -> calibSumZN1[0] %1.0f, calibSumZP1[0] %1.0f,"
1060         "  calibSumZN2[0] %1.0f, calibSumZP2[0] %1.0f, corrADCZEMHG %1.0f\n", 
1061         calibSumZN1[0],calibSumZP1[0],calibSumZN2[0],calibSumZP2[0],corrADCZEMHG);
1062     printf("\t AliZDCReconstructor -> nGenSpecLeft %d nGenSpecRight %d\n", 
1063         nGenSpecLeft, nGenSpecRight);
1064     printf("\t AliZDCReconstructor ->  NpartL %d,  NpartR %d,  b %1.2f fm\n\n",nPartTotLeft, nPartTotRight, impPar);
1065     */
1066   
1067     delete lineC;  delete lineA;
1068
1069   } // ONLY IF fIsCalibrationMB==kFALSE
1070   
1071   AliZDCReco* reco = new AliZDCReco(calibSumZN1, calibSumZP1, calibSumZN2, calibSumZP2, 
1072                   calibTowZN1, calibTowZP1, calibTowZN2, calibTowZP2, 
1073                   calibZEM1, calibZEM2, sPMRef1, sPMRef2,
1074                   nDetSpecNLeft, nDetSpecPLeft, nDetSpecNRight, nDetSpecPRight, 
1075                   nGenSpec, nGenSpecA, nGenSpecC, 
1076                   nPart, nPartA, nPartC, b, bA, bC,
1077                   recoFlag);
1078                     
1079   const Int_t kBufferSize = 4000;
1080   clustersTree->Branch("ZDC", "AliZDCReco", &reco, kBufferSize);
1081   reco->Print("");
1082   // write the output tree
1083   clustersTree->Fill();
1084 }
1085
1086 //_____________________________________________________________________________
1087 void AliZDCReconstructor::BuildRecoParam(Float_t ZDCC, Float_t ZDCA, Float_t ZEM) const
1088 {
1089   // Calculate RecoParam object for Pb-Pb data
1090   (fRecoParam->GethZDCvsZEM())->Fill(ZDCC+ZDCA, ZEM);
1091   (fRecoParam->GethZDCCvsZEM())->Fill(ZDCC, ZEM);
1092   (fRecoParam->GethZDCAvsZEM())->Fill(ZDCA, ZEM);
1093  
1094 }
1095
1096 //_____________________________________________________________________________
1097 void AliZDCReconstructor::FillZDCintoESD(TTree *clustersTree, AliESDEvent* esd) const
1098 {
1099   // fill energies and number of participants to the ESD
1100
1101   if(fIsCalibrationMB==kTRUE) WritePbPbRecoParamInOCDB();
1102   
1103   AliZDCReco reco;
1104   AliZDCReco* preco = &reco;
1105   clustersTree->SetBranchAddress("ZDC", &preco);
1106
1107   clustersTree->GetEntry(0);
1108   //
1109   AliESDZDC * esdzdc = esd->GetESDZDC();
1110   Float_t tZN1Ene[5], tZN2Ene[5], tZP1Ene[5], tZP2Ene[5];
1111   Float_t tZN1EneLR[5], tZN2EneLR[5], tZP1EneLR[5], tZP2EneLR[5];
1112   for(Int_t i=0; i<5; i++){
1113      tZN1Ene[i] = reco.GetZN1HREnTow(i);
1114      tZN2Ene[i] = reco.GetZN2HREnTow(i);
1115      tZP1Ene[i] = reco.GetZP1HREnTow(i);
1116      tZP2Ene[i] = reco.GetZP2HREnTow(i);
1117      //
1118      tZN1EneLR[i] = reco.GetZN1LREnTow(i);
1119      tZN2EneLR[i] = reco.GetZN2LREnTow(i);
1120      tZP1EneLR[i] = reco.GetZP1LREnTow(i);
1121      tZP2EneLR[i] = reco.GetZP2LREnTow(i);
1122   }
1123   //
1124   esdzdc->SetZN1TowerEnergy(tZN1Ene);
1125   esdzdc->SetZN2TowerEnergy(tZN2Ene);
1126   esdzdc->SetZP1TowerEnergy(tZP1Ene);
1127   esdzdc->SetZP2TowerEnergy(tZP2Ene);
1128   //
1129   esdzdc->SetZN1TowerEnergyLR(tZN1EneLR);
1130   esdzdc->SetZN2TowerEnergyLR(tZN2EneLR);
1131   esdzdc->SetZP1TowerEnergyLR(tZP1EneLR);
1132   esdzdc->SetZP2TowerEnergyLR(tZP2EneLR);
1133   // 
1134   Int_t nPart  = reco.GetNParticipants();
1135   Int_t nPartA = reco.GetNPartSideA();
1136   Int_t nPartC = reco.GetNPartSideC();
1137   Double_t b  = reco.GetImpParameter();
1138   Double_t bA = reco.GetImpParSideA();
1139   Double_t bC = reco.GetImpParSideC();
1140   UInt_t recoFlag = reco.GetRecoFlag();
1141   //
1142   esd->SetZDC(reco.GetZN1HREnergy(), reco.GetZP1HREnergy(), reco.GetZEM1HRsignal(), 
1143               reco.GetZEM2HRsignal(), reco.GetZN2HREnergy(), reco.GetZP2HREnergy(), 
1144               nPart, nPartA, nPartC, b, bA, bC, recoFlag);
1145   
1146 }
1147
1148 //_____________________________________________________________________________
1149 AliCDBStorage* AliZDCReconstructor::SetStorage(const char *uri) 
1150 {
1151   // Setting the storage
1152
1153   Bool_t deleteManager = kFALSE;
1154   
1155   AliCDBManager *manager = AliCDBManager::Instance();
1156   AliCDBStorage *defstorage = manager->GetDefaultStorage();
1157   
1158   if(!defstorage || !(defstorage->Contains("ZDC"))){ 
1159      AliWarning("No default storage set or default storage doesn't contain ZDC!");
1160      manager->SetDefaultStorage(uri);
1161      deleteManager = kTRUE;
1162   }
1163  
1164   AliCDBStorage *storage = manager->GetDefaultStorage();
1165
1166   if(deleteManager){
1167     AliCDBManager::Instance()->UnsetDefaultStorage();
1168     defstorage = 0;   // the storage is killed by AliCDBManager::Instance()->Destroy()
1169   }
1170
1171   return storage; 
1172 }
1173
1174 //_____________________________________________________________________________
1175 AliZDCPedestals* AliZDCReconstructor::GetPedestalData() const
1176 {
1177
1178   // Getting pedestal calibration object for ZDC set
1179
1180   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
1181   if(!entry) AliFatal("No calibration data loaded!");  
1182
1183   AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*>  (entry->GetObject());
1184   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
1185
1186   return calibdata;
1187 }
1188
1189 //_____________________________________________________________________________
1190 AliZDCEnCalib* AliZDCReconstructor::GetEnergyCalibData() const
1191 {
1192
1193   // Getting energy and equalization calibration object for ZDC set
1194
1195   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/EnergyCalib");
1196   if(!entry) AliFatal("No calibration data loaded!");  
1197
1198   AliZDCEnCalib *calibdata = dynamic_cast<AliZDCEnCalib*> (entry->GetObject());
1199   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
1200
1201   return calibdata;
1202 }
1203
1204 //_____________________________________________________________________________
1205 AliZDCTowerCalib* AliZDCReconstructor::GetTowerCalibData() const
1206 {
1207
1208   // Getting energy and equalization calibration object for ZDC set
1209
1210   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/TowerCalib");
1211   if(!entry) AliFatal("No calibration data loaded!");  
1212
1213   AliZDCTowerCalib *calibdata = dynamic_cast<AliZDCTowerCalib*> (entry->GetObject());
1214   if(!calibdata)  AliFatal("Wrong calibration object in calibration  file!");
1215
1216   return calibdata;
1217 }
1218
1219 //_____________________________________________________________________________
1220 AliZDCRecoParampp* AliZDCReconstructor::GetppRecoParamFromOCDB() const
1221 {
1222
1223   // Getting reconstruction parameters from OCDB
1224
1225   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/RecoParampp");
1226   if(!entry) AliFatal("No RecoParam data found in OCDB!");  
1227   
1228   AliZDCRecoParampp *param = dynamic_cast<AliZDCRecoParampp*> (entry->GetObject());
1229   if(!param)  AliFatal("No RecoParam object in OCDB entry!");
1230   
1231   return param;
1232
1233 }
1234
1235 //_____________________________________________________________________________
1236 AliZDCRecoParamPbPb* AliZDCReconstructor::GetPbPbRecoParamFromOCDB() const
1237 {
1238
1239   // Getting reconstruction parameters from OCDB
1240
1241   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/RecoParamPbPb");
1242   if(!entry) AliFatal("No RecoParam data found in OCDB!");  
1243   
1244   AliZDCRecoParamPbPb *param = dynamic_cast<AliZDCRecoParamPbPb*> (entry->GetObject());
1245   if(!param)  AliFatal("No RecoParam object in OCDB entry!");
1246   
1247   return param;
1248
1249 }
1250
1251 //_____________________________________________________________________________
1252 void AliZDCReconstructor::WritePbPbRecoParamInOCDB() const
1253 {
1254
1255   // Writing Pb-Pb reconstruction parameters from OCDB
1256
1257   AliCDBManager *man = AliCDBManager::Instance();
1258   AliCDBMetaData *md= new AliCDBMetaData();
1259   md->SetResponsible("Chiara Oppedisano");
1260   md->SetComment("ZDC Pb-Pb reconstruction parameters");
1261   md->SetObjectClassName("AliZDCRecoParamPbPb");
1262   AliCDBId id("ZDC/Calib/RecoParamPbPb",fNRun,AliCDBRunRange::Infinity());
1263   man->Put(fRecoParam, id, md);
1264
1265 }
1266