8 #include "G4Poisson.hh"
9 #include "G4ParticleDefinition.hh"
32 edm::LogInfo(
"HFShower") <<
"HFCherenkov:: initialised with ref_index "
34 << lambda1/cm <<
"/" <<
lambda2/cm
35 <<
" aperture(total/trapped) " <<
aperture <<
"/"
37 << checkSurvive <<
" Gain " <<
gain;
45 double u,
double v,
double w,
46 double step_length,
double zFiber,
47 double dose,
int npe_Dose) {
49 if (pBeta < (1/
ref_index) || step_length < 0.0001) {
return 0;}
51 double uv =
sqrt(u*u + v*v);
54 if (nbOfPhotons < 0) {
56 }
else if (nbOfPhotons > 0) {
58 for (
int i = 0;
i < nbOfPhotons;
i++) {
59 double rand = G4UniformRand();
60 double theta_C = acos(1./(pBeta*
ref_index));
62 double sinTheta =
sin(theta_C);
63 double cosTheta =
cos(theta_C);
64 double cosPhi =
cos(phi_C);
69 w_ph = w * cosTheta - sinTheta * cosPhi * uv;
76 int n_photons = npe_Dose;
81 double u,
double v,
double w,
82 double step_length,
double zFiber,
83 double dose,
int npe_Dose) {
87 if (pBeta < (1/
ref_index) || step_length < 0.0001) {
89 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: pBeta " << pBeta
90 <<
" 1/mu " << (1/
ref_index) <<
" step_length "
96 double uv =
sqrt(u*u + v*v);
99 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: pBeta " << pBeta
100 <<
" u/v/w " << u <<
"/" << v <<
"/" << w
101 <<
" step_length " << step_length <<
" zFib " << zFiber
102 <<
" nbOfPhotons " << nbOfPhotons;
104 if (nbOfPhotons < 0) {
106 }
else if (nbOfPhotons > 0) {
108 for (
int i = 0;
i < nbOfPhotons;
i++) {
109 double rand = G4UniformRand();
110 double theta_C = acos(1./(pBeta*
ref_index));
112 double sinTheta =
sin(theta_C);
113 double cosTheta =
cos(theta_C);
114 double cosPhi =
cos(phi_C);
119 w_ph = w * cosTheta - sinTheta * cosPhi * uv;
121 double r_lambda = G4UniformRand();
124 double lambda = (lambda0/cm) *
pow(
double(10),7);
125 wlini.push_back(lambda);
127 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: " <<
i <<
" lambda "
128 << lambda <<
" w_ph " << w_ph <<
" aperture "
131 if (w_ph > aperture) {
133 double prob_HF = 1.0;
134 double a0_inv = 0.1234;
135 double prob_MX =
exp( - 0.5 * a0_inv );
137 double a_inv = a0_inv + 0.14 *
pow(dose,0.30);
138 double z_meters = zFiber;
139 prob_HF =
exp(-z_meters * a_inv );
141 rand = G4UniformRand();
143 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: probHF " << prob_HF
144 <<
" prob_MX " << prob_MX <<
" Random " << rand
145 <<
" Survive? " << (rand < (prob_HF * prob_MX));
147 if (rand < (prob_HF * prob_MX)) {
149 rand = G4UniformRand();
152 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: w_ph " << w_ph
153 <<
" effHEM " << effHEM <<
" Random " << rand
154 <<
" Survive? " << (w_ph>0.997||(rand<effHEM));
156 if (w_ph>0.997 || (rand<effHEM)) {
157 wlhem.push_back(lambda);
159 rand = G4UniformRand();
161 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: qEffic "
162 << qEffic <<
" Random " << rand
163 <<
" Survive? " <<(rand < qEffic);
167 momZ.push_back(w_ph);
168 wl.push_back(lambda);
178 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPE: npe " << npe;
184 double u,
double v,
double w,
189 if (pBeta < (1/
ref_index) || step_length < 0.0001) {
191 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPEinPMT: pBeta " << pBeta
192 <<
" 1/mu " << (1/
ref_index) <<
" step_length "
198 double uv =
sqrt(u*u + v*v);
201 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPEinPMT: pBeta " << pBeta
202 <<
" u/v/w " << u <<
"/" << v <<
"/" << w
203 <<
" step_length " << step_length
204 <<
" nbOfPhotons " << nbOfPhotons;
206 if (nbOfPhotons < 0) {
208 }
else if (nbOfPhotons > 0) {
210 for (
int i = 0;
i < nbOfPhotons;
i++) {
211 double rand = G4UniformRand();
212 double theta_C = acos(1./(pBeta*
ref_index));
214 double sinTheta =
sin(theta_C);
215 double cosTheta =
cos(theta_C);
216 double cosPhi =
cos(phi_C);
221 w_ph = w * cosTheta - sinTheta * cosPhi * uv;
223 double r_lambda = G4UniformRand();
226 double lambda = (lambda0/cm) *
pow(
double(10),7);
227 wlini.push_back(lambda);
229 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPEinPMT: " <<
i <<
" lambda "
230 << lambda <<
" w_ph " << w_ph <<
" aperture "
233 if (w_ph > aperture) {
235 rand = G4UniformRand();
237 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPEinPMT: Random " << rand
238 <<
" Survive? " << (rand < 1.);
242 rand = G4UniformRand();
244 rand = G4UniformRand();
246 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPEinPMT: qEffic "
247 << qEffic <<
" Random " << rand
248 <<
" Survive? " <<(rand < qEffic);
252 momZ.push_back(w_ph);
253 wl.push_back(lambda);
261 LogDebug(
"HFShower") <<
"HFCherenkov::computeNPEinPMT: npe " << npe_;
267 std::vector<double>
v =
wlini;
272 std::vector<double>
v =
wltrap;
282 std::vector<double>
v =
wlhem;
287 std::vector<double>
v =
wlqeff;
292 std::vector<double>
v =
wl;
297 std::vector<double>
v =
momZ;
304 double alpha = 0.0073;
305 double step_length = stepL;
306 double theta_C = acos(1./(pBeta*
ref_index));
308 double cherenPhPerLength = 2 *
M_PI * alpha * lambdaDiff*cm;
309 double d_NOfPhotons = cherenPhPerLength *
sin(theta_C)*
sin(theta_C) * (step_length/cm);
310 int nbOfPhotons = int(d_NOfPhotons);
312 LogDebug(
"HFShower") <<
"HFCherenkov::computeNbOfPhotons: StepLength "
313 << step_length <<
" theta_C " << theta_C
314 <<
" lambdaDiff " << lambdaDiff
315 <<
" cherenPhPerLength " << cherenPhPerLength
316 <<
" Photons " << d_NOfPhotons <<
" " << nbOfPhotons;
323 double y = (wavelength - 275.) /180.;
324 double func = 1. / (1. + 250.*
pow((y/5.),4));
325 double qE_R7525 = 0.77 * y *
exp(-y) * func ;
326 double qeff = qE_R7525;
328 LogDebug(
"HFShower") <<
"HFCherenkov::computeQEff: wavelength " << wavelength
329 <<
" y/func " << y <<
"/" << func <<
" qeff " << qeff;
337 if (wavelength < 400.) {
339 }
else if (wavelength >= 400. && wavelength < 410.) {
341 hEMEff = (-1322.453 / wavelength) + 4.2056;
342 }
else if (wavelength >= 410.) {
344 if (wavelength > 415. && wavelength < 445.) {
348 }
else if (wavelength > 550. && wavelength < 600.) {
352 }
else if (wavelength > 565. && wavelength <= 635.) {
354 hEMEff = (701.7268 / wavelength) - 0.186;
358 LogDebug(
"HFShower") <<
"HFCherenkov::computeHEMEff: wavelength "
359 << wavelength <<
" hEMEff " << hEMEff;
368 for (
int i = 0;
i < npe; ++
i) {
369 double val = G4Poisson(
gain);
370 pe += (val/
gain) + 0.001*G4UniformRand();
374 LogDebug(
"HFShower") <<
"HFCherenkov::smearNPE: npe " << npe <<
" pe " << pe;
391 bool tmp = (aParticleType->GetPDGCharge() != 0);
393 LogDebug(
"HFShower") <<
"HFCherenkov::isApplicable: aParticleType "
394 << aParticleType->GetParticleName() <<
" PDGCharge "
395 << aParticleType->GetPDGCharge() <<
" Result " <<
tmp;
T getParameter(std::string const &) const
double computeQEff(double wavelength)
Sin< T >::type sin(const T &t)
HFCherenkov(edm::ParameterSet const &p)
int computeNPE(G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose)
std::vector< double > wlqeff
int computeNPEinPMT(G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
std::vector< double > wlini
double computeHEMEff(double wavelength)
double smearNPE(G4int Npe)
std::vector< double > getWL()
std::vector< double > getWLQEff()
std::vector< double > wltrap
Cos< T >::type cos(const T &t)
std::vector< double > wlhem
std::vector< double > getMom()
std::vector< double > getWLIni()
std::vector< std::vector< double > > tmp
int computeNbOfPhotons(double pBeta, double step_length)
std::vector< double > getWLTrap()
std::vector< double > getWLAtten()
bool isApplicable(const G4ParticleDefinition *aParticleType)
std::vector< double > wlatten
std::vector< double > momZ
int computeNPhTrapped(double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose)
Power< A, B >::type pow(const A &a, const B &b)
std::vector< double > getWLHEM()