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