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