]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ZDC/AliZDCDigitizer.cxx
Adding parameters for LHC13d and LHC13e (muon_pass2) (Diego)
[u/mrichter/AliRoot.git] / ZDC / AliZDCDigitizer.cxx
CommitLineData
8309c1ab 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// ZDC digitizer class //
21// //
22///////////////////////////////////////////////////////////////////////////////
23
3d7e0d0c 24// --- Standard libraries
8309c1ab 25#include <stdlib.h>
3d7e0d0c 26#include <stdio.h>
8309c1ab 27
28// --- ROOT system
29#include <TTree.h>
30#include <TFile.h>
31#include <TNtuple.h>
32#include <TRandom.h>
33
34// --- AliRoot header files
8309c1ab 35#include "AliRun.h"
36#include "AliHeader.h"
37#include "AliGenHijingEventHeader.h"
07e38cef 38#include "AliGenCocktailEventHeader.h"
f21fc003 39#include "AliDigitizationInput.h"
8309c1ab 40#include "AliRunLoader.h"
c93255fe 41#include "AliLoader.h"
12434381 42#include "AliGRPObject.h"
78d18275 43#include "AliCDBManager.h"
48642b09 44#include "AliCDBEntry.h"
8309c1ab 45#include "AliZDCSDigit.h"
46#include "AliZDCDigit.h"
47#include "AliZDCFragment.h"
a18a0543 48#include "AliZDCv3.h"
8309c1ab 49#include "AliZDCDigitizer.h"
8a2624cc 50
51class AliCDBStorage;
6024ec85 52class AliZDCPedestals;
8309c1ab 53
54ClassImp(AliZDCDigitizer)
55
56
57//____________________________________________________________________________
a718c993 58AliZDCDigitizer::AliZDCDigitizer() :
59 fIsCalibration(0),
fd9afd60 60 fIsSignalInADCGate(kFALSE),
61 fFracLostSignal(0.),
a718c993 62 fPedData(0),
2d9c70ab 63 fSpectators2Track(kFALSE),
f2e2aa97 64 fBeamEnergy(0.),
5d24ec95 65 fBeamType(""),
66 fIspASystem(kFALSE),
67 fIsRELDISgen(kFALSE)
8309c1ab 68{
698e394d 69 // Default constructor
83534754 70 for(Int_t i=0; i<2; i++) fADCRes[i]=0.;
8309c1ab 71}
72
73//____________________________________________________________________________
f21fc003 74AliZDCDigitizer::AliZDCDigitizer(AliDigitizationInput* digInput):
75 AliDigitizer(digInput),
a718c993 76 fIsCalibration(0), //By default the simulation doesn't create calib. data
fd9afd60 77 fIsSignalInADCGate(kFALSE),
78 fFracLostSignal(0.),
a718c993 79 fPedData(GetPedData()),
2d9c70ab 80 fSpectators2Track(kFALSE),
f2e2aa97 81 fBeamEnergy(0.),
5d24ec95 82 fBeamType(""),
83 fIspASystem(kFALSE),
84 fIsRELDISgen(kFALSE)
8309c1ab 85{
78d18275 86 // Get calibration data
d27e20dd 87 if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
edb4e893 88 for(Int_t i=0; i<5; i++){
698e394d 89 for(Int_t j=0; j<5; j++)
90 fPMGain[i][j] = 0.;
91 }
8309c1ab 92}
93
94//____________________________________________________________________________
95AliZDCDigitizer::~AliZDCDigitizer()
96{
97// Destructor
a718c993 98// Not implemented
8309c1ab 99}
100
101
cc2abffd 102//____________________________________________________________________________
103AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
a718c993 104 AliDigitizer(),
105 fIsCalibration(digitizer.fIsCalibration),
fd9afd60 106 fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
107 fFracLostSignal(digitizer.fFracLostSignal),
a718c993 108 fPedData(digitizer.fPedData),
2d9c70ab 109 fSpectators2Track(digitizer.fSpectators2Track),
f2e2aa97 110 fBeamEnergy(digitizer.fBeamEnergy),
5d24ec95 111 fBeamType(digitizer.fBeamType),
112 fIspASystem(digitizer.fIspASystem),
113 fIsRELDISgen(digitizer.fIsRELDISgen)
cc2abffd 114{
115 // Copy constructor
116
edb4e893 117 for(Int_t i=0; i<5; i++){
cc2abffd 118 for(Int_t j=0; j<5; j++){
119 fPMGain[i][j] = digitizer.fPMGain[i][j];
120 }
121 }
122 for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
cc2abffd 123
12434381 124
cc2abffd 125}
126
8309c1ab 127//____________________________________________________________________________
128Bool_t AliZDCDigitizer::Init()
129{
48642b09 130 // Initialize the digitizer
12434381 131
5d24ec95 132
133 //printf(" **** Initializing AliZDCDigitizer fBeamEnergy = %1.0f GeV\n\n", fBeamEnergy);
12434381 134 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
83534754 135 if(!entry) AliFatal("No calibration data loaded!");
12434381 136 AliGRPObject* grpData = 0x0;
137 if(entry){
138 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
139 if(m){
140 //m->Print();
141 grpData = new AliGRPObject();
142 grpData->ReadValuesFromMap(m);
143 }
144 else{
145 grpData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
146 }
147 entry->SetOwner(0);
148 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
149 }
83534754 150 if(!grpData){
151 AliError("No GRP entry found in OCDB! \n ");
152 return kFALSE;
153 }
12434381 154
5d24ec95 155 fBeamType = grpData->GetBeamType();
156 if(fBeamType==AliGRPObject::GetInvalidString()){
12434381 157 AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
158 }
159
5d24ec95 160 // If beam energy is not set from Config.C (RELDIS / spectator generators)
161 if(fBeamEnergy<0.01){
162 fBeamEnergy = grpData->GetBeamEnergy();
163 if(fBeamEnergy==AliGRPObject::GetInvalidFloat()){
164 AliWarning("GRP/GRP/Data entry: missing value for the beam energy ! Using 0.");
165 AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
166 fBeamEnergy = 0.;
a736fd93 167 }
12434381 168 }
f2e2aa97 169
f6c61b40 170 /*if(fIspASystem){
5d24ec95 171 fBeamType = "p-A";
172 AliInfo(" AliZDCDigitizer -> Manually setting beam type to p-A\n");
f6c61b40 173 }*/
5d24ec95 174
175 // Setting beam type for spectator generator and RELDIS generator
176 if(((fBeamType.CompareTo("UNKNOWN")) == 0) || fIsRELDISgen){
177 fBeamType = "A-A";
178 AliInfo(" AliZDCDigitizer -> Manually setting beam type to A-A\n");
179 }
180 printf("\t AliZDCDigitizer -> beam type %s - beam energy = %f GeV\n", fBeamType.Data(), fBeamEnergy);
181
db247fc3 182 if(fBeamEnergy>0.1){
183 ReadPMTGains();
f6c61b40 184 //CalculatePMTGains();
db247fc3 185 }
186 else{
187 AliWarning("\n Beam energy is 0 -> ZDC PMT gains can't be set -> NO ZDC DIGITS!!!\n");
188 }
189
698e394d 190 // ADC Caen V965
7f73eb6b 191 fADCRes[0] = 0.0000008; // ADC Resolution high gain: 200 fC/adcCh
192 fADCRes[1] = 0.0000064; // ADC Resolution low gain: 25 fC/adcCh
8309c1ab 193
194 return kTRUE;
195}
196
197//____________________________________________________________________________
f21fc003 198void AliZDCDigitizer::Digitize(Option_t* /*option*/)
8309c1ab 199{
48642b09 200 // Execute digitization
8309c1ab 201
d79f8d50 202 // ------------------------------------------------------------
f91d1e35 203 // !!! 2nd ZDC set added
204 // *** 1st 3 arrays are digits from REAL (simulated) hits
205 // *** last 2 are copied from simulated digits
fd9afd60 206 // --- pm[0][...] = light in ZN side C [C, Q1, Q2, Q3, Q4]
207 // --- pm[1][...] = light in ZP side C [C, Q1, Q2, Q3, Q4]
48642b09 208 // --- pm[2][...] = light in ZEM [x, 1, 2, x, x]
fd9afd60 209 // --- pm[3][...] = light in ZN side A [C, Q1, Q2, Q3, Q4] ->NEW!
210 // --- pm[4][...] = light in ZP side A [C, Q1, Q2, Q3, Q4] ->NEW!
d79f8d50 211 // ------------------------------------------------------------
212 Float_t pm[5][5];
8574b308 213 for(Int_t iSector1=0; iSector1<5; iSector1++)
214 for(Int_t iSector2=0; iSector2<5; iSector2++){
8309c1ab 215 pm[iSector1][iSector2] = 0;
216 }
d79f8d50 217
218 // ------------------------------------------------------------
219 // ### Out of time ADC added (22 channels)
220 // --- same codification as for signal PTMs (see above)
221 // ------------------------------------------------------------
3d7ab74d 222 // Float_t pmoot[5][5];
223 // for(Int_t iSector1=0; iSector1<5; iSector1++)
224 // for(Int_t iSector2=0; iSector2<5; iSector2++){
225 // pmoot[iSector1][iSector2] = 0;
226 // }
8309c1ab 227
228 // impact parameter and number of spectators
131b5d31 229 Float_t impPar = 0;
fd9afd60 230 Int_t specNTarg = 0, specPTarg = 0;
231 Int_t specNProj = 0, specPProj = 0;
a18a0543 232 Float_t signalTime0 = 0.;
8309c1ab 233
234 // loop over input streams
f21fc003 235 for(Int_t iInput = 0; iInput<fDigInput->GetNinputs(); iInput++){
8309c1ab 236
237 // get run loader and ZDC loader
238 AliRunLoader* runLoader =
f21fc003 239 AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
8309c1ab 240 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
8574b308 241 if(!loader) continue;
8309c1ab 242
243 // load sdigits
244 loader->LoadSDigits();
245 TTree* treeS = loader->TreeS();
8574b308 246 if(!treeS) continue;
8309c1ab 247 AliZDCSDigit sdigit;
248 AliZDCSDigit* psdigit = &sdigit;
249 treeS->SetBranchAddress("ZDC", &psdigit);
250
251 // loop over sdigits
8574b308 252 for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
8309c1ab 253 treeS->GetEntry(iSDigit);
7f73eb6b 254 //
8574b308 255 if(!psdigit) continue;
256 if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
8309c1ab 257 AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
258 sdigit.GetSector(0), sdigit.GetSector(1)));
259 continue;
260 }
fd9afd60 261 // Checking if signal is inside ADC gate
262 if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
263 else{
264 // Assuming a signal lenght of 20 ns, signal is in gate if
265 // signal ENDS in signalTime0+50. -> sdigit.GetTrackTime()+20<=signalTime0+50.
266 if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
267 if(sdigit.GetTrackTime()>signalTime0+30.){
268 fIsSignalInADCGate = kFALSE;
5d24ec95 269 // Vedi quaderno per spiegazione approx. usata nel calcolo della fraz. di segnale perso
fd9afd60 270 fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
271 }
272 }
698e394d 273
fd9afd60 274 Float_t sdSignal = sdigit.GetLightPM();
275 if(fIsSignalInADCGate == kFALSE){
7e858b1d 276 AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
0fba8107 277 sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
fd9afd60 278 sdSignal = (1-fFracLostSignal)*sdSignal;
279 }
280
48642b09 281 pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
698e394d 282 //Ch. debug
283 /*printf("\t Detector %d, Tower %d -> pm[%d][%d] = %.0f \n",
48642b09 284 sdigit.GetSector(0), sdigit.GetSector(1),sdigit.GetSector(0)-1,
698e394d 285 sdigit.GetSector(1), pm[sdigit.GetSector(0)-1][sdigit.GetSector(1)]);
878fa0b3 286 */
fd9afd60 287
8309c1ab 288 }
289
8309c1ab 290 loader->UnloadSDigits();
291
292 // get the impact parameter and the number of spectators in case of hijing
8574b308 293 if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
d3b3a3b2 294 AliHeader* header = runLoader->GetHeader();
8574b308 295 if(!header) continue;
8309c1ab 296 AliGenEventHeader* genHeader = header->GenEventHeader();
8574b308 297 if(!genHeader) continue;
07e38cef 298 AliGenHijingEventHeader *hijingHeader = 0;
299 if(genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (genHeader);
300 else if(genHeader->InheritsFrom(AliGenCocktailEventHeader::Class())){
301 TList* listOfHeaders = ((AliGenCocktailEventHeader*) genHeader)->GetHeaders();
c47d6938 302 if(listOfHeaders){
303 hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (listOfHeaders->FindObject("Hijing"));
304 if(!hijingHeader) hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (listOfHeaders->FindObject("HijingpPbPb"));
305 if(!hijingHeader) hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (listOfHeaders->FindObject("HijingpPb"));
306 if(!hijingHeader) hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (listOfHeaders->FindObject("HijingpPb_0"));
307 }
01e01248 308 else{
309 AliWarning(" No list of headers from generator -> skipping event\n");
310 continue;
311 }
07e38cef 312 }
c47d6938 313 if(!hijingHeader){
c47d6938 314 printf(" No HIJING header found in list of headers from generator\n");
315 //continue;
316 }
07e38cef 317
8fced614 318 if(fSpectators2Track==kTRUE){
319 impPar = hijingHeader->ImpactParameter();
07e38cef 320 specNProj = hijingHeader->ProjSpectatorsn();
321 specPProj = hijingHeader->ProjSpectatorsp();
322 specNTarg = hijingHeader->TargSpectatorsn();
323 specPTarg = hijingHeader->TargSpectatorsp();
5d24ec95 324 AliInfo(Form("\t AliZDCDigitizer: b = %1.2f fm\n"
8fced614 325 " \t PROJECTILE: #spectator n %d, #spectator p %d\n"
326 " \t TARGET: #spectator n %d, #spectator p %d\n",
5d24ec95 327 impPar, specNProj, specPProj, specNTarg, specPTarg));
a18a0543 328
5d24ec95 329 // Applying fragmentation algorithm and adding spectator signal
330 Int_t freeSpecNProj=0, freeSpecPProj=0;
331 if(specNProj!=0 || specPProj!=0) Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
332 Int_t freeSpecNTarg=0, freeSpecPTarg=0;
333 if(specNTarg!=0 || specPTarg!=0) Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
334 if(freeSpecNProj!=0) SpectatorSignal(1, freeSpecNProj, pm);
335 if(freeSpecPProj!=0) SpectatorSignal(2, freeSpecPProj, pm);
336 AliInfo(Form("\t AliZDCDigitizer -> Adding spectator signal for PROJECTILE: %d free n and %d free p\n",freeSpecNProj,freeSpecPProj));
337 if(freeSpecNTarg!=0) SpectatorSignal(3, freeSpecNTarg, pm);
338 if(freeSpecPTarg!=0) SpectatorSignal(4, freeSpecPTarg, pm);
339 AliInfo(Form("\t AliZDCDigitizer -> Adding spectator signal for TARGET: %d free n and %d free p\n",freeSpecNTarg,freeSpecPTarg));
340 }
8309c1ab 341 }
342
8a2624cc 343
8309c1ab 344 // get the output run loader and loader
345 AliRunLoader* runLoader =
f21fc003 346 AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
8309c1ab 347 AliLoader* loader = runLoader->GetLoader("ZDCLoader");
8574b308 348 if(!loader) {
8309c1ab 349 AliError("no ZDC loader found");
350 return;
351 }
352
353 // create the output tree
354 const char* mode = "update";
8574b308 355 if(runLoader->GetEventNumber() == 0) mode = "recreate";
8309c1ab 356 loader->LoadDigits(mode);
357 loader->MakeTree("D");
358 TTree* treeD = loader->TreeD();
359 AliZDCDigit digit;
360 AliZDCDigit* pdigit = &digit;
361 const Int_t kBufferSize = 4000;
362 treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
363
364 // Create digits
80e87581 365 Int_t sector[2];
366 Int_t digi[2], digioot[2];
367 for(sector[0]=1; sector[0]<6; sector[0]++){
abf60186 368 for(sector[1]=0; sector[1]<5; sector[1]++){
369 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
8574b308 370 for(Int_t res=0; res<2; res++){
abf60186 371 digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
48642b09 372 + Pedestal(sector[0], sector[1], res);
7f73eb6b 373 }
abf60186 374 new(pdigit) AliZDCDigit(sector, digi);
8309c1ab 375 treeD->Fill();
698e394d 376
377 //Ch. debug
378 //printf("\t DIGIT added -> det %d quad %d - digi[0,1] = [%d, %d]\n",
379 // sector[0], sector[1], digi[0], digi[1]); // Chiara debugging!
380
8309c1ab 381 }
83347831 382 } // Loop over detector
383 // Adding in-time digits for 2 reference PTM signals (after signal ch.)
384 // (for the moment the ref. signal is completely invented assuming a PMgain of 5*10^4!)
385 Int_t sectorRef[2];
386 sectorRef[1] = 5;
9d38ced8 387 Int_t sigRef[2];
80e87581 388 // Reference signal are set to 100 (high gain chain) and 800 (low gain chain)
9d38ced8 389 if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
80e87581 390 else {sigRef[0]=0; sigRef[1]=0;} // calibration -> simulation of pedestal values
9d38ced8 391 //
83347831 392 for(Int_t iref=0; iref<2; iref++){
393 sectorRef[0] = 3*iref+1;
394 for(Int_t res=0; res<2; res++){
395 sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
396 }
83347831 397 new(pdigit) AliZDCDigit(sectorRef, sigRef);
398 treeD->Fill();
698e394d 399
400 //Ch. debug
401 //printf("\t RefDigit added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
402 // sectorRef[0], sectorRef[1], sigRef[0], sigRef[1]); // Chiara debugging!
abf60186 403 }
d79f8d50 404 //
405 // --- Adding digits for out-of-time channels after signal digits
80e87581 406 for(sector[0]=1; sector[0]<6; sector[0]++){
d79f8d50 407 for(sector[1]=0; sector[1]<5; sector[1]++){
408 if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
8574b308 409 for(Int_t res=0; res<2; res++){
d79f8d50 410 digioot[res] = Pedestal(sector[0], sector[1], res); // out-of-time ADCs
411 }
d79f8d50 412 new(pdigit) AliZDCDigit(sector, digioot);
413 treeD->Fill();
698e394d 414
415 //Ch. debug
416 //printf("\t DIGIToot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
417 // sector[0], sector[1], digioot[0], digioot[1]); // Chiara debugging!
d79f8d50 418 }
419 }
83347831 420 // Adding out-of-time digits for 2 reference PTM signals (after out-of-time ch.)
421 Int_t sigRefoot[2];
422 for(Int_t iref=0; iref<2; iref++){
423 sectorRef[0] = 3*iref+1;
424 for(Int_t res=0; res<2; res++){
425 sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
426 }
83347831 427 new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
428 treeD->Fill();
698e394d 429 //Ch. debug
430 //printf("\t RefDigitoot added -> det = %d, quad = %d - digi[0,1] = [%d, %d]\n",
431 // sectorRef[0], sectorRef[1], sigRefoot[0], sigRefoot[1]); // Chiara debugging!
432
83347831 433 }
434 //printf("\t AliZDCDigitizer -> TreeD has %d entries\n",(Int_t) treeD->GetEntries());
435
8309c1ab 436 // write the output tree
437 loader->WriteDigits("OVERWRITE");
438 loader->UnloadDigits();
439}
440
441
25175f38 442//_____________________________________________________________________________
443void AliZDCDigitizer::ReadPMTGains()
444{
445// Read PMT gain from an external file
446
447 char *fname = gSystem->ExpandPathName("$ALICE_ROOT/ZDC/PMTGainsdata.txt");
448 FILE *fdata = fopen(fname,"r");
449 if(fdata==NULL){
450 AliWarning(" Can't open file $ALICE_ROOT/ZDC/PMTGainsdata.txt to read ZDC PMT Gains\n");
451 AliWarning(" -> ZDC signal will be pedestal!!!!!!!!!!!!\n\n");
452 return;
453 }
3d7e0d0c 454 int read=1;
25175f38 455 Float_t data[5];
456 Int_t beam[12], det[12];
457 Float_t gain[12], aEne[12], bEne[12];
458 for(int ir=0; ir<12; ir++){
459 for(int ic=0; ic<5; ic++){
3d7e0d0c 460 read = fscanf(fdata,"%f ",&data[ic]);
461 if(read==0) AliDebug(3, " Error in reading PMT gains from external file ");
25175f38 462 }
463 beam[ir] = data[0];
464 det[ir] = data[1];
465 gain[ir] = data[2];
466 aEne[ir] = data[3];
467 bEne[ir] = data[4];
468 }
3d7ab74d 469 fclose(fdata);
25175f38 470
471 if(((fBeamType.CompareTo("P-P")) == 0)){
472 for(int i=0; i<12; i++){
473 if(beam[i]==0 && fBeamEnergy!=0.){
474 if(det[i]!=31 && det[i]!=32){
475 for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
476 }
477 else if(det[i] == 31) fPMGain[2][1] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
478 else if(det[i] == 32) fPMGain[2][2] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
479 }
480 }
481 //
0e9277fd 482 printf("\n AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for p-p @ %1.0f+%1.0f GeV: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
483 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
25175f38 484 }
485 else if(((fBeamType.CompareTo("A-A")) == 0)){
486 for(int i=0; i<12; i++){
487 if(beam[i]==1){
488 Float_t scalGainFactor = fBeamEnergy/2760.;
489 if(det[i]!=31 && det[i]!=32){
490 for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]/(aEne[i]*scalGainFactor);
491 }
492 else{
493 for(int iq=1; iq<3; iq++) fPMGain[2][iq] = gain[i]/(aEne[i]*scalGainFactor);
494 }
495 }
496 }
497 //
0e9277fd 498 printf("\n AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for Pb-Pb @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
499 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
25175f38 500 }
0e9277fd 501 else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A")) == 0)
502 || ((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P")) == 0)){
25175f38 503 for(int i=0; i<12; i++){
504 if(beam[i]==0 && fBeamEnergy!=0.){
505 if(det[i]==1 || det[i]==2){
506 for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
507 }
508 }
509 if(beam[i]==1){
510 Float_t scalGainFactor = fBeamEnergy/2760.;
511 Float_t npartScalingFactor = 208./15.;
512 if(det[i]==4 || det[i]==5){
513 for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
514 }
515 else if(det[i]==31 || det[i]==32){
516 for(int iq=1; iq<3; iq++) fPMGain[2][iq] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
517 }
518 }
519 }
0e9277fd 520 printf("\n AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for p-Pb: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
521 fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
25175f38 522 }
523}
524
5d24ec95 525//_____________________________________________________________________________
526void AliZDCDigitizer::CalculatePMTGains()
527{
528// Calculate PMT gain according to beam type and beam energy
58b2fa8a 529 if( ((fBeamType.CompareTo("P-P")) == 0) || ((fBeamType.CompareTo("p-p"))) ){
5d24ec95 530 // PTM gains rescaled to beam energy for p-p
531 // New correction coefficients for PMT gains needed
532 // to reproduce experimental spectra (from Grazia Jul 2010)
533 if(fBeamEnergy != 0){
534 for(Int_t j = 0; j < 5; j++){
535 fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
536 fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
537 fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;
538 fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
539 }
540 fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
541 fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
542 //
0e9277fd 543 printf("\n AliZDCDigitizer::CalculatePMTGains -> ZDC PMT gains for p-p @ %1.0f+%1.0f GeV: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
544 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
5d24ec95 545
546 }
547 }
548 else if(((fBeamType.CompareTo("A-A")) == 0)){
549 // PTM gains for Pb-Pb @ 2.7+2.7 A TeV ***************
550 // rescaled for Pb-Pb @ 1.38+1.38 A TeV ***************
551 // Values corrected after 2010 Pb-Pb data taking (7/2/2011 - Ch.)
552 // Experimental data compared to EMD simulation for single nucleon peaks:
553 // ZN gains must be divided by 4, ZP gains by 10!
554 Float_t scalGainFactor = fBeamEnergy/2760.;
555 for(Int_t j = 0; j < 5; j++){
556 fPMGain[0][j] = 50000./(4*scalGainFactor); // ZNC
557 fPMGain[1][j] = 100000./(5*scalGainFactor); // ZPC
558 fPMGain[2][j] = 100000./scalGainFactor; // ZEM
559 fPMGain[3][j] = 50000./(4*scalGainFactor); // ZNA
560 fPMGain[4][j] = 100000./(5*scalGainFactor); // ZPA
561 }
0e9277fd 562 printf("\n AliZDCDigitizer::CalculatePMTGains -> ZDC PMT gains for Pb-Pb @ %1.0f+%1.0f A GeV: ZN(%1.0f), ZP(%1.0f), ZEM(%1.0f)\n",
563 fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
5d24ec95 564 }
f6c61b40 565 else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A"))) ){
5d24ec95 566 // PTM gains for Pb-Pb @ 1.38+1.38 A TeV on side A
567 // PTM gains rescaled to beam energy for p-p on side C
568 // WARNING! Energies are set by hand for 2011 pA RUN!!!
312b8ed7 569 Float_t scalGainFactor = fBeamEnergy/2760.;
f6c61b40 570 Float_t npartScalingFactor = 208./15.;
5d24ec95 571
572 for(Int_t j = 0; j < 5; j++){
573 fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000; //ZNC (p)
574 fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000; //ZPC (p)
f6c61b40 575 fPMGain[2][j] = npartScalingFactor*100000./scalGainFactor; // ZEM (Pb)
312b8ed7 576 // Npart max scales from 400 in Pb-Pb to ~8 in pPb -> *40.
f6c61b40 577 fPMGain[3][j] = npartScalingFactor*50000/(4*scalGainFactor); // ZNA (Pb)
578 fPMGain[4][j] = npartScalingFactor*100000/(5*scalGainFactor); // ZPA (Pb)
5d24ec95 579 }
0e9277fd 580 printf("\n AliZDCDigitizer::CalculatePMTGains -> ZDC PMT gains for p-Pb: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
581 fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
5d24ec95 582 }
583}
584
8309c1ab 585//_____________________________________________________________________________
586void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
587 Int_t &freeSpecN, Int_t &freeSpecP) const
588{
589// simulate fragmentation of spectators
590
8309c1ab 591 AliZDCFragment frag(impPar);
8309c1ab 592
593 // Fragments generation
92339a90 594 frag.GenerateIMF();
595 Int_t nAlpha = frag.GetNalpha();
8309c1ab 596
597 // Attach neutrons
98711e01 598 frag.AttachNeutrons();
92339a90 599 Int_t ztot = frag.GetZtot();
600 Int_t ntot = frag.GetNtot();
98711e01 601
602 // Removing fragments and alpha pcs
8309c1ab 603 freeSpecN = specN-ntot-2*nAlpha;
604 freeSpecP = specP-ztot-2*nAlpha;
98711e01 605
3848c666 606 // Removing deuterons
a0607c1f 607 Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
3848c666 608 freeSpecN -= ndeu;
98711e01 609 freeSpecP -= ndeu;
3848c666 610 //
8309c1ab 611 if(freeSpecN<0) freeSpecN=0;
612 if(freeSpecP<0) freeSpecP=0;
613 AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
614}
615
616//_____________________________________________________________________________
617void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents,
fd9afd60 618 Float_t pm[5][5]) const
8309c1ab 619{
620// add signal of the spectators
fd9afd60 621
2d9c70ab 622 TFile *specSignalFile = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
623 if(!specSignalFile || !specSignalFile->IsOpen()) {
624 AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
625 return;
8309c1ab 626 }
2d9c70ab 627
e7eabc21 628 TNtuple* zdcSignal=0x0;
2d9c70ab 629
630 Float_t sqrtS = 2*fBeamEnergy;
631 //
632 if(TMath::Abs(sqrtS-5500) < 100.){
5d24ec95 633 AliInfo(" Extracting signal from SpectatorSignal/energy5500 ");
2d9c70ab 634 specSignalFile->cd("energy5500");
635 //
636 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
637 specSignalFile->GetObject("energy5500/ZNCSignal;1",zdcSignal);
638 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
639 }
640 else if(SpecType == 2) { // --- Signal for projectile spectator protons
641 specSignalFile->GetObject("energy5500/ZPCSignal;1",zdcSignal);
642 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
643 }
644 else if(SpecType == 3) { // --- Signal for target spectator neutrons
645 specSignalFile->GetObject("energy5500/ZNASignal;1",zdcSignal);
646 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
647 }
648 else if(SpecType == 4) { // --- Signal for target spectator protons
649 specSignalFile->GetObject("energy5500/ZPASignal;1",zdcSignal);
650 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
651 }
fd9afd60 652 }
2d9c70ab 653 else if(TMath::Abs(sqrtS-2760) < 100.){
5d24ec95 654 AliInfo(" Extracting signal from SpectatorSignal/energy2760 ");
2d9c70ab 655 specSignalFile->cd("energy2760");
656 //
657 if(SpecType == 1) { // --- Signal for projectile spectator neutrons
658 specSignalFile->GetObject("energy2760/ZNCSignal;1",zdcSignal);
659 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
660 }
661 else if(SpecType == 2) { // --- Signal for projectile spectator protons
662 specSignalFile->GetObject("energy2760/ZPCSignal;1",zdcSignal);
663 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
664 }
665 else if(SpecType == 3) { // --- Signal for target spectator neutrons
666 specSignalFile->GetObject("energy2760/ZNASignal;1",zdcSignal);
667 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
668 }
669 else if(SpecType == 4) { // --- Signal for target spectator protons
670 specSignalFile->GetObject("energy2760/ZPASignal;1",zdcSignal);
671 if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
672 }
fd9afd60 673 }
674
83534754 675 if(!zdcSignal){
676 printf("\n No spectator signal available for ZDC digitization\n");
677 return;
678 }
679
8309c1ab 680 Int_t nentries = (Int_t) zdcSignal->GetEntries();
681
fd9afd60 682 Float_t *entry;
683 Int_t pl, i, k, iev=0, rnd[125], volume[2];
48642b09 684 for(pl=0;pl<125;pl++) rnd[pl] = 0;
8574b308 685 if(numEvents > 125) {
fea9b334 686 AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
8309c1ab 687 numEvents = 125;
688 }
689 for(pl=0;pl<numEvents;pl++){
690 rnd[pl] = (Int_t) (9999*gRandom->Rndm());
85a96223 691 if(rnd[pl] >= 9999) rnd[pl] = 9998;
8309c1ab 692 //printf(" rnd[%d] = %d\n",pl,rnd[pl]);
693 }
694 // Sorting vector in ascending order with C function QSORT
695 qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
696 do{
697 for(i=0; i<nentries; i++){
698 zdcSignal->GetEvent(i);
699 entry = zdcSignal->GetArgs();
700 if(entry[0] == rnd[iev]){
701 for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
48642b09 702 //
fd9afd60 703 Float_t lightQ = entry[7];
704 Float_t lightC = entry[8];
705 //
acf9550c 706 if(volume[0] != 3) { // ZN or ZP
8309c1ab 707 pm[volume[0]-1][0] += lightC;
708 pm[volume[0]-1][volume[1]] += lightQ;
48642b09 709 //printf("\n pm[%d][0] = %.0f, pm[%d][%d] = %.0f\n",(volume[0]-1),pm[volume[0]-1][0],
710 // (volume[0]-1),volume[1],pm[volume[0]-1][volume[1]]);
711 }
712 else {
8574b308 713 if(volume[1] == 1) pm[2][1] += lightC; // ZEM 1
acf9550c 714 else pm[2][2] += lightQ; // ZEM 2
48642b09 715 //printf("\n pm[2][1] = %.0f, pm[2][2] = %.0f\n",pm[2][1],pm[2][2]);
716 }
8309c1ab 717 }
718 else if(entry[0] > rnd[iev]){
719 iev++;
720 continue;
721 }
722 }
723 }while(iev<numEvents);
724
2d9c70ab 725 specSignalFile->Close();
726 delete specSignalFile;
8309c1ab 727}
728
729
730//_____________________________________________________________________________
731Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
732 Int_t Res) const
733{
48642b09 734 // Evaluation of the ADC channel corresponding to the light yield Light
8a2624cc 735 Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
698e394d 736 // Ch. debug
737 //printf("\t Phe2ADCch -> det %d quad %d - PMGain[%d][%d] %1.0f phe %1.2f ADC %d\n",
738 // Det,Quad,Det-1,Quad,fPMGain[Det-1][Quad],Light,vADCch);
7f73eb6b 739
8a2624cc 740 return vADCch;
48642b09 741}
8309c1ab 742
48642b09 743//_____________________________________________________________________________
744Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
745{
8a2624cc 746 // Returns a pedestal for detector det, PM quad, channel with res.
747 //
65b16e87 748 Float_t pedValue;
abf60186 749 // Normal run
750 if(fIsCalibration == 0){
c35ed519 751 Int_t index=0, kNch=24;
83347831 752 if(Quad!=5){
c35ed519 753 if(Det==1) index = Quad+kNch*Res; // ZNC
754 else if(Det==2) index = (Quad+5)+kNch*Res; // ZPC
755 else if(Det==3) index = (Quad+9)+kNch*Res; // ZEM
756 else if(Det==4) index = (Quad+12)+kNch*Res; // ZNA
757 else if(Det==5) index = (Quad+17)+kNch*Res; // ZPA
83347831 758 }
c35ed519 759 else index = (Det-1)/3+22+kNch*Res; // Reference PMs
83347831 760 //
c35ed519 761 Float_t meanPed = fPedData->GetMeanPed(index);
762 Float_t pedWidth = fPedData->GetMeanPedWidth(index);
65b16e87 763 pedValue = gRandom->Gaus(meanPed,pedWidth);
abf60186 764 //
58e12fee 765 /*printf("\t AliZDCDigitizer::Pedestal -> det %d quad %d res %d - Ped[%d] = %d\n",
65b16e87 766 Det, Quad, Res, index,(Int_t) pedValue); // Chiara debugging!
abf60186 767 */
768 }
7f73eb6b 769 // To create calibration object
83347831 770 else{
65b16e87 771 if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm())); //High gain
772 else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm())); //Low gain
83347831 773 }
dbc8d7fc 774
65b16e87 775 return (Int_t) pedValue;
8309c1ab 776}
777
778//_____________________________________________________________________________
78d18275 779AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
8309c1ab 780{
8309c1ab 781
78d18275 782 Bool_t deleteManager = kFALSE;
48642b09 783
78d18275 784 AliCDBManager *manager = AliCDBManager::Instance();
785 AliCDBStorage *defstorage = manager->GetDefaultStorage();
48642b09 786
78d18275 787 if(!defstorage || !(defstorage->Contains("ZDC"))){
788 AliWarning("No default storage set or default storage doesn't contain ZDC!");
789 manager->SetDefaultStorage(uri);
790 deleteManager = kTRUE;
791 }
792
793 AliCDBStorage *storage = manager->GetDefaultStorage();
794
795 if(deleteManager){
796 AliCDBManager::Instance()->UnsetDefaultStorage();
797 defstorage = 0; // the storage is killed by AliCDBManager::Instance()->Destroy()
798 }
799
800 return storage;
801}
48642b09 802
78d18275 803//_____________________________________________________________________________
6024ec85 804AliZDCPedestals* AliZDCDigitizer::GetPedData() const
805{
806
807 // Getting pedestal calibration object for ZDC set
808
809 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
810 if(!entry) AliFatal("No calibration data loaded!");
811
812 AliZDCPedestals *calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
813 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
814
815 return calibdata;
816}
f91d1e35 817