15 #include "G4VPhysicalVolume.hh"
18 #include "G4NavigationHistory.hh"
19 #include "Randomize.hh"
21 #include "CLHEP/Units/GlobalPhysicalConstants.h"
22 #include "CLHEP/Units/GlobalSystemOfUnits.h"
37 bool useShowerLibrary = m_HF.
getParameter<
bool>(
"UseShowerLibrary");
42 double lambdaMean = m_HF.
getParameter<
double>(
"LambdaMean");
46 edm::LogInfo(
"HFShower") <<
"HFShowerParam::Use of shower library is set to "
47 << useShowerLibrary <<
" Use of Gflash is set to "
48 << useGflash <<
" P.E. per GeV " << pePerGeV
49 <<
", ref. index of fibre " << ref_index
50 <<
", Track EM Flag " << trackEM <<
", edMin "
51 << edMin <<
" GeV, use of Short fibre info in"
52 <<
" shower library set to " << !(
onlyLong)
53 <<
", use of parametrization for last part set to "
54 << parametrizeLast <<
", Mean lambda " << lambdaMean
55 <<
", aperture (cutoff) " << aperture
62 LogDebug(
"HFShower") <<
"HFShowerParam::Save histos in directory "
63 <<
"ProfileFromParam";
65 hzvem = showerDir.
make<TH1F>(
"hzvem",
"Longitudinal Profile (EM Part);Number of PE",330,0.0,1650.0);
66 hzvhad = showerDir.
make<TH1F>(
"hzvhad",
"Longitudinal Profile (Had Part);Number of PE",330,0.0,1650.0);
67 em_2d_1 = showerDir.
make<TH2F>(
"em_2d_1",
"Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
68 em_long_1 = showerDir.
make<TH1F>(
"em_long_1",
"Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
69 em_long_1_tuned = showerDir.
make<TH1F>(
"em_long_1_tuned",
"Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
70 em_lateral_1 = showerDir.
make<TH1F>(
"em_lateral_1",
"Lateral Profile;cm;Events",100,50.0,150.0);
71 em_2d_2 = showerDir.
make<TH2F>(
"em_2d_2",
"Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
72 em_long_2 = showerDir.
make<TH1F>(
"em_long_2",
"Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
73 em_lateral_2 = showerDir.
make<TH1F>(
"em_lateral_2",
"Lateral Profile;cm;Events",100,50.0,150.0);
74 em_long_gflash = showerDir.
make<TH1F>(
"em_long_gflash",
"Longitudinal Profile From GFlash;cm;Number of Spots",800,800.0,1600.0);
75 em_long_sl = showerDir.
make<TH1F>(
"em_long_sl",
"Longitudinal Profile From Shower Library;cm;Number of Spots",800,800.0,1600.0);
78 edm::LogInfo(
"HFShower") <<
"HFShowerParam::No file is available for "
79 <<
"saving histos so the flag is set to false";
87 edm::LogInfo(
"HFShower") <<
"att. length used for (lambda=" << lambdaMean
99 emPDG = theParticleTable->FindParticle(
"e-")->GetPDGEncoding();
100 epPDG = theParticleTable->FindParticle(
"e+")->GetPDGEncoding();
101 gammaPDG = theParticleTable->FindParticle(
"gamma")->GetPDGEncoding();
111 edm::LogInfo(
"HFShower") <<
"HFShowerParam: " <<
gpar.size() <<
" gpar (cm)";
112 for (
unsigned int ig=0; ig<
gpar.size(); ig++)
113 edm::LogInfo(
"HFShower") <<
"HFShowerParam: gpar[" << ig <<
"] = "
114 <<
gpar[ig]/cm <<
" cm";
119 G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
120 G4Track * track = aStep->GetTrack();
121 G4ThreeVector hitPoint = preStepPoint->GetPosition();
122 G4int particleCode = track->GetDefinition()->GetPDGEncoding();
124 G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(),hitPoint.y(),zv);
126 double pin = (preStepPoint->GetTotalEnergy())/
GeV;
127 double zint = hitPoint.z();
132 << track->GetDefinition()->GetParticleName()
133 <<
" of energy " << pin <<
" GeV"
134 <<
" Pos x,y,z = " << hitPoint.x() <<
","
135 << hitPoint.y() <<
"," << zint <<
" (" << zz <<
","
136 << localPoint.z() <<
", "
137 << (localPoint.z()+0.5*
gpar[1]) <<
") Local "
140 std::vector<HFShowerParam::Hit> hits;
146 double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
147 double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).
z();
148 if (hitPoint.z() < 0) dirz *= -1.;
150 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits Momentum "
151 <<track->GetDynamicParticle()->GetMomentumDirection()
152 <<
" HitPoint " << hitPoint <<
" dirz " << dirz;
155 if (track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1/
ref_index) &&
156 aStep->GetTotalEnergyDeposit() > 0) other =
true;
160 if (particleCode ==
emPDG || particleCode ==
epPDG ||
161 particleCode ==
gammaPDG || other) {
168 }
else if ((track->GetDefinition()->GetPDGCharge() != 0) &&
170 edep = (aStep->GetTotalEnergyDeposit())/
GeV;
174 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits edep = " << edep
175 <<
" weight " << weight <<
" final " <<edep*weight
176 <<
", Kill = " << kill <<
", pin = " << pin
177 <<
", edMin = " <<
edMin <<
" Other " << other;
184 for (
unsigned int i=0;
i<hitSL.size();
i++) {
191 if (npmt <= 0) ok =
false;
195 hit.
depth = hitSL[
i].depth;
196 hit.
time = hitSL[
i].time;
206 if (hit.
depth == 1) {
210 }
else if (hit.
depth == 2) {
218 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth "
219 << hit.
depth <<
" with edep " << hit.
edep
220 <<
" Time " << hit.
time;
226 for (
unsigned int i=0;
i<hitSL.size(); ++
i) {
228 G4ThreeVector pe_effect(hitSL[
i].
position.x(), hitSL[
i].position.y(),
229 hitSL[
i].position.z());
234 if (zv < 0. || zv >
gpar[1]) {
236 std::cout<<
"-#Zcut-HFShowerParam::getHits:z="<<zv<<
",m="<<gpar[1]<<std::endl;
243 edm::LogInfo(
"HFShower") <<
"HFShowerParam::getHits:#PMT= "
244 << npmt <<
",z = " << zv;
248 edm::LogInfo(
"HFShower") <<
"-#PMT=0 cut-HFShowerParam::"
249 <<
"getHits: npmt = " << npmt;
252 }
else if (npmt > 24) {
257 edm::LogInfo(
"HFShower") <<
"-SHORT cut-HFShowerParam::"
258 <<
"getHits:zMin=" << gpar[0];
264 edm::LogInfo(
"HFShower") <<
"HFShowerParam: npmt " << npmt
265 <<
" zv " <<
std::abs(pe_effect.z())
266 <<
":" << gpar[4] <<
":" << zv <<
":"
267 << gpar[0] <<
" ok " << ok <<
" depth "
271 if (G4UniformRand() > 0.5) depth = 2;
272 if (depth == 2 && zv < gpar[0]) ok =
false;
276 double r1 = G4UniformRand();
278 edm::LogInfo(
"HFShower") <<
"HFShowerParam:Distance to PMT (" <<npmt
279 <<
") " << dist <<
", exclusion flag "
284 double r2 = G4UniformRand();
286 edm::LogInfo(
"HFShower") <<
"HFShowerParam:Extra exclusion "
287 << r2 <<
">" << weight <<
" "
295 hit.
time = time + hitSL[
i].time;
304 if (hit.
depth == 1) {
308 }
else if (hit.
depth == 2) {
316 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth "
317 << hit.
depth <<
" with edep "
318 << hit.
edep <<
" Time " << hit.
time;
327 double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
332 if (npmt <= 0) ok =
false;
335 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits hitPoint " << hitPoint <<
" flag " <<
ok;
339 hit.
time = tSlice+time;
343 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth 1 with edep "
344 << edep <<
" Time " << tSlice <<
":" << time
356 hit.
time = tSlice+time;
359 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth 2 with edep "
360 << edep <<
" Time " << tSlice <<
":" << time
372 for (
unsigned int ii=0;
ii<hits.size(); ++
ii) {
374 if (zv > 12790)
edm::LogInfo(
"HFShower")<<
"HFShowerParam: Abnormal hit along "
376 << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
377 <<
" at " << hits[
ii].position <<
" zz "
378 << zv <<
" Edep " << edep <<
" due to "
379 <<track->GetDefinition()->GetParticleName()
380 <<
" time " << hit.
time;
384 track->SetTrackStatus(fStopAndKill);
385 G4TrackVector tv = *(aStep->GetSecondary());
386 for (
unsigned int kk=0;
kk<tv.size(); ++
kk) {
387 if (tv[
kk]->
GetVolume() == preStepPoint->GetPhysicalVolume())
388 tv[
kk]->SetTrackStatus(fStopAndKill);
392 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits kill (" << kill
393 <<
") track " << track->GetTrackID()
394 <<
" at " << hitPoint
395 <<
" and deposit " << edep <<
" " << hits.size()
396 <<
" times" <<
" ZZ " << zz <<
" " <<
gpar[0];
407 LogDebug(
"HFShower") <<
"HFShowerParam:getDDDArray called for " << str;
415 const std::vector<double> & fvec = value.
doubles();
416 int nval = fvec.size();
419 <<
" bins " << nval <<
" < 2 ==> illegal";
420 throw cms::Exception(
"Unknown",
"HFShowerParam") <<
"nval < 2 for array "
425 edm::LogError(
"HFShower") <<
"HFShowerParam : cannot get array " << str;
426 throw cms::Exception(
"Unknown",
"HFShowerParam") <<
"cannot get array "
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
std::vector< Hit > gfParameterization(G4Step *aStep, bool &ok, bool onlyLong=false)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
void initRun(HcalDDDSimConstants *)
std::vector< Hit > getHits(G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
double attLength(double lambda)
std::vector< double > gpar
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
HFShowerLibrary * showerLibrary
HFShowerParam(std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
type of data representation of DDCompactView
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
static int PMTNumber(const G4ThreeVector &pe_effect)
Cos< T >::type cos(const T &t)
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Abs< T >::type abs(const T &t)
T * make(const Args &...args) const
make new ROOT object
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
const std::vector< double > & getGparHF() const
static const G4LogicalVolume * GetVolume(const std::string &name)
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
std::vector< Hit > getHits(G4Step *aStep, double weight)
void initRun(G4ParticleTable *, HcalDDDSimConstants *)
static int position[264][3]
volatile std::atomic< bool > shutdown_flag false
void initRun(G4ParticleTable *, HcalDDDSimConstants *)
Power< A, B >::type pow(const A &a, const B &b)