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