+
/**************************************************************************
* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
* *
* provided "as is" without express or implied warranty. *
**************************************************************************/
-/*
-$Log$
-Revision 1.11 2001/11/14 10:50:46 cblume
-Changes in digits IO. Add merging of summable digits
-
-Revision 1.10 2001/05/31 16:53:26 alibrary
-Correction to the destructor
-
-Revision 1.9 2001/05/21 16:45:47 hristov
-Last minute changes (C.Blume)
-
-Revision 1.8 2001/01/26 19:56:57 hristov
-Major upgrade of AliRoot code
-
-Revision 1.7 2000/12/20 13:00:45 cblume
-Modifications for the HP-compiler
-
-Revision 1.6 2000/12/12 10:20:10 cblume
-Initialize fSepctrum = 0 in ctors
-
-Revision 1.5 2000/10/15 23:40:01 cblume
-Remove AliTRDconst
-
-Revision 1.4 2000/10/06 16:49:46 cblume
-Made Getters const
-
-Revision 1.3.2.1 2000/09/18 13:45:30 cblume
-New class AliTRDsim that simulates TR photons
-
-Revision 1.2 1999/09/29 09:24:35 fca
-Introduction of the Copyright and cvs Log
-
-*/
+/* $Id$ */
///////////////////////////////////////////////////////////////////////////////
// //
// 17.07.1998 - A.Andronic //
// 08.12.1998 - simplified version //
// 11.07.2000 - Adapted code to aliroot environment (C.Blume) //
+// 04.06.2004 - Momentum dependent parameters implemented (CBL) //
// //
///////////////////////////////////////////////////////////////////////////////
fSpectrum = 0;
fSigma = 0;
+ fNFoils = 0;
+ fNFoilsUp = 0;
Init();
fSpectrum = 0;
fSigma = 0;
+ fNFoils = 0;
+ fNFoilsUp = 0;
Init();
}
//_____________________________________________________________________________
-AliTRDsim::AliTRDsim(const AliTRDsim &s)
+AliTRDsim::AliTRDsim(const AliTRDsim &s):TObject(s)
{
//
// AliTRDsim copy constructor
// if (fSpectrum) delete fSpectrum;
if (fSigma) delete [] fSigma;
+ if (fNFoils) delete [] fNFoils;
+ if (fNFoilsUp) delete [] fNFoilsUp;
}
}
//_____________________________________________________________________________
-void AliTRDsim::Copy(TObject &s)
+void AliTRDsim::Copy(TObject &s) const
{
//
// Copy function
//
- ((AliTRDsim &) s).fNFoils = fNFoils;
((AliTRDsim &) s).fFoilThick = fFoilThick;
((AliTRDsim &) s).fFoilDens = fFoilDens;
((AliTRDsim &) s).fFoilOmega = fFoilOmega;
((AliTRDsim &) s).fSpLower = fSpLower;
((AliTRDsim &) s).fSpUpper = fSpUpper;
- if (((AliTRDsim &) s).fSigma) delete [] ((AliTRDsim &) s).fSigma;
+ if (((AliTRDsim &) s).fNFoils) delete [] ((AliTRDsim &) s).fNFoils;
+ ((AliTRDsim &) s).fNFoils = new Int_t[fNFoilsDim];
+ for (Int_t iFoil = 0; iFoil < fNFoilsDim; iFoil++) {
+ ((AliTRDsim &) s).fNFoils[iFoil] = fNFoils[iFoil];
+ }
+
+ if (((AliTRDsim &) s).fNFoilsUp) delete [] ((AliTRDsim &) s).fNFoilsUp;
+ ((AliTRDsim &) s).fNFoilsUp = new Double_t[fNFoilsDim];
+ for (Int_t iFoil = 0; iFoil < fNFoilsDim; iFoil++) {
+ ((AliTRDsim &) s).fNFoilsUp[iFoil] = fNFoilsUp[iFoil];
+ }
+
+ if (((AliTRDsim &) s).fSigma) delete [] ((AliTRDsim &) s).fSigma;
((AliTRDsim &) s).fSigma = new Double_t[fSpNBins];
for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
((AliTRDsim &) s).fSigma[iBin] = fSigma[iBin];
{
//
// Initialization
- // The default radiator are 100 prolypropilene foils of 13 mu thickness
- // with gaps of 60 mu filled with CO2.
+ // The default radiator are prolypropilene foils of 10 mu thickness
+ // with gaps of 80 mu filled with N2.
//
- fNFoils = 100;
+ fNFoilsDim = 7;
+
+ if (fNFoils) delete [] fNFoils;
+ fNFoils = new Int_t[fNFoilsDim];
+ fNFoils[0] = 170;
+ fNFoils[1] = 225; //250;
+ fNFoils[2] = 275; //310;
+ fNFoils[3] = 305; //380;
+ fNFoils[4] = 325; //430;
+ fNFoils[5] = 340; //490;
+ fNFoils[6] = 350; //550;
+
+ if (fNFoilsUp) delete [] fNFoilsUp;
+ fNFoilsUp = new Double_t[fNFoilsDim];
+ fNFoilsUp[0] = 1.25;
+ fNFoilsUp[1] = 1.75;
+ fNFoilsUp[2] = 2.50;
+ fNFoilsUp[3] = 3.50;
+ fNFoilsUp[4] = 4.50;
+ fNFoilsUp[5] = 5.50;
+ fNFoilsUp[6] = 10000.0;
fFoilThick = 0.0013;
fFoilDens = 0.92;
fFoilOmega = Omega(fFoilDens,fFoilZ,fFoilA);
fGapThick = 0.0060;
- fGapDens = 0.001977;
- fGapZ = 7.45455;
- fGapA = 14.9091;
+ fGapDens = 0.00125;
+ fGapZ = 7.0;
+ fGapA = 14.00674;
fGapOmega = Omega(fGapDens ,fGapZ ,fGapA );
fTemp = 293.16;
break;
};
- // Calculate gamma
- Double_t gamma = TMath::Sqrt(p*p + mass*mass) / mass;
-
// Calculate the TR photons
- return TrPhotons(gamma, nPhoton, ePhoton);
+ return TrPhotons(p, mass, nPhoton, ePhoton);
}
//_____________________________________________________________________________
-Int_t AliTRDsim::TrPhotons(Double_t gamma, Int_t &nPhoton, Float_t *ePhoton)
+Int_t AliTRDsim::TrPhotons(Float_t p, Float_t mass
+ , Int_t &nPhoton, Float_t *ePhoton)
{
//
// Produces TR photons.
//
const Double_t kAlpha = 0.0072973;
- const Int_t kSumMax = 10;
+ const Int_t kSumMax = 30;
+
+ Double_t tau = fGapThick / fFoilThick;
- Double_t kappa = fGapThick / fFoilThick;
+ // Calculate gamma
+ Double_t gamma = TMath::Sqrt(p*p + mass*mass) / mass;
+
+ // Select the number of foils corresponding to momentum
+ Int_t foils = SelectNFoils(p);
fSpectrum->Reset();
// The TR spectrum
- Double_t stemp = 0;
- for (Int_t iBin = 0; iBin < fSpNBins; iBin++) {
+ Double_t csi1,csi2,rho1,rho2;
+ Double_t fSigma,Sum,Nequ,theta_n,aux;
+ Double_t energyeV, energykeV;
+ for (Int_t iBin = 1; iBin <= fSpNBins; iBin++) {
+ energykeV = fSpectrum->GetBinCenter(iBin);
+ energyeV = energykeV * 1.e3;
- // keV -> eV
- Double_t energyeV = (fSpBinWidth * iBin + 1.0) * 1e3;
+ fSigma = Sigma(energykeV);
- Double_t csFoil = fFoilOmega / energyeV;
- Double_t csGap = fGapOmega / energyeV;
+ csi1 = fFoilOmega / energyeV;
+ csi2 = fGapOmega / energyeV;
- Double_t rho1 = energyeV * fFoilThick * 1e4 * 2.5
- * (1.0 / (gamma*gamma) + csFoil*csFoil);
- Double_t rho2 = energyeV * fFoilThick * 1e4 * 2.5
- * (1.0 / (gamma*gamma) + csGap *csGap);
+ rho1 = 2.5 * energyeV * fFoilThick * 1.E4
+ * (1. / (gamma*gamma) + csi1*csi1);
+ rho2 = 2.5 * energyeV * fFoilThick * 1.E4
+ * (1.0 / (gamma*gamma) + csi2 *csi2);
// Calculate the sum
- Double_t sum = 0;
- for (Int_t iSum = 0; iSum < kSumMax; iSum++) {
- Double_t tetan = (TMath::Pi() * 2.0 * (iSum+1) - (rho1 + kappa * rho2))
- / (kappa + 1.0);
- if (tetan < 0.0) tetan = 0.0;
- Double_t aux = 1.0 / (rho1 + tetan) - 1.0 / (rho2 + tetan);
- sum += tetan * (aux*aux) * (1.0 - TMath::Cos(rho1 + tetan));
+ Sum = 0;
+ for (Int_t n = 1; n <= kSumMax; n++) {
+ theta_n = (TMath::Pi() * 2.0 * n - (rho1 + tau * rho2)) / (1.+ tau);
+ if (theta_n < 0.) theta_n = 0.0;
+ aux = 1. / (rho1 + theta_n) - 1. / (rho2 + theta_n);
+ Sum += theta_n * (aux*aux) * (1.0 - TMath::Cos(rho1 + theta_n));
}
- // Absorbtion
- Double_t conv = 1.0 - TMath::Exp(-fNFoils * fSigma[iBin]);
+ // Equivalent number of foils
+ Nequ = (1. - TMath::Exp(-foils * fSigma)) / (1.- TMath::Exp(-fSigma));
- // eV -> keV
- Float_t energykeV = energyeV * 0.001;
// dN / domega
- Double_t wn = kAlpha * 4.0 / (fSigma[iBin] * (kappa + 1.0))
- * conv * sum / energykeV;
- fSpectrum->SetBinContent(iBin,wn);
-
- stemp += wn;
-
+ fSpectrum->SetBinContent(iBin,4. * kAlpha * Nequ * Sum / (energykeV * (1. + tau)));
}
// <nTR> (binsize corr.)
- Float_t ntr = stemp * fSpBinWidth;
+ Float_t ntr = fSpBinWidth*fSpectrum->Integral();
// Number of TR photons from Poisson distribution with mean <ntr>
nPhoton = gRandom->Poisson(ntr);
// Energy of the TR photons
return 0;
}
+
+//_____________________________________________________________________________
+Int_t AliTRDsim::SelectNFoils(Float_t p)
+{
+ //
+ // Selects the number of foils corresponding to the momentum
+ //
+
+ Int_t foils = fNFoils[fNFoilsDim-1];
+
+ for (Int_t iFoil = 0; iFoil < fNFoilsDim; iFoil++) {
+ if (p < fNFoilsUp[iFoil]) {
+ foils = fNFoils[iFoil];
+ break;
+ }
+ }
+
+ return foils;
+
+}