#include <stdlib.h>
#include <stdio.h>
#include <TTree.h>
#include <TFile.h>
#include <TNtuple.h>
#include <TRandom.h>
#include "AliRun.h"
#include "AliHeader.h"
#include "AliGenHijingEventHeader.h"
#include "AliGenCocktailEventHeader.h"
#include "AliDigitizationInput.h"
#include "AliRunLoader.h"
#include "AliLoader.h"
#include "AliGRPObject.h"
#include "AliCDBManager.h"
#include "AliCDBEntry.h"
#include "AliZDCSDigit.h"
#include "AliZDCDigit.h"
#include "AliZDCFragment.h"
#include "AliZDCv3.h"
#include "AliZDCDigitizer.h"
class AliCDBStorage;
class AliZDCPedestals;
ClassImp(AliZDCDigitizer)
AliZDCDigitizer::AliZDCDigitizer() :
fIsCalibration(0),
fIsSignalInADCGate(kFALSE),
fFracLostSignal(0.),
fPedData(0),
fSpectators2Track(kFALSE),
fBeamEnergy(0.),
fBeamType(""),
fIspASystem(kFALSE),
fIsRELDISgen(kFALSE),
fSpectatorData(0)
{
for(Int_t i=0; i<2; i++) fADCRes[i]=0.;
}
AliZDCDigitizer::AliZDCDigitizer(AliDigitizationInput* digInput):
AliDigitizer(digInput),
fIsCalibration(0),
fIsSignalInADCGate(kFALSE),
fFracLostSignal(0.),
fPedData(GetPedData()),
fSpectators2Track(kFALSE),
fBeamEnergy(0.),
fBeamType(""),
fIspASystem(kFALSE),
fIsRELDISgen(kFALSE),
fSpectatorData(NULL)
{
if(fIsCalibration!=0) printf("\n\t AliZDCDigitizer -> Creating calibration data (pedestals)\n");
for(Int_t i=0; i<5; i++){
for(Int_t j=0; j<5; j++)
fPMGain[i][j] = 0.;
}
}
AliZDCDigitizer::~AliZDCDigitizer()
{
if(fSpectatorData!=NULL){
fSpectatorData->Close();
delete fSpectatorData;
}
}
AliZDCDigitizer::AliZDCDigitizer(const AliZDCDigitizer &digitizer):
AliDigitizer(),
fIsCalibration(digitizer.fIsCalibration),
fIsSignalInADCGate(digitizer.fIsSignalInADCGate),
fFracLostSignal(digitizer.fFracLostSignal),
fPedData(digitizer.fPedData),
fSpectators2Track(digitizer.fSpectators2Track),
fBeamEnergy(digitizer.fBeamEnergy),
fBeamType(digitizer.fBeamType),
fIspASystem(digitizer.fIspASystem),
fIsRELDISgen(digitizer.fIsRELDISgen),
fSpectatorData(digitizer.fSpectatorData)
{
for(Int_t i=0; i<5; i++){
for(Int_t j=0; j<5; j++){
fPMGain[i][j] = digitizer.fPMGain[i][j];
}
}
for(Int_t i=0; i<2; i++) fADCRes[i] = digitizer.fADCRes[i];
}
Bool_t AliZDCDigitizer::Init()
{
AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
if(!entry) AliFatal("No calibration data loaded!");
AliGRPObject* grpData = 0x0;
if(entry){
TMap* m = dynamic_cast<TMap*>(entry->GetObject());
if(m){
grpData = new AliGRPObject();
grpData->ReadValuesFromMap(m);
}
else{
grpData = dynamic_cast<AliGRPObject*>(entry->GetObject());
}
entry->SetOwner(0);
AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
}
if(!grpData){
AliError("No GRP entry found in OCDB! \n ");
return kFALSE;
}
fBeamType = grpData->GetBeamType();
if(fBeamType==AliGRPObject::GetInvalidString()){
AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
}
if(fBeamEnergy<0.01){
fBeamEnergy = grpData->GetBeamEnergy();
if(fBeamEnergy==AliGRPObject::GetInvalidFloat()){
AliWarning("GRP/GRP/Data entry: missing value for the beam energy ! Using 0.");
AliError("\t UNKNOWN beam type from GRP obj -> PMT gains not set in ZDC digitizer!!!\n");
fBeamEnergy = 0.;
}
}
if(((fBeamType.CompareTo("UNKNOWN")) == 0) || fIsRELDISgen){
fBeamType = "A-A";
AliInfo(" AliZDCDigitizer -> Manually setting beam type to A-A\n");
}
printf("\n\t AliZDCDigitizer -> beam type %s - beam energy = %f GeV\n", fBeamType.Data(), fBeamEnergy);
if(fSpectators2Track) printf("\n\t AliZDCDigitizer -> spectator signal added at digit level\n\n");
if(fBeamEnergy>0.1){
ReadPMTGains();
}
else{
AliWarning("\n Beam energy is 0 -> ZDC PMT gains can't be set -> NO ZDC DIGITS!!!\n");
}
fADCRes[0] = 0.0000008;
fADCRes[1] = 0.0000064;
return kTRUE;
}
void AliZDCDigitizer::Digitize(Option_t* )
{
Float_t pm[5][5];
for(Int_t iSector1=0; iSector1<5; iSector1++)
for(Int_t iSector2=0; iSector2<5; iSector2++){
pm[iSector1][iSector2] = 0;
}
Float_t impPar = 0;
Int_t specNTarg = 0, specPTarg = 0;
Int_t specNProj = 0, specPProj = 0;
Float_t signalTime0 = 0.;
for(Int_t iInput = 0; iInput<fDigInput->GetNinputs(); iInput++){
AliRunLoader* runLoader =
AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
AliLoader* loader = runLoader->GetLoader("ZDCLoader");
if(!loader) continue;
loader->LoadSDigits();
TTree* treeS = loader->TreeS();
if(!treeS) continue;
AliZDCSDigit sdigit;
AliZDCSDigit* psdigit = &sdigit;
treeS->SetBranchAddress("ZDC", &psdigit);
for(Int_t iSDigit=0; iSDigit<treeS->GetEntries(); iSDigit++){
treeS->GetEntry(iSDigit);
if(!psdigit) continue;
if((sdigit.GetSector(1) < 0) || (sdigit.GetSector(1) > 4)){
AliError(Form("\nsector[0] = %d, sector[1] = %d\n",
sdigit.GetSector(0), sdigit.GetSector(1)));
continue;
}
if(iSDigit==0) signalTime0 = sdigit.GetTrackTime();
else{
if(sdigit.GetTrackTime()<=signalTime0+30.) fIsSignalInADCGate = kTRUE;
if(sdigit.GetTrackTime()>signalTime0+30.){
fIsSignalInADCGate = kFALSE;
fFracLostSignal = (sdigit.GetTrackTime()-30)*(sdigit.GetTrackTime()-30)/280.;
}
}
Float_t sdSignal = sdigit.GetLightPM();
if(fIsSignalInADCGate == kFALSE){
AliDebug(2,Form("\t Signal time %f -> fraction %f of ZDC signal for det.(%d, %d) out of ADC gate\n",
sdigit.GetTrackTime(),fFracLostSignal,sdigit.GetSector(0),sdigit.GetSector(1)));
sdSignal = (1-fFracLostSignal)*sdSignal;
}
pm[(sdigit.GetSector(0))-1][sdigit.GetSector(1)] += sdigit.GetLightPM();
}
loader->UnloadSDigits();
if(!runLoader->GetAliRun()) runLoader->LoadgAlice();
AliHeader* header = runLoader->GetHeader();
if(!header) continue;
AliGenEventHeader* genHeader = header->GenEventHeader();
if(!genHeader) continue;
AliGenHijingEventHeader *hijingHeader = 0;
if(genHeader->InheritsFrom(AliGenHijingEventHeader::Class())) hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (genHeader);
else if(genHeader->InheritsFrom(AliGenCocktailEventHeader::Class())){
TList* listOfHeaders = ((AliGenCocktailEventHeader*) genHeader)->GetHeaders();
if(listOfHeaders){
for(Int_t iH = 0; iH < listOfHeaders->GetEntries(); ++iH) {
AliGenEventHeader *currHeader = dynamic_cast <AliGenEventHeader *> (listOfHeaders->At(iH));
if (currHeader && currHeader->InheritsFrom(AliGenHijingEventHeader::Class())) {
hijingHeader = dynamic_cast <AliGenHijingEventHeader*> (currHeader);
break;
}
}
}
else{
printf(" No list of headers from generator \n");
}
}
if(!hijingHeader){
printf(" No HIJING header found in list of headers from generator\n");
}
if(hijingHeader && fSpectators2Track==kTRUE){
impPar = hijingHeader->ImpactParameter();
specNProj = hijingHeader->ProjSpectatorsn();
specPProj = hijingHeader->ProjSpectatorsp();
specNTarg = hijingHeader->TargSpectatorsn();
specPTarg = hijingHeader->TargSpectatorsp();
Int_t freeSpecNProj=0, freeSpecPProj=0;
if(specNProj!=0 || specPProj!=0) Fragmentation(impPar, specNProj, specPProj, freeSpecNProj, freeSpecPProj);
Int_t freeSpecNTarg=0, freeSpecPTarg=0;
if(specNTarg!=0 || specPTarg!=0) Fragmentation(impPar, specNTarg, specPTarg, freeSpecNTarg, freeSpecPTarg);
if(freeSpecNProj!=0) SpectatorSignal(1, freeSpecNProj, pm);
if(freeSpecPProj!=0) SpectatorSignal(2, freeSpecPProj, pm);
if(freeSpecNTarg!=0) SpectatorSignal(3, freeSpecNTarg, pm);
if(freeSpecPTarg!=0) SpectatorSignal(4, freeSpecPTarg, pm);
}
}
AliRunLoader* runLoader =
AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
AliLoader* loader = runLoader->GetLoader("ZDCLoader");
if(!loader) {
AliError("no ZDC loader found");
return;
}
const char* mode = "update";
if(runLoader->GetEventNumber() == 0) mode = "recreate";
loader->LoadDigits(mode);
loader->MakeTree("D");
TTree* treeD = loader->TreeD();
AliZDCDigit digit;
AliZDCDigit* pdigit = &digit;
const Int_t kBufferSize = 4000;
treeD->Branch("ZDC", "AliZDCDigit", &pdigit, kBufferSize);
Int_t sector[2];
Int_t digi[2], digioot[2];
for(sector[0]=1; sector[0]<6; sector[0]++){
for(sector[1]=0; sector[1]<5; sector[1]++){
if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
for(Int_t res=0; res<2; res++){
digi[res] = Phe2ADCch(sector[0], sector[1], pm[sector[0]-1][sector[1]], res)
+ Pedestal(sector[0], sector[1], res);
}
new(pdigit) AliZDCDigit(sector, digi);
treeD->Fill();
}
}
Int_t sectorRef[2];
sectorRef[1] = 5;
Int_t sigRef[2];
if(fIsCalibration==0) {sigRef[0]=100; sigRef[1]=800;}
else {sigRef[0]=0; sigRef[1]=0;}
for(Int_t iref=0; iref<2; iref++){
sectorRef[0] = 3*iref+1;
for(Int_t res=0; res<2; res++){
sigRef[res] += Pedestal(sectorRef[0], sectorRef[1], res);
}
new(pdigit) AliZDCDigit(sectorRef, sigRef);
treeD->Fill();
}
for(sector[0]=1; sector[0]<6; sector[0]++){
for(sector[1]=0; sector[1]<5; sector[1]++){
if((sector[0]==3) && ((sector[1]<1) || (sector[1]>2))) continue;
for(Int_t res=0; res<2; res++){
digioot[res] = Pedestal(sector[0], sector[1], res);
}
new(pdigit) AliZDCDigit(sector, digioot);
treeD->Fill();
}
}
Int_t sigRefoot[2];
for(Int_t iref=0; iref<2; iref++){
sectorRef[0] = 3*iref+1;
for(Int_t res=0; res<2; res++){
sigRefoot[res] = Pedestal(sectorRef[0], sectorRef[1], res);
}
new(pdigit) AliZDCDigit(sectorRef, sigRefoot);
treeD->Fill();
}
loader->WriteDigits("OVERWRITE");
loader->UnloadDigits();
}
void AliZDCDigitizer::ReadPMTGains()
{
char *fname = gSystem->ExpandPathName("$ALICE_ROOT/ZDC/PMTGainsdata.txt");
FILE *fdata = fopen(fname,"r");
if(fdata==NULL){
AliWarning(" Can't open file $ALICE_ROOT/ZDC/PMTGainsdata.txt to read ZDC PMT Gains\n");
AliWarning(" -> ZDC signal will be pedestal!!!!!!!!!!!!\n\n");
return;
}
int read=1;
Float_t data[5];
Int_t beam[12], det[12];
Float_t gain[12], aEne[12], bEne[12];
for(int ir=0; ir<12; ir++){
for(int ic=0; ic<5; ic++){
read = fscanf(fdata,"%f ",&data[ic]);
if(read==0) AliDebug(3, " Error in reading PMT gains from external file ");
}
beam[ir] = (int) (data[0]);
det[ir] = (int) (data[1]);
gain[ir] = data[2];
aEne[ir] = data[3];
bEne[ir] = data[4];
}
fclose(fdata);
if(((fBeamType.CompareTo("P-P")) == 0) || ((fBeamType.CompareTo("p-p")) == 0)){
for(int i=0; i<12; i++){
if(beam[i]==0 && fBeamEnergy!=0.){
if(det[i]!=31 && det[i]!=32){
for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
}
else if(det[i] == 31) fPMGain[2][1] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
else if(det[i] == 32) fPMGain[2][2] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
}
}
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",
fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
}
else if(((fBeamType.CompareTo("A-A")) == 0)){
for(int i=0; i<12; i++){
if(beam[i]==1){
Float_t scalGainFactor = fBeamEnergy/2760.;
if(det[i]!=31 && det[i]!=32){
for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]/(aEne[i]*scalGainFactor);
}
else{
for(int iq=1; iq<3; iq++) fPMGain[2][iq] = gain[i]/(aEne[i]*scalGainFactor);
}
}
}
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",
fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
}
else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A")) == 0)){
for(int i=0; i<12; i++){
if(beam[i]==0 && fBeamEnergy!=0.){
if(det[i]==1 || det[i]==2){
for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
}
}
if(beam[i]==1){
Float_t scalGainFactor = fBeamEnergy/2760.;
Float_t npartScalingFactor = 208./15.;
if(det[i]==4 || det[i]==5){
for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
}
else if(det[i]==31 || det[i]==32){
for(int iq=1; iq<3; iq++) fPMGain[2][iq] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
}
}
}
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",
fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
}
else if(((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P")) == 0)){
for(int i=0; i<12; i++){
if(beam[i]==0 && fBeamEnergy!=0.){
if(det[i]==4 || det[i]==5){
for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = gain[i]*(aEne[i]/fBeamEnergy+bEne[i]);
}
else if(det[i] == 31) fPMGain[2][1] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
else if(det[i] == 32) fPMGain[2][2] = gain[i]*(aEne[i]-fBeamEnergy*bEne[i]);
}
if(beam[i]==1){
Float_t scalGainFactor = fBeamEnergy/2760.;
Float_t npartScalingFactor = 208./15.;
if(det[i]==1 || det[i]==2){
for(Int_t j=0; j<5; j++) fPMGain[det[i]-1][j] = npartScalingFactor*gain[i]/(aEne[i]*scalGainFactor);
}
}
}
printf("\n AliZDCDigitizer::ReadPMTGains -> ZDC PMT gains for Pb-p: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
}
}
void AliZDCDigitizer::CalculatePMTGains()
{
if( ((fBeamType.CompareTo("P-P")) == 0) || ((fBeamType.CompareTo("p-p"))) ){
if(fBeamEnergy != 0){
for(Int_t j = 0; j < 5; j++){
fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;
fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
}
fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
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",
fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
}
}
else if(((fBeamType.CompareTo("A-A")) == 0)){
Float_t scalGainFactor = fBeamEnergy/2760.;
for(Int_t j = 0; j < 5; j++){
fPMGain[0][j] = 50000./(4*scalGainFactor);
fPMGain[1][j] = 100000./(5*scalGainFactor);
fPMGain[2][j] = 100000./scalGainFactor;
fPMGain[3][j] = 50000./(4*scalGainFactor);
fPMGain[4][j] = 100000./(5*scalGainFactor);
}
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",
fBeamEnergy, fBeamEnergy, fPMGain[0][0], fPMGain[1][0], fPMGain[2][1]);
}
else if(((fBeamType.CompareTo("p-A")) == 0) || ((fBeamType.CompareTo("P-A"))) ){
Float_t scalGainFactor = fBeamEnergy/2760.;
Float_t npartScalingFactor = 208./15.;
for(Int_t j = 0; j < 5; j++){
fPMGain[0][j] = 1.515831*(661.444/fBeamEnergy+0.000740671)*10000000;
fPMGain[1][j] = 0.674234*(864.350/fBeamEnergy+0.00234375)*10000000;
if(j<2) fPMGain[2][j] = npartScalingFactor*100000./scalGainFactor;
fPMGain[3][j] = npartScalingFactor*50000/(4*scalGainFactor);
fPMGain[4][j] = npartScalingFactor*100000/(5*scalGainFactor);
}
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",
fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
}
else if(((fBeamType.CompareTo("A-p")) == 0) || ((fBeamType.CompareTo("A-P"))) ){
Float_t scalGainFactor = fBeamEnergy/2760.;
Float_t npartScalingFactor = 208./15.;
for(Int_t j = 0; j < 5; j++){
fPMGain[3][j] = 1.350938*(661.444/fBeamEnergy+0.000740671)*10000000;
fPMGain[4][j] = 0.678597*(864.350/fBeamEnergy+0.00234375)*10000000;
fPMGain[1][j] = npartScalingFactor*50000/(4*scalGainFactor);
fPMGain[2][j] = npartScalingFactor*100000/(5*scalGainFactor);
}
fPMGain[2][1] = 0.869654*(1.32312-0.000101515*fBeamEnergy)*10000000;
fPMGain[2][2] = 1.030883*(1.32312-0.000101515*fBeamEnergy)*10000000;
printf("\n AliZDCDigitizer::CalculatePMTGains -> ZDC PMT gains for Pb-p: ZNC(%1.0f), ZPC(%1.0f), ZEM(%1.0f), ZNA(%1.0f) ZPA(%1.0f)\n",
fPMGain[0][0], fPMGain[1][0], fPMGain[2][1], fPMGain[3][0], fPMGain[4][0]);
}
}
void AliZDCDigitizer::Fragmentation(Float_t impPar, Int_t specN, Int_t specP,
Int_t &freeSpecN, Int_t &freeSpecP) const
{
AliZDCFragment frag(impPar);
frag.GenerateIMF();
Int_t nAlpha = frag.GetNalpha();
frag.AttachNeutrons();
Int_t ztot = frag.GetZtot();
Int_t ntot = frag.GetNtot();
freeSpecN = specN-ntot-2*nAlpha;
freeSpecP = specP-ztot-2*nAlpha;
Int_t ndeu = (Int_t) (freeSpecN*frag.DeuteronNumber());
freeSpecN -= ndeu;
freeSpecP -= ndeu;
if(freeSpecN<0) freeSpecN=0;
if(freeSpecP<0) freeSpecP=0;
AliDebug(2, Form("FreeSpn = %d, FreeSpp = %d", freeSpecN, freeSpecP));
}
void AliZDCDigitizer::SpectatorSignal(Int_t SpecType, Int_t numEvents, Float_t pm[5][5])
{
if(!fSpectatorData) fSpectatorData = TFile::Open("$ALICE_ROOT/ZDC/SpectatorSignal.root");
if(!fSpectatorData || !fSpectatorData->IsOpen()) {
AliError((" Opening file $ALICE_ROOT/ZDC/SpectatorSignal.root failed\n"));
return;
}
TNtuple* zdcSignal=0x0;
Float_t sqrtS = 2*fBeamEnergy;
if(TMath::Abs(sqrtS-5500) < 100.){
AliInfo(" Extracting signal from SpectatorSignal/energy5500 ");
fSpectatorData->cd("energy5500");
if(SpecType == 1) {
fSpectatorData->GetObject("energy5500/ZNCSignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
}
else if(SpecType == 2) {
fSpectatorData->GetObject("energy5500/ZPCSignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
}
else if(SpecType == 3) {
fSpectatorData->GetObject("energy5500/ZNASignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
}
else if(SpecType == 4) {
fSpectatorData->GetObject("energy5500/ZPASignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
}
}
else if(TMath::Abs(sqrtS-2760) < 100.){
AliInfo(" Extracting signal from SpectatorSignal/energy2760 ");
fSpectatorData->cd("energy2760");
if(SpecType == 1) {
fSpectatorData->GetObject("energy2760/ZNCSignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNCSignal from SpectatorSignal.root file");
}
else if(SpecType == 2) {
fSpectatorData->GetObject("energy2760/ZPCSignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPCSignal from SpectatorSignal.root file");
}
else if(SpecType == 3) {
fSpectatorData->GetObject("energy2760/ZNASignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZNASignal from SpectatorSignal.root file");
}
else if(SpecType == 4) {
fSpectatorData->GetObject("energy2760/ZPASignal;1",zdcSignal);
if(!zdcSignal) AliError(" PROBLEM!!! Can't retrieve ZPASignal from SpectatorSignal.root file");
}
}
if(!zdcSignal){
printf("\n No spectator signal available for ZDC digitization\n");
return;
}
Int_t nentries = (Int_t) zdcSignal->GetEntries();
Float_t *entry;
Int_t pl, i, k, iev=0, rnd[125], volume[2];
for(pl=0;pl<125;pl++) rnd[pl] = 0;
if(numEvents > 125) {
AliDebug(2,Form("numEvents (%d) is larger than 125", numEvents));
numEvents = 125;
}
for(pl=0;pl<numEvents;pl++){
rnd[pl] = (Int_t) (9999*gRandom->Rndm());
if(rnd[pl] >= 9999) rnd[pl] = 9998;
}
qsort((void*)rnd,numEvents,sizeof(Int_t),comp);
do{
for(i=0; i<nentries; i++){
zdcSignal->GetEvent(i);
entry = zdcSignal->GetArgs();
if(entry[0] == rnd[iev]){
for(k=0; k<2; k++) volume[k] = (Int_t) entry[k+1];
Float_t lightQ = entry[7];
Float_t lightC = entry[8];
if(volume[0] != 3) {
pm[volume[0]-1][0] += lightC;
pm[volume[0]-1][volume[1]] += lightQ;
}
else {
if(volume[1] == 1) pm[2][1] += lightC;
else pm[2][2] += lightQ;
}
}
else if(entry[0] > rnd[iev]){
iev++;
continue;
}
}
}while(iev<numEvents);
}
Int_t AliZDCDigitizer::Phe2ADCch(Int_t Det, Int_t Quad, Float_t Light,
Int_t Res) const
{
Int_t vADCch = (Int_t) (Light * fPMGain[Det-1][Quad] * fADCRes[Res]);
return vADCch;
}
Int_t AliZDCDigitizer::Pedestal(Int_t Det, Int_t Quad, Int_t Res) const
{
Float_t pedValue=0.;
if(fIsCalibration == 0){
Int_t index=0, kNch=24;
if(Quad!=5){
if(Det==1) index = Quad+kNch*Res;
else if(Det==2) index = (Quad+5)+kNch*Res;
else if(Det==3) index = (Quad+9)+kNch*Res;
else if(Det==4) index = (Quad+12)+kNch*Res;
else if(Det==5) index = (Quad+17)+kNch*Res;
}
else index = (Det-1)/3+22+kNch*Res;
if(fPedData){
Float_t meanPed = fPedData->GetMeanPed(index);
Float_t pedWidth = fPedData->GetMeanPedWidth(index);
pedValue = gRandom->Gaus(meanPed,pedWidth);
}
else{
printf(" AliZDCDigitizer::Pedestal -> No valid pedestal calibration object loaded!\n\n");
}
}
else{
if(Res == 0) pedValue = gRandom->Gaus((35.+10.*gRandom->Rndm()),(0.5+0.2*gRandom->Rndm()));
else pedValue = gRandom->Gaus((250.+100.*gRandom->Rndm()),(3.5+2.*gRandom->Rndm()));
}
return (Int_t) pedValue;
}
AliCDBStorage* AliZDCDigitizer::SetStorage(const char *uri)
{
Bool_t deleteManager = kFALSE;
AliCDBManager *manager = AliCDBManager::Instance();
AliCDBStorage *defstorage = manager->GetDefaultStorage();
if(!defstorage || !(defstorage->Contains("ZDC"))){
AliWarning("No default storage set or default storage doesn't contain ZDC!");
manager->SetDefaultStorage(uri);
deleteManager = kTRUE;
}
AliCDBStorage *storage = manager->GetDefaultStorage();
if(deleteManager){
AliCDBManager::Instance()->UnsetDefaultStorage();
defstorage = 0;
}
return storage;
}
AliZDCPedestals* AliZDCDigitizer::GetPedData() const
{
AliZDCPedestals *calibdata = 0x0;
AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
if(!entry) AliFatal("No calibration data loaded!");
else{
calibdata = dynamic_cast<AliZDCPedestals*> (entry->GetObject());
if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
}
return calibdata;
}