class TFile;
#include <TF1.h>
#include <TFolder.h>
#include <TGeoGlobalMagField.h>
#include <TH1F.h>
#include <TRandom.h>
#include <TTree.h>
#include <TVirtualMC.h>
#include "AliMagF.h"
#include "AliPHOS.h"
#include "AliPHOSLoader.h"
#include "AliRun.h"
#include "AliRawReader.h"
#include "AliPHOSDigitizer.h"
#include "AliPHOSSDigitizer.h"
#include "AliPHOSDigit.h"
#include "AliAltroBuffer.h"
#include "AliAltroMapping.h"
#include "AliCaloAltroMapping.h"
#include "AliLog.h"
#include "AliCDBManager.h"
#include "AliCDBEntry.h"
#include "AliCDBStorage.h"
#include "AliPHOSCalibData.h"
#include "AliPHOSPulseGenerator.h"
#include "AliDAQ.h"
#include "AliPHOSRawFitterv0.h"
#include "AliPHOSCalibData.h"
#include "AliPHOSRawDigiProducer.h"
#include "AliPHOSQAChecker.h"
#include "AliPHOSRecoParam.h"
#include "AliPHOSSimParam.h"
ClassImp(AliPHOS)
AliPHOS:: AliPHOS() : AliDetector(),fgCalibData(0)
{
fName = "PHOS" ;
}
AliPHOS::AliPHOS(const char* name, const char* title): AliDetector(name, title),
fgCalibData(0)
{
}
AliPHOS::~AliPHOS()
{
if(fgCalibData) delete fgCalibData ;
}
AliDigitizer* AliPHOS::CreateDigitizer(AliDigitizationInput* digInput) const
{
return new AliPHOSDigitizer(digInput);
}
void AliPHOS::CreateMaterials()
{
Float_t aX[3] = {207.19, 183.85, 16.0} ;
Float_t zX[3] = {82.0, 74.0, 8.0} ;
Float_t wX[3] = {1.0, 1.0, 4.0} ;
Float_t dX = 8.28 ;
AliMixture(0, "PbWO4$", aX, zX, dX, -3, wX) ;
Float_t aP[2] = {12.011, 1.00794} ;
Float_t zP[2] = {6.0, 1.0} ;
Float_t wP[2] = {1.0, 1.0} ;
Float_t dP = 1.032 ;
AliMixture(1, "Polystyrene$", aP, zP, dP, -2, wP) ;
AliMaterial(2, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
Float_t aT[2] = {12.011, 1.00794} ;
Float_t zT[2] = {6.0, 1.0} ;
Float_t wT[2] = {1.0, 2.0} ;
Float_t dT = 0.331 ;
AliMixture(3, "Tyvek$", aT, zT, dT, -2, wT) ;
Float_t aF[2] = {12.011, 1.00794} ;
Float_t zF[2] = {6.0, 1.0} ;
Float_t wF[2] = {1.0, 1.0} ;
Float_t dF = 0.12 ;
AliMixture(4, "Foam$", aF, zF, dF, -2, wF) ;
Float_t aTIT[3] = {47.88, 26.98, 54.94} ;
Float_t zTIT[3] = {22.0, 13.0, 25.0} ;
Float_t wTIT[3] = {69.0, 6.0, 1.0} ;
Float_t dTIT = 4.5 ;
AliMixture(5, "Titanium$", aTIT, zTIT, dTIT, -3, wTIT);
AliMaterial(6, "Si$", 28.0855, 14., 2.33, 9.36, 42.3, 0, 0) ;
Float_t aTI[2] = {12.011, 1.00794} ;
Float_t zTI[2] = {6.0, 1.0} ;
Float_t wTI[2] = {1.0, 1.0} ;
Float_t dTI = 0.04 ;
AliMixture(7, "Thermo Insul.$", aTI, zTI, dTI, -2, wTI) ;
Float_t aTX[4] = {16.0, 28.09, 12.011, 1.00794} ;
Float_t zTX[4] = {8.0, 14.0, 6.0, 1.0} ;
Float_t wTX[4] = {292.0, 68.0, 462.0, 736.0} ;
Float_t dTX = 1.75 ;
AliMixture(8, "Textolit$", aTX, zTX, dTX, -4, wTX) ;
Float_t aFR[4] = {16.0, 28.09, 12.011, 1.00794} ;
Float_t zFR[4] = {8.0, 14.0, 6.0, 1.0} ;
Float_t wFR[4] = {292.0, 68.0, 462.0, 736.0} ;
Float_t dFR = 1.8 ;
AliMixture(9, "FR4$", aFR, zFR, dFR, -4, wFR) ;
Float_t aCM[2] = {12.01, 1.} ;
Float_t zCM[2] = {6., 1.} ;
Float_t wCM[2] = {1., 2.} ;
Float_t dCM = 0.935 ;
AliMixture(10, "Compo Mat$", aCM, zCM, dCM, -2, wCM) ;
AliMaterial(11, "Cu$", 63.546, 29, 8.96, 1.43, 14.8, 0, 0) ;
Float_t aG10[4] = { 12., 1., 16., 28.} ;
Float_t zG10[4] = { 6., 1., 8., 14.} ;
Float_t wG10[4] = { .259, .288, .248, .205} ;
Float_t dG10 = 1.7 ;
AliMixture(12, "G10$", aG10, zG10, dG10, -4, wG10);
AliMaterial(13, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
Float_t aCO[2] = {12.0, 16.0} ;
Float_t zCO[2] = {6.0, 8.0} ;
Float_t wCO[2] = {1.0, 2.0} ;
Float_t dCO = 0.001977 ;
AliMixture(14, "CO2$", aCO, zCO, dCO, -2, wCO);
Float_t dAr = 0.001782 ;
AliMaterial(15, "Ar$", 39.948, 18.0, dAr, 14.0, 0., 0, 0) ;
Float_t arContent = 0.80 ;
Float_t aArCO[3] = {39.948, 12.0, 16.0} ;
Float_t zArCO[3] = {18.0 , 6.0, 8.0} ;
Float_t wArCO[3];
wArCO[0] = arContent;
wArCO[1] = (1-arContent)*1;
wArCO[2] = (1-arContent)*2;
Float_t dArCO = arContent*dAr + (1-arContent)*dCO ;
AliMixture(16, "ArCO2$", aArCO, zArCO, dArCO, -3, wArCO) ;
AliMaterial(17, "Steel$", 55.845, 26, 7.87, 1.76, 0., 0, 0) ;
Float_t aFG[4] = {16.0, 28.09, 12.011, 1.00794} ;
Float_t zFG[4] = {8.0, 14.0, 6.0, 1.0} ;
Float_t wFG[4] = {292.0, 68.0, 462.0, 736.0} ;
Float_t dFG = 1.9 ;
AliMixture(18, "Fibergla$", aFG, zFG, dFG, -4, wFG) ;
Float_t aCA[4] = { 1.,12.,55.8,63.5 };
Float_t zCA[4] = { 1.,6.,26.,29. };
Float_t wCA[4] = { .014,.086,.42,.48 };
Float_t dCA = 0.8 ;
AliMixture(19, "Cables $", aCA, zCA, dCA, -4, wCA) ;
Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
Float_t zAir[4]={6.,7.,8.,18.};
Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
Float_t dAir = 1.20479E-3;
AliMixture(99, "Air$", aAir, zAir, dAir, 4, wAir) ;
Int_t isxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ() ;
Float_t sxmgmx = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max() ;
AliMedium(0, "PHOS Xtal $", 0, 1,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(1, "CPV scint. $", 1, 1,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(2, "Al parts $", 2, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
AliMedium(3, "Tyvek wrapper$", 3, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
AliMedium(4, "Polyst. foam $", 4, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(5, "Titan. cover $", 5, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.0001, 0.0001, 0, 0) ;
AliMedium(6, "Si PIN $", 6, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.01, 0.01, 0, 0) ;
AliMedium(7, "Thermo Insul.$", 7, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(8, "Textolit $", 8, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(9, "FR4 $", 9, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
AliMedium(10, "CompoMat $", 10, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(11, "Copper $", 11, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
AliMedium(12, "G10 $", 12, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.01, 0, 0) ;
AliMedium(13, "Lead $", 13, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(16, "ArCo2 $", 16, 1,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.01, 0, 0) ;
AliMedium(17, "Steel $", 17, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.0001, 0, 0) ;
AliMedium(18, "Fiberglass$", 18, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(19, "Cables $", 19, 0,
isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
AliMedium(99, "Air $", 99, 0,
isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
}
void AliPHOS::Init()
{
}
void AliPHOS::Digits2Raw()
{
if(AliPHOSSimParam::GetInstance()->IsEDigitizationOn()){
AliError("Energy digitization should be OFF if use Digits2Raw") ;
}
AliPHOSLoader * loader = static_cast<AliPHOSLoader*>(fLoader) ;
loader->LoadDigits();
TClonesArray* digits = loader->Digits() ;
if (!digits) {
AliError(Form("No digits found !"));
return;
}
AliPHOSGeometry* geom = GetGeometry();
if (!geom) {
AliError(Form("No geometry found !"));
return;
}
const TObjArray* maps = AliPHOSRecoParam::GetMappings();
if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
const Float_t kThreshold = 1.;
const Int_t kAdcThreshold = 1;
Int_t prevDDL = -1;
if(fgCalibData==0)
fgCalibData= new AliPHOSCalibData(-1) ;
AliPHOSPulseGenerator *pulse = new AliPHOSPulseGenerator();
pulse->SetTimeStep(TMath::Abs(fgCalibData->GetSampleTimeStep())) ;
Int_t *adcValuesLow = new Int_t[pulse->GetRawFormatTimeBins()];
Int_t *adcValuesHigh= new Int_t[pulse->GetRawFormatTimeBins()];
const Int_t maxDDL = 20;
AliAltroBuffer *buffer[maxDDL];
AliAltroMapping *mapping[maxDDL];
for(Int_t jDDL=0; jDDL<maxDDL; jDDL++) {
buffer[jDDL]=0;
mapping[jDDL]=0;
}
Int_t modMax=-111;
Int_t colMax=-111;
Int_t rowMax=-111;
Float_t eMax=-333;
for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
AliPHOSDigit* digit = static_cast<AliPHOSDigit *>(digits->At(iDigit)) ;
if (digit->GetEnergy() < kThreshold)
continue;
Int_t relId[4];
geom->AbsToRelNumbering(digit->GetId(), relId);
Int_t module = 5-relId[0];
if (relId[1] != 0)
continue;
Int_t row = relId[2]-1;
Int_t col = relId[3]-1;
Int_t iRCU = -111;
if(0<=row&&row<16 && 0<=col&&col<56) iRCU=0;
if(16<=row&&row<32 && 0<=col&&col<56) iRCU=1;
if(32<=row&&row<48 && 0<=col&&col<56) iRCU=2;
if(48<=row&&row<64 && 0<=col&&col<56) iRCU=3;
Int_t iDDL = 4 * module + iRCU;
if (iDDL != prevDDL) {
if (buffer[iDDL] == 0) {
TString fileName = AliDAQ::DdlFileName("PHOS",iDDL);
mapping[iDDL] = (AliAltroMapping*)maps->At(iDDL);
buffer[iDDL] = new AliAltroBuffer(fileName.Data(),mapping[iDDL]);
buffer[iDDL]->WriteDataHeader(kTRUE, kFALSE);
}
prevDDL = iDDL;
}
AliDebug(2,Form("digit E=%.4f GeV, t=%g s, (mod,col,row)=(%d,%d,%d)\n",
digit->GetEnergy(),digit->GetTimeR(),
relId[0]-1,relId[3]-1,relId[2]-1));
if (digit->GetTimeR() > pulse->GetRawFormatTimeMax()*0.5 ) {
AliDebug(2,"Signal is out of time range.\n");
buffer[iDDL]->FillBuffer(0);
buffer[iDDL]->FillBuffer(pulse->GetRawFormatTimeBins() );
buffer[iDDL]->FillBuffer(3);
buffer[iDDL]->WriteTrailer(3, relId[3]-1, relId[2]-1, 0);
} else {
Double_t energy = 0 ;
if (digit->GetId() <= geom->GetNModules() * geom->GetNCristalsInModule()) {
energy=digit->GetEnergy();
if(energy>eMax) {eMax=energy; modMax=relId[0]; colMax=col; rowMax=row;}
}
else {
energy = 0;
}
pulse->SetAmplitude(energy);
pulse->SetTZero(digit->GetTimeR());
Double_t r =fgCalibData->GetHighLowRatioEmc(relId[0],relId[3],relId[2]) ;
pulse->SetHG2LGRatio(r) ;
pulse->MakeSamples();
pulse->GetSamples(adcValuesHigh, adcValuesLow) ;
buffer[iDDL]->WriteChannel(relId[3]-1, relId[2]-1, 0,
pulse->GetRawFormatTimeBins(), adcValuesLow , kAdcThreshold);
buffer[iDDL]->WriteChannel(relId[3]-1, relId[2]-1, 1,
pulse->GetRawFormatTimeBins(), adcValuesHigh, kAdcThreshold);
}
}
delete [] adcValuesLow;
delete [] adcValuesHigh;
for (Int_t iDDL=0; iDDL<maxDDL; iDDL++) {
if (buffer[iDDL]) {
buffer[iDDL]->Flush();
buffer[iDDL]->WriteDataHeader(kFALSE, kFALSE);
delete buffer[iDDL];
}
}
AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n",
modMax,colMax,rowMax,eMax));
delete pulse;
loader->UnloadDigits();
}
void AliPHOS::Hits2SDigits()
{
AliPHOSSDigitizer phosDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
phosDigitizer.SetEventRange(0, -1) ;
phosDigitizer.Digitize("all") ;
}
AliLoader* AliPHOS::MakeLoader(const char* topfoldername)
{
fLoader = new AliPHOSLoader(GetName(),topfoldername);
return fLoader;
}
void AliPHOS::SetTreeAddress()
{
TBranch *branch;
char branchname[20];
snprintf(branchname,20,"%s",GetName());
TTree *treeH = fLoader->TreeH();
if (treeH) {
branch = treeH->GetBranch(branchname);
if (branch)
{
if (fHits == 0x0) fHits= new TClonesArray("AliPHOSHit",1000);
branch->SetAddress(&fHits);
}
}
}
Bool_t AliPHOS::Raw2SDigits(AliRawReader* rawReader)
{
AliPHOSLoader * loader = static_cast<AliPHOSLoader*>(fLoader) ;
TTree * tree = 0 ;
tree = loader->TreeS() ;
if ( !tree ) {
loader->MakeTree("S");
tree = loader->TreeS() ;
}
TClonesArray * sdigits = loader->SDigits() ;
if(!sdigits) {
loader->MakeSDigitsArray();
sdigits = loader->SDigits();
}
sdigits->Clear();
rawReader->Reset() ;
const TObjArray* maps = AliPHOSRecoParam::GetMappings();
if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
AliAltroMapping *mapping[20];
for(Int_t i = 0; i < 20; i++) {
mapping[i] = (AliAltroMapping*)maps->At(i);
}
AliPHOSRawFitterv0 fitter;
fitter.SubtractPedestals(AliPHOSSimParam::GetInstance()->EMCSubtractPedestals());
fitter.SetAmpOffset(AliPHOSSimParam::GetInstance()->GetGlobalAltroOffset());
fitter.SetAmpThreshold(AliPHOSSimParam::GetInstance()->GetGlobalAltroThreshold());
AliPHOSRawDigiProducer pr(rawReader,mapping);
pr.SetSampleQualityCut(AliPHOSSimParam::GetInstance()->GetEMCSampleQualityCut());
pr.MakeDigits(sdigits,&fitter);
Int_t bufferSize = 32000 ;
tree->Branch("PHOS",&sdigits,bufferSize);
tree->Fill();
fLoader->WriteSDigits("OVERWRITE");
return kTRUE;
}