CMS 3D CMS Logo

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