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