1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for ZDC reconstruction //
22 ///////////////////////////////////////////////////////////////////////////////
27 #include "AliRunLoader.h"
28 #include "AliRawReader.h"
30 #include "AliZDCDigit.h"
31 #include "AliZDCRawStream.h"
32 #include "AliZDCReco.h"
33 #include "AliZDCReconstructor.h"
34 #include "AliZDCCalibData.h"
37 ClassImp(AliZDCReconstructor)
40 //_____________________________________________________________________________
41 AliZDCReconstructor:: AliZDCReconstructor()
43 // **** Default constructor
46 // --- Number of generated spectator nucleons and impact parameter
47 // --------------------------------------------------------------------------------------------------
48 // [1] ### Results from a new production -> 0<b<18 fm (Apr 2002)
49 // Fit results for neutrons (Nspectator n true vs. EZN)
50 fZNCen = new TF1("fZNCen",
51 "(-2.287920+sqrt(2.287920*2.287920-4*(-0.007629)*(11.921710-x)))/(2*(-0.007629))",0.,164.);
52 fZNPer = new TF1("fZNPer",
53 "(-37.812280-sqrt(37.812280*37.812280-4*(-0.190932)*(-1709.249672-x)))/(2*(-0.190932))",0.,164.);
54 // Fit results for protons (Nspectator p true vs. EZP)
55 fZPCen = new TF1("fZPCen",
56 "(-1.321353+sqrt(1.321353*1.321353-4*(-0.007283)*(3.550697-x)))/(2*(-0.007283))",0.,60.);
57 fZPPer = new TF1("fZPPer",
58 "(-42.643308-sqrt(42.643308*42.643308-4*(-0.310786)*(-1402.945615-x)))/(2*(-0.310786))",0.,60.);
59 // Fit results for total number of spectators (Nspectators true vs. EZDC)
60 fZDCCen = new TF1("fZDCCen",
61 "(-1.934991+sqrt(1.934991*1.934991-4*(-0.004080)*(15.111124-x)))/(2*(-0.004080))",0.,225.);
62 fZDCPer = new TF1("fZDCPer",
63 "(-34.380639-sqrt(34.380639*34.380639-4*(-0.104251)*(-2612.189017-x)))/(2*(-0.104251))",0.,225.);
64 // --------------------------------------------------------------------------------------------------
65 // Fit results for b (b vs. EZDC)
66 // [2] ### Results from a new production -> 0<b<18 fm (Apr 2002)
67 fbCen = new TF1("fbCen","-0.056923+0.079703*x-0.0004301*x*x+0.000001366*x*x*x",0.,220.);
68 fbPer = new TF1("fbPer","17.943998-0.046846*x+0.000074*x*x",0.,220.);
69 // --------------------------------------------------------------------------------------------------
70 // Evaluating Nspectators and b from ZEM energy
71 // [2] ### Results from a new production -> 0<b<18 fm (Apr 2002)
72 fZEMn = new TF1("fZEMn","126.2-0.05399*x+0.000005679*x*x",0.,4000.);
73 fZEMp = new TF1("fZEMp","82.49-0.03611*x+0.00000385*x*x",0.,4000.);
74 fZEMsp = new TF1("fZEMsp","208.7-0.09006*x+0.000009526*x*x",0.,4000.);
75 fZEMb = new TF1("fZEMb","16.06-0.01633*x+1.44e-5*x*x-6.778e-9*x*x*x+1.438e-12*x*x*x*x-1.112e-16*x*x*x*x*x",0.,4000.);
78 fStorage = SetStorage("local://$ALICE_ROOT");
80 // Get calibration data
82 fCalibData = GetCalibData(runNumber);
85 //_____________________________________________________________________________
86 AliZDCReconstructor::AliZDCReconstructor(const AliZDCReconstructor&
88 AliReconstructor(reconstructor)
92 Fatal("AliZDCReconstructor", "copy constructor not implemented");
95 //_____________________________________________________________________________
96 AliZDCReconstructor& AliZDCReconstructor::operator =
97 (const AliZDCReconstructor& /*reconstructor*/)
99 // assignment operator
101 Fatal("operator =", "assignment operator not implemented");
105 //_____________________________________________________________________________
106 AliZDCReconstructor::~AliZDCReconstructor()
125 //_____________________________________________________________________________
126 void AliZDCReconstructor::Reconstruct(AliRunLoader* runLoader) const
128 // *** Local ZDC reconstruction for digits
131 for(Int_t jj=0; jj<47; jj++) meanPed[jj] = fCalibData->GetMeanPed(jj);
133 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
135 loader->LoadDigits("read");
136 loader->LoadRecPoints("recreate");
138 AliZDCDigit* pdigit = &digit;
141 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
142 runLoader->GetEvent(iEvent);
145 loader->LoadDigits();
146 TTree* treeD = loader->TreeD();
147 if (!treeD) continue;
148 treeD->SetBranchAddress("ZDC", &pdigit);
151 Float_t zn1corr=0, zp1corr=0, zn2corr=0, zp2corr=0, zemcorr=0;
152 for (Int_t iDigit = 0; iDigit < treeD->GetEntries(); iDigit++) {
153 treeD->GetEntry(iDigit);
154 if (!pdigit) continue;
156 if(digit.GetSector(0) == 1)
157 zn1corr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)]); // high gain ZN1 ADCs
158 else if(digit.GetSector(0) == 2)
159 zp1corr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+10]); // high gain ZP1 ADCs
160 else if(digit.GetSector(0) == 3){
161 if(digit.GetSector(1)==1)
162 zemcorr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+20]); // high gain ZEM1 ADCs
163 else if(digit.GetSector(1)==2)
164 zemcorr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+22]); // high gain ZEM2 ADCs
166 else if(digit.GetSector(0) == 4)
167 zn2corr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+24]); // high gain ZN2 ADCs
168 else if(digit.GetSector(0) == 5)
169 zp2corr += (Float_t) (digit.GetADCValue(0)-meanPed[digit.GetSector(1)+34]); // high gain ZP2 ADCs
171 if(zn1corr<0) zn1corr=0;
172 if(zp1corr<0) zp1corr=0;
173 if(zn2corr<0) zn2corr=0;
174 if(zp2corr<0) zp2corr=0;
175 if(zemcorr<0) zemcorr=0;
177 // reconstruct the event
178 //printf("\n \t ZDCReco from digits-> Ev.#%d ZN = %.0f, ZP = %.0f, ZEM = %.0f\n",iEvent,zncorr,zpcorr,zemcorr);
179 ReconstructEvent(loader, zn1corr, zp1corr, zemcorr, zn2corr, zp2corr);
182 loader->UnloadDigits();
183 loader->UnloadRecPoints();
186 //_____________________________________________________________________________
187 void AliZDCReconstructor::Reconstruct(AliRunLoader* runLoader,
188 AliRawReader* rawReader) const
190 // *** Local ZDC reconstruction for raw data
193 for(Int_t jj=0; jj<47; jj++) meanPed[jj] = fCalibData->GetMeanPed(jj);
195 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
197 loader->LoadRecPoints("recreate");
200 while (rawReader->NextEvent()) {
201 runLoader->GetEvent(iEvent++);
203 // loop over raw data digits
204 Float_t zn1corr=0, zp1corr=0, zn2corr=0, zp2corr=0,zemcorr=0;
205 AliZDCRawStream digit(rawReader);
206 while (digit.Next()) {
207 if(digit.IsADCDataWord()){
208 if(digit.GetADCGain() == 0){
209 if(digit.GetSector(0) == 1)
210 zn1corr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)]); // high gain ZN1 ADCs;
211 else if(digit.GetSector(0) == 2)
212 zp1corr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+10]); // high gain ZP1 ADCs;
213 else if(digit.GetSector(0) == 3)
214 if(digit.GetSector(1)==1)
215 zemcorr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+20]); // high gain ZEM1 ADCs
216 else if(digit.GetSector(1)==2)
217 zemcorr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+22]); // high gain ZEM2 ADCs
218 else if(digit.GetSector(0) == 4)
219 zn2corr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+24]); // high gain ZN2 ADCs;
220 else if(digit.GetSector(0) == 5)
221 zp2corr += (Float_t) (digit.GetADCValue()-meanPed[digit.GetSector(1)+34]); // high gain ZP2 ADCs;
225 if(zn1corr<0) zn1corr=0;
226 if(zp1corr<0) zp1corr=0;
227 if(zn2corr<0) zn2corr=0;
228 if(zp2corr<0) zp2corr=0;
229 if(zemcorr<0) zemcorr=0;
231 // reconstruct the event
232 //printf("\n\t ZDCReco from raw-> Ev.#%d ZN = %.0f, ZP = %.0f, ZEM = %.0f\n",iEvent,zncorr,zpcorr,zemcorr);
233 ReconstructEvent(loader, zn1corr, zp1corr, zemcorr, zn2corr, zp2corr);
236 loader->UnloadRecPoints();
239 //_____________________________________________________________________________
240 void AliZDCReconstructor::ReconstructEvent(AliLoader* loader, Float_t zn1corr,
241 Float_t zp1corr, Float_t zemcorr, Float_t zn2corr, Float_t zp2corr) const
243 // ***** Reconstruct one event
245 // --- ADCchannel -> photoelectrons
246 // NB-> PM gain = 10^(5), ADC resolution = 6.4*10^(-7)
247 // Move to V965 (E.S.,15/09/04) NB-> PM gain = 10^(5), ADC resolution = 8*10^(-7)
248 Float_t zn1phe, zp1phe, zemphe, zn2phe, zp2phe, convFactor = 0.08;
249 zn1phe = zn1corr/convFactor;
250 zp1phe = zp1corr/convFactor;
251 zemphe = zemcorr/convFactor;
252 zn2phe = zn2corr/convFactor;
253 zp2phe = zp2corr/convFactor;
254 //if AliDebug(1,Form("\n znphe = %f, zpphe = %f, zemphe = %f\n",znphe, zpphe, zemphe);
256 // --- Energy calibration
257 // Conversion factors for hadronic ZDCs goes from phe yield to TRUE
258 // incident energy (conversion from GeV to TeV is included); while for EM
259 // calos conversion is from light yield to detected energy calculated by
260 // GEANT NB -> ZN and ZP conversion factors are constant since incident
261 // spectators have all the same energy, ZEM energy is obtained through a
262 // fit over the whole range of incident particle energies
263 // (obtained with full HIJING simulations)
264 Float_t zn1energy, zp1energy, zemenergy, zdc1energy, zn2energy, zp2energy, zdc2energy;
265 Float_t zn1phexTeV=329., zp1phexTeV=369., zn2phexTeV=329., zp2phexTeV=369.;
266 zn1energy = zn1phe/zn1phexTeV;
267 zp1energy = zp1phe/zp1phexTeV;
268 zdc1energy = zn1energy+zp1energy;
269 zn2energy = zn2phe/zn2phexTeV;
270 zp2energy = zp2phe/zp2phexTeV;
271 zdc2energy = zn2energy+zp2energy;
272 zemenergy = -4.81+0.3238*zemphe;
273 if(zemenergy<0) zemenergy=0;
274 // if AliDebug(1,Form(" znenergy = %f TeV, zpenergy = %f TeV, zdcenergy = %f GeV, "
275 // "\n zemenergy = %f TeV\n", znenergy, zpenergy,
276 // zdcenergy, zemenergy);
278 // if AliDebug(1,Form("\n\n ### ATTENZIONE!!! -> ev# %d: znenergy = %f TeV, zpenergy = %f TeV, zdcenergy = %f GeV, "
279 // " zemenergy = %f TeV\n\n", fMerger->EvNum(), znenergy, zpenergy, zdcenergy, zemenergy);
281 // --- Number of detected spectator nucleons
282 // *** N.B. -> It works only in Pb-Pb
283 Int_t nDetSpecNLeft, nDetSpecPLeft, nDetSpecNRight, nDetSpecPRight;
284 nDetSpecNLeft = (Int_t) (zn1energy/2.760);
285 nDetSpecPLeft = (Int_t) (zp1energy/2.760);
286 nDetSpecNRight = (Int_t) (zn2energy/2.760);
287 nDetSpecPRight = (Int_t) (zp2energy/2.760);
289 // --- Number of generated spectator nucleons (from HIJING parameterization)
290 // *** N.B. -> Only one side!!!
291 Int_t nGenSpecN=0, nGenSpecP=0, nGenSpec=0;
293 // Cut value for Ezem (GeV)
294 // ### Results from production -> 0<b<18 fm (Apr 2002)
295 Float_t eZEMCut = 420.;
296 Float_t deltaEZEMSup = 690.;
297 Float_t deltaEZEMInf = 270.;
298 if(zemenergy > (eZEMCut+deltaEZEMSup)){
299 nGenSpecN = (Int_t) (fZNCen->Eval(zn1energy));
300 nGenSpecP = (Int_t) (fZPCen->Eval(zp1energy));
301 nGenSpec = (Int_t) (fZDCCen->Eval(zdc1energy));
302 impPar = fbCen->Eval(zdc1energy);
304 else if(zemenergy < (eZEMCut-deltaEZEMInf)){
305 nGenSpecN = (Int_t) (fZNPer->Eval(zn1energy));
306 nGenSpecP = (Int_t) (fZPPer->Eval(zp1energy));
307 nGenSpec = (Int_t) (fZDCPer->Eval(zdc1energy));
308 impPar = fbPer->Eval(zdc1energy);
310 else if(zemenergy >= (eZEMCut-deltaEZEMInf) && zemenergy <= (eZEMCut+deltaEZEMSup)){
311 nGenSpecN = (Int_t) (fZEMn->Eval(zemenergy));
312 nGenSpecP = (Int_t) (fZEMp->Eval(zemenergy));
313 nGenSpec = (Int_t)(fZEMsp->Eval(zemenergy));
314 impPar = fZEMb->Eval(zemenergy);
316 // ### Results from production -> 0<b<18 fm (Apr 2002)
317 if(zn1energy>162.) nGenSpecN = (Int_t) (fZEMn->Eval(zemenergy));
318 if(zp1energy>59.75) nGenSpecP = (Int_t) (fZEMp->Eval(zemenergy));
319 if(zdc1energy>221.5) nGenSpec = (Int_t)(fZEMsp->Eval(zemenergy));
320 if(zdc1energy>220.) impPar = fZEMb->Eval(zemenergy);
322 if(nGenSpecN>125) nGenSpecN=125;
323 else if(nGenSpecN<0) nGenSpecN=0;
324 if(nGenSpecP>82) nGenSpecP=82;
325 else if(nGenSpecP<0) nGenSpecP=0;
326 if(nGenSpec>207) nGenSpec=207;
327 else if(nGenSpec<0) nGenSpec=0;
329 // --- Number of generated participants (from HIJING parameterization)
330 Int_t nPart, nPartTot;
331 nPart = 207-nGenSpecN-nGenSpecP;
332 nPartTot = 207-nGenSpec;
333 //printf("\t ZDCeventReco-> ZNEn = %.0f GeV, ZPEn = %.0f GeV, ZEMEn = %.0f GeV\n",
334 // znenergy, zpenergy, zemenergy);
336 // create the output tree
337 loader->MakeTree("R");
338 TTree* treeR = loader->TreeR();
339 AliZDCReco reco(zn1energy, zp1energy, zdc1energy, zemenergy,
340 zn2energy, zp2energy, zdc2energy,
341 nDetSpecNLeft, nDetSpecPLeft, nDetSpecNRight, nDetSpecPRight,
342 nGenSpecN, nGenSpecP, nGenSpec,nPartTot, impPar);
343 AliZDCReco* preco = &reco;
344 const Int_t kBufferSize = 4000;
345 treeR->Branch("ZDC", "AliZDCReco", &preco, kBufferSize);
347 // write the output tree
349 loader->WriteRecPoints("OVERWRITE");
352 //_____________________________________________________________________________
353 void AliZDCReconstructor::FillESD(AliRunLoader* runLoader,
356 // fill energies and number of participants to the ESD
358 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
360 loader->LoadRecPoints();
362 TTree* treeR = loader->TreeR();
365 AliZDCReco* preco = &reco;
366 treeR->SetBranchAddress("ZDC", &preco);
369 esd->SetZDC(reco.GetZN1energy(), reco.GetZP1energy(), reco.GetZEMenergy(),
370 reco.GetZN2energy(), reco.GetZP2energy(), reco.GetNPart());
372 loader->UnloadRecPoints();
375 //_____________________________________________________________________________
376 AliCDBStorage* AliZDCReconstructor::SetStorage(const char *uri)
378 //printf("\n\t AliZDCReconstructor::SetStorage \n");
380 Bool_t deleteManager = kFALSE;
382 AliCDBManager *manager = AliCDBManager::Instance();
383 AliCDBStorage *defstorage = manager->GetDefaultStorage();
385 if(!defstorage || !(defstorage->Contains("ZDC"))){
386 AliWarning("No default storage set or default storage doesn't contain ZDC!");
387 manager->SetDefaultStorage(uri);
388 deleteManager = kTRUE;
391 AliCDBStorage *storage = manager->GetDefaultStorage();
394 AliCDBManager::Instance()->UnsetDefaultStorage();
395 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
401 //_____________________________________________________________________________
402 AliZDCCalibData* AliZDCReconstructor::GetCalibData(int runNumber) const
405 //printf("\n\t AliZDCReconstructor::GetCalibData \n");
407 AliCDBEntry *entry = fStorage->Get("ZDC/Calib/Data",runNumber);
408 AliZDCCalibData *calibdata = (AliZDCCalibData*) entry->GetObject();
410 if (!calibdata) AliWarning("No calibration data from calibration database !");