CMS 3D CMS Logo

HFShowerPMT.cc
Go to the documentation of this file.
1 // File: HFShowerPMT.cc
3 // Description: Parametrized version of HF hits
5 
10 
12 
13 #include "G4NavigationHistory.hh"
14 #include "G4VPhysicalVolume.hh"
15 #include "G4Step.hh"
16 #include "G4Track.hh"
17 #include "CLHEP/Units/GlobalPhysicalConstants.h"
18 #include "CLHEP/Units/GlobalSystemOfUnits.h"
19 
20 //#define DebugLog
21 
23  edm::ParameterSet const & p) : cherenkov(nullptr) {
24 
25  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShowerPMT");
26  pePerGeV = m_HF.getParameter<double>("PEPerGeVPMT");
27 
28  //Special Geometry parameters
29  std::string attribute = "Volume";
30  std::string value = "HFPMT";
31  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
32  DDFilteredView fv1(cpv,filter1);
33  if (fv1.firstChild()) {
35  std::vector<double> neta;
36  neta = getDDDArray("indexPMTR",sv1);
37  for (unsigned int ii=0; ii<neta.size(); ii++) {
38  int index = static_cast<int>(neta[ii]);
39  int ir=-1, ifib=-1;
40  if (index >= 0) {
41  ir = index/10; ifib = index%10;
42  }
43  pmtR1.push_back(ir);
44  pmtFib1.push_back(ifib);
45  }
46  neta = getDDDArray("indexPMTL",sv1);
47  for (unsigned int ii=0; ii<neta.size(); ii++) {
48  int index = static_cast<int>(neta[ii]);
49  int ir=-1, ifib=-1;
50  if (index >= 0) {
51  ir = index/10; ifib = index%10;
52  }
53  pmtR2.push_back(ir);
54  pmtFib2.push_back(ifib);
55  }
56  edm::LogInfo("HFShower") << "HFShowerPMT: gets the Index matches for "
57  << neta.size() << " PMTs";
58  for (unsigned int ii=0; ii<neta.size(); ii++)
59  edm::LogInfo("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = "
60  << pmtR1[ii] << " fibreR[" << ii << "] = "
61  << pmtFib1[ii] << " rIndexL[" << ii << "] = "
62  << pmtR2[ii] << " fibreL[" << ii << "] = "
63  << pmtFib2[ii];
64  } else {
65  edm::LogWarning("HFShower") << "HFShowerPMT: cannot get filtered "
66  << " view for " << attribute << " matching "
67  << value;
68  }
69 
70  cherenkov = new HFCherenkov(m_HF);
71 }
72 
74  if (cherenkov) delete cherenkov;
75 }
76 
77 void HFShowerPMT::initRun(G4ParticleTable *, HcalDDDSimConstants* hcons) {
78 
79  // Special Geometry parameters
80  rTable = hcons->getRTableHF();
81  edm::LogInfo("HFShower") << "HFShowerPMT: " << rTable.size()
82  << " rTable (cm)";
83  for (unsigned int ig=0; ig<rTable.size(); ig++)
84  edm::LogInfo("HFShower") << "HFShowerPMT: rTable[" << ig << "] = "
85  << rTable[ig]/cm << " cm";
86 }
87 
88 double HFShowerPMT::getHits(G4Step * aStep) {
89 
90  indexR = indexF = -1;
91 
92  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
93  const G4VTouchable* touch = preStepPoint->GetTouchable();
94  int boxNo = touch->GetReplicaNumber(2);
95  int pmtNo = touch->GetReplicaNumber(1);
96  if (boxNo <= 1) {
97  indexR = pmtR1[pmtNo-1];
98  indexF = pmtFib1[pmtNo-1];
99  } else {
100  indexR = pmtR2[pmtNo-1];
101  indexF = pmtFib2[pmtNo-1];
102  }
103 
104 #ifdef DebugLog
105  double edep = aStep->GetTotalEnergyDeposit();
106  LogDebug("HFShower") << "HFShowerPMT: Box " << boxNo << " PMT "
107  << pmtNo << " Mapped Indices " << indexR << ", "
108  << indexF << " Edeposit " << edep/MeV << " MeV; PE "
109  << edep*pePerGeV/GeV;
110 #endif
111 
112  double photons = 0;
113  if (indexR >= 0 && indexF > 0) {
114  G4Track *aTrack = aStep->GetTrack();
115  G4ParticleDefinition *particleDef = aTrack->GetDefinition();
116  double stepl = aStep->GetStepLength();
117  double beta = preStepPoint->GetBeta();
118  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
119  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
120  GetTopTransform().TransformAxis(pDir);
121  photons = cherenkov->computeNPEinPMT(particleDef, beta, localMom.x(),
122  localMom.y(), localMom.z(), stepl);
123 #ifdef DebugLog
124  LogDebug("HFShower") << "HFShowerPMT::getHits: for particle "
125  << particleDef->GetParticleName() << " Step " << stepl
126  << " Beta " << beta << " Direction " << pDir
127  << " Local " << localMom << " p.e. " << photons;
128 #endif
129 
130  }
131  return photons;
132 }
133 
135 
136  double r = 0.;
137  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
138  r = 0.5*(rTable[indexR]+rTable[indexR+1]);
139 #ifdef DebugLog
140  else
141  LogDebug("HFShower") << "HFShowerPMT::getRadius: R " << indexR
142  << " F " << indexF;
143 #endif
144  if (indexF == 2) r =-r;
145 #ifdef DebugLog
146  LogDebug("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF
147  << ") " << r;
148 #endif
149  return r;
150 }
151 
152 std::vector<double> HFShowerPMT::getDDDArray(const std::string & str,
153  const DDsvalues_type & sv) {
154 
155 #ifdef DebugLog
156  LogDebug("HFShower") << "HFShowerPMT:getDDDArray called for " << str;
157 #endif
158  DDValue value(str);
159  if (DDfetch(&sv,value)) {
160 #ifdef DebugLog
161  LogDebug("HFShower") << value;
162 #endif
163  const std::vector<double> & fvec = value.doubles();
164  int nval = fvec.size();
165  if (nval < 2) {
166  edm::LogError("HFShower") << "HFShowerPMT: # of " << str
167  << " bins " << nval << " < 2 ==> illegal";
168  throw cms::Exception("Unknown", "HFShowerPMT")
169  << "nval < 2 for array " << str << "\n";
170  }
171 
172  return fvec;
173  } else {
174  edm::LogError("HFShower") << "HFShowerPMT: cannot get array " << str;
175  throw cms::Exception("Unknown", "HFShowerPMT")
176  << "cannot get array " << str << "\n";
177  }
178 }
#define LogDebug(id)
const double beta
T getParameter(std::string const &) const
double getRadius()
Definition: HFShowerPMT.cc:134
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:140
const double GeV
Definition: MathUtil.h:16
std::vector< double > rTable
Definition: HFShowerPMT.h:40
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:41
#define nullptr
type of data representation of DDCompactView
Definition: DDCompactView.h:90
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
virtual ~HFShowerPMT()
Definition: HFShowerPMT.cc:73
int computeNPEinPMT(G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
Definition: HFCherenkov.cc:229
void initRun(G4ParticleTable *, HcalDDDSimConstants *)
Definition: HFShowerPMT.cc:77
const double MeV
double pePerGeV
Definition: HFShowerPMT.h:38
std::vector< int > pmtR2
Definition: HFShowerPMT.h:42
const std::vector< double > & getRTableHF() const
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::...
Definition: DDsvalues.h:20
Definition: value.py:1
ii
Definition: cuy.py:588
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: HFShowerPMT.cc:152
DDsvalues_type mergedSpecifics() const
double getHits(G4Step *aStep)
Definition: HFShowerPMT.cc:88
bool firstChild()
set the current node to the first child ...
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:37
HFShowerPMT(std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
Definition: HFShowerPMT.cc:22
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:42
std::vector< int > pmtR1
Definition: HFShowerPMT.h:41