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";
83 G4String attribute =
"ReadOutName";
95 edm::LogInfo(
"HFShower") <<
"HFShowerParam: " <<gpar.size() <<
" gpar (cm)";
96 for (
unsigned int ig=0; ig<gpar.size(); ig++)
97 edm::LogInfo(
"HFShower") <<
"HFShowerParam: gpar[" << ig <<
"] = "
98 << gpar[ig]/cm <<
" cm";
100 edm::LogError(
"HFShower") <<
"HFShowerParam: cannot get filtered "
101 <<
" view for " << attribute <<
" matching " <<
name;
102 throw cms::Exception(
"Unknown",
"HFShowerParam") <<
"cannot match " << attribute
103 <<
" to " << name <<
"\n";
110 edm::LogInfo(
"HFShower") <<
"att. length used for (lambda=" << lambdaMean
121 emPDG = theParticleTable->FindParticle(
"e-")->GetPDGEncoding();
122 epPDG = theParticleTable->FindParticle(
"e+")->GetPDGEncoding();
123 gammaPDG = theParticleTable->FindParticle(
"gamma")->GetPDGEncoding();
133 G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
134 G4Track * track = aStep->GetTrack();
135 G4ThreeVector hitPoint = preStepPoint->GetPosition();
136 G4int particleCode = track->GetDefinition()->GetPDGEncoding();
138 G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(),hitPoint.y(),zv);
140 double pin = (preStepPoint->GetTotalEnergy())/GeV;
141 double zint = hitPoint.z();
146 << track->GetDefinition()->GetParticleName()
147 <<
" of energy " << pin <<
" GeV"
148 <<
" Pos x,y,z = " << hitPoint.x() <<
","
149 << hitPoint.y() <<
"," << zint <<
" (" << zz <<
","
150 << localPoint.z() <<
", "
151 << (localPoint.z()+0.5*
gpar[1]) <<
") Local "
154 std::vector<HFShowerParam::Hit> hits;
160 double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
161 double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).
z();
162 if (hitPoint.z() < 0) dirz *= -1.;
164 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits Momentum "
165 <<track->GetDynamicParticle()->GetMomentumDirection()
166 <<
" HitPoint " << hitPoint <<
" dirz " << dirz;
169 if (track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1/
ref_index) &&
170 aStep->GetTotalEnergyDeposit() > 0) other =
true;
174 if (particleCode ==
emPDG || particleCode ==
epPDG ||
175 particleCode ==
gammaPDG || other) {
182 }
else if ((track->GetDefinition()->GetPDGCharge() != 0) &&
184 edep = (aStep->GetTotalEnergyDeposit())/GeV;
188 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits edep = " << edep
189 <<
" weight " << weight <<
" final " <<edep*weight
190 <<
", Kill = " << kill <<
", pin = " << pin
191 <<
", edMin = " <<
edMin <<
" Other " << other;
198 for (
unsigned int i=0;
i<hitSL.size();
i++) {
205 if (npmt <= 0) ok =
false;
209 hit.
depth = hitSL[
i].depth;
210 hit.
time = hitSL[
i].time;
220 if (hit.
depth == 1) {
224 }
else if (hit.
depth == 2) {
232 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth "
233 << hit.
depth <<
" with edep " << hit.
edep
234 <<
" Time " << hit.
time;
240 for (
unsigned int i=0;
i<hitSL.size(); ++
i) {
242 G4ThreeVector pe_effect(hitSL[
i].
position.x(), hitSL[
i].position.y(),
243 hitSL[
i].position.z());
248 if (zv < 0. || zv >
gpar[1]) {
250 std::cout<<
"-#Zcut-HFShowerParam::getHits:z="<<zv<<
",m="<<gpar[1]<<std::endl;
257 edm::LogInfo(
"HFShower") <<
"HFShowerParam::getHits:#PMT= "
258 << npmt <<
",z = " << zv;
262 edm::LogInfo(
"HFShower") <<
"-#PMT=0 cut-HFShowerParam::"
263 <<
"getHits: npmt = " << npmt;
266 }
else if (npmt > 24) {
271 edm::LogInfo(
"HFShower") <<
"-SHORT cut-HFShowerParam::"
272 <<
"getHits:zMin=" << gpar[0];
278 edm::LogInfo(
"HFShower") <<
"HFShowerParam: npmt " << npmt
279 <<
" zv " <<
std::abs(pe_effect.z())
280 <<
":" << gpar[4] <<
":" << zv <<
":"
281 << gpar[0] <<
" ok " << ok <<
" depth "
285 if (G4UniformRand() > 0.5) depth = 2;
286 if (depth == 2 && zv < gpar[0]) ok =
false;
290 double r1 = G4UniformRand();
292 edm::LogInfo(
"HFShower") <<
"HFShowerParam:Distance to PMT (" <<npmt
293 <<
") " << dist <<
", exclusion flag "
298 double r2 = G4UniformRand();
300 edm::LogInfo(
"HFShower") <<
"HFShowerParam:Extra exclusion "
301 << r2 <<
">" << weight <<
" "
309 hit.
time = time + hitSL[
i].time;
318 if (hit.
depth == 1) {
322 }
else if (hit.
depth == 2) {
330 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth "
331 << hit.
depth <<
" with edep "
332 << hit.
edep <<
" Time " << hit.
time;
341 double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
346 if (npmt <= 0) ok =
false;
349 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits hitPoint " << hitPoint <<
" flag " <<
ok;
357 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth 1 with edep "
358 << edep <<
" Time " << tSlice <<
":" << time
373 edm::LogInfo(
"HFShower") <<
"HFShowerParam: Hit at depth 2 with edep "
374 << edep <<
" Time " << tSlice <<
":" << time
386 for (
unsigned int ii=0;
ii<hits.size(); ++
ii) {
388 if (zv > 12790)
edm::LogInfo(
"HFShower")<<
"HFShowerParam: Abnormal hit along "
390 << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
391 <<
" at " << hits[
ii].position <<
" zz "
392 << zv <<
" Edep " << edep <<
" due to "
393 <<track->GetDefinition()->GetParticleName()
394 <<
" time " << hit.
time;
398 track->SetTrackStatus(fStopAndKill);
399 G4TrackVector tv = *(aStep->GetSecondary());
400 for (
unsigned int kk=0;
kk<tv.size(); ++
kk) {
401 if (tv[
kk]->
GetVolume() == preStepPoint->GetPhysicalVolume())
402 tv[
kk]->SetTrackStatus(fStopAndKill);
406 edm::LogInfo(
"HFShower") <<
"HFShowerParam: getHits kill (" << kill
407 <<
") track " << track->GetTrackID()
408 <<
" at " << hitPoint
409 <<
" and deposit " << edep <<
" " << hits.size()
410 <<
" times" <<
" ZZ " << zz <<
" " <<
gpar[0];
421 LogDebug(
"HFShower") <<
"HFShowerParam:getDDDArray called for " << str;
429 const std::vector<double> & fvec = value.
doubles();
430 int nval = fvec.size();
433 <<
" bins " << nval <<
" < 2 ==> illegal";
434 throw cms::Exception(
"Unknown",
"HFShowerParam") <<
"nval < 2 for array "
439 edm::LogError(
"HFShower") <<
"HFShowerParam : cannot get array " << str;
440 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 addFilter(const DDFilter &, log_op op=AND)
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 position[TOTALCHAMBERS][3]
static int PMTNumber(const G4ThreeVector &pe_effect)
tuple path
else: Piece not in the list, fine.
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 &)
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)
DDsvalues_type mergedSpecifics() const
void initRun(G4ParticleTable *theParticleTable)
void initRun(G4ParticleTable *)
bool firstChild()
set the current node to the first child ...
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
volatile std::atomic< bool > shutdown_flag false
Power< A, B >::type pow(const A &a, const B &b)
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.