#include <HFCherenkov.h>
Public Member Functions | |
void | clearVectors () |
int | computeNPE (G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose) |
int | computeNPEinPMT (G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length) |
int | computeNPhTrapped (double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose) |
std::vector< double > | getMom () |
std::vector< double > | getWL () |
std::vector< double > | getWLAtten () |
std::vector< double > | getWLHEM () |
std::vector< double > | getWLIni () |
std::vector< double > | getWLQEff () |
std::vector< double > | getWLTrap () |
HFCherenkov (edm::ParameterSet const &p) | |
double | smearNPE (G4int Npe) |
virtual | ~HFCherenkov () |
Private Member Functions | |
double | computeHEMEff (double wavelength) |
int | computeNbOfPhotons (double pBeta, double step_length) |
double | computeQEff (double wavelength) |
bool | isApplicable (const G4ParticleDefinition *aParticleType) |
Private Attributes | |
double | aperture |
double | apertureTrap |
bool | checkSurvive |
double | gain |
double | lambda1 |
double | lambda2 |
std::vector< double > | momZ |
G4ThreeVector | phMom |
double | ref_index |
std::vector< double > | wl |
std::vector< double > | wlatten |
std::vector< double > | wlhem |
std::vector< double > | wlini |
std::vector< double > | wlqeff |
std::vector< double > | wltrap |
Definition at line 18 of file HFCherenkov.h.
HFCherenkov::HFCherenkov | ( | edm::ParameterSet const & | p | ) |
Definition at line 13 of file HFCherenkov.cc.
References aperture, apertureTrap, checkSurvive, clearVectors(), funct::cos(), gain, edm::ParameterSet::getParameter(), lambda1, lambda2, funct::pow(), and ref_index.
{ //Simple configurables //static SimpleConfigurable<double> p1(1.459, "HFCherenkov:RefIndex"); //static SimpleConfigurable<double> p2(280.0, "HFCherenkov:Lambda1"); //static SimpleConfigurable<double> p3(700.0, "HFCherenkov:Lambda2"); //static SimpleConfigurable<double> p4(0.33, "HFCherenkov:Aperture"); //static SimpleConfigurable<double> p5(0.22, "HFCherenkov:ApertureTrapped"); //static SimpleConfigurable<double> p6(0.33, "HFCherenkov:Gain"); //static SimpleConfigurable<bool> p7(false, "HFCherenkov:CheckSurvive"); ref_index = m_HF.getParameter<double>("RefIndex"); lambda1 = ((m_HF.getParameter<double>("Lambda1"))/pow(double(10),7))*cm; lambda2 = ((m_HF.getParameter<double>("Lambda2"))/pow(double(10),7))*cm; aperture = cos(asin(m_HF.getParameter<double>("Aperture"))); apertureTrap = cos(asin(m_HF.getParameter<double>("ApertureTrapped"))); gain = m_HF.getParameter<double>("Gain"); checkSurvive = m_HF.getParameter<bool>("CheckSurvive"); edm::LogInfo("HFShower") << "HFCherenkov:: initialised with ref_index " << ref_index << " lambda1/lambda2 (cm) " << lambda1/cm << "/" << lambda2/cm << " aperture(total/trapped) " << aperture << "/" << apertureTrap << " Check photon survival in HF " << checkSurvive << " Gain " << gain; clearVectors(); }
HFCherenkov::~HFCherenkov | ( | ) | [virtual] |
Definition at line 42 of file HFCherenkov.cc.
{}
void HFCherenkov::clearVectors | ( | ) |
Definition at line 394 of file HFCherenkov.cc.
References momZ, wl, wlatten, wlhem, wlini, wlqeff, and wltrap.
Referenced by computeNPE(), computeNPEinPMT(), and HFCherenkov().
double HFCherenkov::computeHEMEff | ( | double | wavelength | ) | [private] |
Definition at line 349 of file HFCherenkov.cc.
References LogDebug.
Referenced by computeNPE().
{ double hEMEff = 0; if (wavelength < 400.) { hEMEff = 0.; } else if (wavelength >= 400. && wavelength < 410.) { //hEMEff = .99 * (wavelength - 400.) / 10.; hEMEff = (-1322.453 / wavelength) + 4.2056; } else if (wavelength >= 410.) { hEMEff = 0.99; if (wavelength > 415. && wavelength < 445.) { //abs(wavelength - 430.) < 15. //hEMEff = 0.95; hEMEff = 0.97; } else if (wavelength > 550. && wavelength < 600.) { // abs(wavelength - 575.) < 25.) //hEMEff = 0.96; hEMEff = 0.98; } else if (wavelength > 565. && wavelength <= 635.) { // added later // abs(wavelength - 600.) < 35.) hEMEff = (701.7268 / wavelength) - 0.186; } } #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength " << wavelength << " hEMEff " << hEMEff; #endif return hEMEff; }
int HFCherenkov::computeNbOfPhotons | ( | double | pBeta, |
double | step_length | ||
) | [private] |
Definition at line 316 of file HFCherenkov.cc.
References alpha, beta, lambda1, lambda2, LogDebug, M_PI, ref_index, and funct::sin().
Referenced by computeNPE(), computeNPEinPMT(), and computeNPhTrapped().
{ double pBeta = beta; double alpha = 0.0073; double step_length = stepL; double theta_C = acos(1./(pBeta*ref_index)); double lambdaDiff = (1./lambda1 - 1./lambda2); double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm; double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) * (step_length/cm); int nbOfPhotons = int(d_NOfPhotons); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength " << step_length << " theta_C " << theta_C << " lambdaDiff " << lambdaDiff << " cherenPhPerLength " << cherenPhPerLength << " Photons " << d_NOfPhotons << " " << nbOfPhotons; #endif return nbOfPhotons; }
int HFCherenkov::computeNPE | ( | G4ParticleDefinition * | pDef, |
double | pBeta, | ||
double | u, | ||
double | v, | ||
double | w, | ||
double | step_length, | ||
double | zFiber, | ||
double | Dose, | ||
int | Npe_Dose | ||
) |
Definition at line 85 of file HFCherenkov.cc.
References aperture, checkSurvive, clearVectors(), computeHEMEff(), computeNbOfPhotons(), computeQEff(), funct::cos(), funct::exp(), i, isApplicable(), lambda1, lambda2, LogDebug, M_PI, momZ, funct::pow(), ref_index, funct::sin(), mathSSE::sqrt(), wl, wlatten, wlhem, wlini, wlqeff, and wltrap.
Referenced by HFShower::getHits().
{ clearVectors(); if (!isApplicable(pDef)) {return 0;} if (pBeta < (1/ref_index) || step_length < 0.0001) { #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta << " 1/mu " << (1/ref_index) << " step_length " << step_length; #endif return 0; } double uv = sqrt(u*u + v*v); int nbOfPhotons = computeNbOfPhotons(pBeta, step_length); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta << " u/v/w " << u << "/" << v << "/" << w << " step_length " << step_length << " zFib " << zFiber << " nbOfPhotons " << nbOfPhotons; #endif if (nbOfPhotons < 0) { return 0; } else if (nbOfPhotons > 0) { double u_ph, v_ph, w_ph=0; for (int i = 0; i < nbOfPhotons; i++) { double rand = G4UniformRand(); double theta_C = acos(1./(pBeta*ref_index)); double phi_C = 2*M_PI*rand; double sinTheta = sin(theta_C); double cosTheta = cos(theta_C); double sinPhi = sin(phi_C); double cosPhi = cos(phi_C); //photon momentum if (uv < 0.001) { // aligned with z-axis u_ph = sinTheta * cosPhi; v_ph = sinTheta * sinPhi; w_ph = cosTheta; } else { // general case u_ph = u * cosTheta + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv; v_ph = v * cosTheta + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv; w_ph = w * cosTheta - sinTheta * cosPhi * uv; } double r_lambda = G4UniformRand(); double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda * (lambda2 - lambda1)); double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm wlini.push_back(lambda); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: " << i << " lambda " << lambda << " w_ph " << w_ph << " aperture " << aperture; #endif if (w_ph > aperture) { // phton trapped inside fiber wltrap.push_back(lambda); double prob_HF = 1.0; //photon survived in HF double a0_inv = 0.1234; //meter^-1 double prob_MX = exp( - 0.5 * a0_inv ); //light mixer if (checkSurvive) { double a_inv = a0_inv + 0.14 * pow(dose,0.30); double z_meters = zFiber; prob_HF = exp(-z_meters * a_inv ); //photon survived in HF } rand = G4UniformRand(); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: probHF " << prob_HF << " prob_MX " << prob_MX << " Random " << rand << " Survive? " << (rand < (prob_HF * prob_MX)); #endif if (rand < (prob_HF * prob_MX)) { // survived and sent to light mixer wlatten.push_back(lambda); rand = G4UniformRand(); double effHEM = computeHEMEff(lambda); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: w_ph " << w_ph << " effHEM " << effHEM << " Random " << rand << " Survive? " << (w_ph>0.997||(rand<effHEM)); #endif if (w_ph>0.997 || (rand<effHEM)) { // survived HEM wlhem.push_back(lambda); double qEffic = computeQEff(lambda); rand = G4UniformRand(); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: qEffic " << qEffic << " Random " << rand << " Survive? " <<(rand < qEffic); #endif if (rand < qEffic) { // made photoelectron npe_Dose += 1; momZ.push_back(w_ph); wl.push_back(lambda); wlqeff.push_back(lambda); } // made pe } // passed HEM } // passed fiber } // end of if(w_ph < w_aperture), trapped inside fiber }// end of ++NbOfPhotons } // end of if(NbOfPhotons)} int npe = npe_Dose; // Nb of photoelectrons #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPE: npe " << npe; #endif return npe; }
int HFCherenkov::computeNPEinPMT | ( | G4ParticleDefinition * | pDef, |
double | pBeta, | ||
double | u, | ||
double | v, | ||
double | w, | ||
double | step_length | ||
) |
Definition at line 193 of file HFCherenkov.cc.
References aperture, clearVectors(), computeNbOfPhotons(), computeQEff(), funct::cos(), i, isApplicable(), lambda1, lambda2, LogDebug, M_PI, momZ, funct::pow(), ref_index, funct::sin(), mathSSE::sqrt(), wl, wlatten, wlini, wlqeff, and wltrap.
Referenced by HFShowerFibreBundle::getHits(), and HFShowerPMT::getHits().
{ clearVectors(); int npe_ = 0; if (!isApplicable(pDef)) {return 0;} if (pBeta < (1/ref_index) || step_length < 0.0001) { #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta << " 1/mu " << (1/ref_index) << " step_length " << step_length; #endif return 0; } double uv = sqrt(u*u + v*v); int nbOfPhotons = computeNbOfPhotons(pBeta, step_length); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta << " u/v/w " << u << "/" << v << "/" << w << " step_length " << step_length << " nbOfPhotons " << nbOfPhotons; #endif if (nbOfPhotons < 0) { return 0; } else if (nbOfPhotons > 0) { double u_ph, v_ph, w_ph=0; for (int i = 0; i < nbOfPhotons; i++) { double rand = G4UniformRand(); double theta_C = acos(1./(pBeta*ref_index)); double phi_C = 2*M_PI*rand; double sinTheta = sin(theta_C); double cosTheta = cos(theta_C); double sinPhi = sin(phi_C); double cosPhi = cos(phi_C); //photon momentum if (uv < 0.001) { // aligned with z-axis u_ph = sinTheta * cosPhi; v_ph = sinTheta * sinPhi; w_ph = cosTheta; } else { // general case u_ph = u * cosTheta + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv; v_ph = v * cosTheta + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv; w_ph = w * cosTheta - sinTheta * cosPhi * uv; } double r_lambda = G4UniformRand(); double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda * (lambda2 - lambda1)); double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm wlini.push_back(lambda); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: " <<i <<" lambda " << lambda << " w_ph " << w_ph << " aperture " << aperture; #endif if (w_ph > aperture) { // phton trapped inside PMT glass wltrap.push_back(lambda); rand = G4UniformRand(); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: Random " << rand << " Survive? " << (rand < 1.); #endif if (rand < 1.0) { // survived all the times and sent to photo-cathode wlatten.push_back(lambda); rand = G4UniformRand(); double qEffic = computeQEff(lambda);//Quantum efficiency of the PMT rand = G4UniformRand(); #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: qEffic " << qEffic << " Random " << rand << " Survive? " <<(rand < qEffic); #endif if (rand < qEffic) { // made photoelectron npe_ += 1; momZ.push_back(w_ph); wl.push_back(lambda); wlqeff.push_back(lambda); } // made pe } // accepted all Cherenkov photons } // end of if(w_ph < w_aperture), trapped inside glass }// end of ++NbOfPhotons } // end of if(NbOfPhotons)} #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: npe " << npe_; #endif return npe_; }
int HFCherenkov::computeNPhTrapped | ( | double | pBeta, |
double | u, | ||
double | v, | ||
double | w, | ||
double | step_length, | ||
double | zFiber, | ||
double | Dose, | ||
int | Npe_Dose | ||
) |
Definition at line 44 of file HFCherenkov.cc.
References apertureTrap, computeNbOfPhotons(), funct::cos(), i, M_PI, ref_index, funct::sin(), and mathSSE::sqrt().
{ if (pBeta < (1/ref_index) || step_length < 0.0001) {return 0;} double uv = sqrt(u*u + v*v); int nbOfPhotons = computeNbOfPhotons(pBeta, step_length); if (nbOfPhotons < 0) { return 0; } else if (nbOfPhotons > 0) { double u_ph, v_ph, w_ph=0; for (int i = 0; i < nbOfPhotons; i++) { double rand = G4UniformRand(); double theta_C = acos(1./(pBeta*ref_index)); double phi_C = 2*M_PI*rand; double sinTheta = sin(theta_C); double cosTheta = cos(theta_C); double sinPhi = sin(phi_C); double cosPhi = cos(phi_C); //photon momentum if (uv < 0.001) { // aligned with z-axis u_ph = sinTheta * cosPhi; v_ph = sinTheta * sinPhi; w_ph = cosTheta; } else { // general case u_ph = u * cosTheta + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv; v_ph = v * cosTheta + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv; w_ph = w * cosTheta - sinTheta * cosPhi * uv; } if (w_ph > apertureTrap) { // phton trapped inside fiber npe_Dose += 1; } } } int n_photons = npe_Dose; return n_photons; }
double HFCherenkov::computeQEff | ( | double | wavelength | ) | [private] |
Definition at line 336 of file HFCherenkov.cc.
References funct::exp(), LogDebug, funct::pow(), and ExpressReco_HICollisions_FallBack::y.
Referenced by computeNPE(), and computeNPEinPMT().
{ double y = (wavelength - 275.) /180.; double func = 1. / (1. + 250.*pow((y/5.),4)); double qE_R7525 = 0.77 * y * exp(-y) * func ; double qeff = qE_R7525; #ifdef DebugLog LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength << " y/func " << y << "/" << func << " qeff " << qeff; #endif return qeff; }
std::vector< double > HFCherenkov::getMom | ( | ) |
Definition at line 311 of file HFCherenkov.cc.
Referenced by HFShower::getHits().
std::vector< double > HFCherenkov::getWL | ( | ) |
Definition at line 306 of file HFCherenkov.cc.
Referenced by HFShower::getHits().
std::vector< double > HFCherenkov::getWLAtten | ( | ) |
std::vector< double > HFCherenkov::getWLHEM | ( | ) |
std::vector< double > HFCherenkov::getWLIni | ( | ) |
std::vector< double > HFCherenkov::getWLQEff | ( | ) |
std::vector< double > HFCherenkov::getWLTrap | ( | ) |
bool HFCherenkov::isApplicable | ( | const G4ParticleDefinition * | aParticleType | ) | [private] |
Definition at line 405 of file HFCherenkov.cc.
Referenced by computeNPE(), and computeNPEinPMT().
double HFCherenkov::smearNPE | ( | G4int | Npe | ) |
Definition at line 379 of file HFCherenkov.cc.
double HFCherenkov::aperture [private] |
Definition at line 61 of file HFCherenkov.h.
Referenced by computeNPE(), computeNPEinPMT(), and HFCherenkov().
double HFCherenkov::apertureTrap [private] |
Definition at line 61 of file HFCherenkov.h.
Referenced by computeNPhTrapped(), and HFCherenkov().
bool HFCherenkov::checkSurvive [private] |
Definition at line 63 of file HFCherenkov.h.
Referenced by computeNPE(), and HFCherenkov().
double HFCherenkov::gain [private] |
Definition at line 62 of file HFCherenkov.h.
Referenced by HFCherenkov(), and smearNPE().
double HFCherenkov::lambda1 [private] |
Definition at line 60 of file HFCherenkov.h.
Referenced by computeNbOfPhotons(), computeNPE(), computeNPEinPMT(), and HFCherenkov().
double HFCherenkov::lambda2 [private] |
Definition at line 60 of file HFCherenkov.h.
Referenced by computeNbOfPhotons(), computeNPE(), computeNPEinPMT(), and HFCherenkov().
std::vector<double> HFCherenkov::momZ [private] |
Definition at line 67 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getMom().
G4ThreeVector HFCherenkov::phMom [private] |
Definition at line 65 of file HFCherenkov.h.
double HFCherenkov::ref_index [private] |
Definition at line 59 of file HFCherenkov.h.
Referenced by computeNbOfPhotons(), computeNPE(), computeNPEinPMT(), computeNPhTrapped(), and HFCherenkov().
std::vector<double> HFCherenkov::wl [private] |
Definition at line 66 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWL().
std::vector<double> HFCherenkov::wlatten [private] |
Definition at line 70 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLAtten().
std::vector<double> HFCherenkov::wlhem [private] |
Definition at line 71 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), and getWLHEM().
std::vector<double> HFCherenkov::wlini [private] |
Definition at line 68 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLIni().
std::vector<double> HFCherenkov::wlqeff [private] |
Definition at line 72 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLQEff().
std::vector<double> HFCherenkov::wltrap [private] |
Definition at line 69 of file HFCherenkov.h.
Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLTrap().