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 
20 //#define DebugLog
21 
23  const DDCompactView & cpv,
24  edm::ParameterSet const & p) {
25 
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::LogInfo("HFShower") << "HFShowerFibreBundle intialized with factors: "
33  << facTube << " for the straight portion and "
34  << facCone << " for the curved portion";
35 
36  //Special Geometry parameters
37  std::string attribute = "Volume";
38  std::string value = "HFPMT";
39  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
40  DDFilteredView fv1(cpv,filter1);
41  if (fv1.firstChild()) {
43  std::vector<double> neta;
44  neta = getDDDArray("indexPMTR",sv1);
45  for (unsigned int ii=0; ii<neta.size(); ii++) {
46  int index = static_cast<int>(neta[ii]);
47  int ir=-1, ifib=-1;
48  if (index >= 0) {
49  ir = index/10; 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; ifib = index%10;
60  }
61  pmtR2.push_back(ir);
62  pmtFib2.push_back(ifib);
63  }
64  edm::LogInfo("HFShower") << "HFShowerFibreBundle: gets the Index matches "
65  << "for " << neta.size() << " PMTs";
66  for (unsigned int ii=0; ii<neta.size(); ii++)
67  edm::LogInfo("HFShower") << "HFShowerFibreBundle: rIndexR[" << ii
68  << "] = " << pmtR1[ii] << " fibreR[" << ii
69  << "] = " << pmtFib1[ii] << " rIndexL[" << ii
70  << "] = " << pmtR2[ii] << " fibreL[" << ii
71  << "] = " << pmtFib2[ii];
72  } else {
73  edm::LogWarning("HFShower") << "HFShowerFibreBundle: cannot get filtered "
74  << " view for " << attribute << " matching "
75  << value;
76  }
77 
78 }
79 
81  delete cherenkov1;
82  delete cherenkov2;
83 }
84 
85 void HFShowerFibreBundle::initRun(G4ParticleTable *, HcalDDDSimConstants* hcons) {
86 
87  // Special Geometry parameters
88  rTable = hcons->getRTableHF();
89  edm::LogInfo("HFShower") << "HFShowerFibreBundle: " << rTable.size()
90  << " rTable (cm)";
91  for (unsigned int ig=0; ig<rTable.size(); ig++)
92  edm::LogInfo("HFShower") << "HFShowerFibreBundle: rTable[" << ig << "] = "
93  << rTable[ig]/cm << " cm";
94 }
95 
96 double HFShowerFibreBundle::getHits(G4Step * aStep, bool type) {
97 
98  indexR = indexF = -1;
99 
100  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
101  const G4VTouchable* touch = preStepPoint->GetTouchable();
102  int boxNo = touch->GetReplicaNumber(1);
103  int pmtNo = touch->GetReplicaNumber(0);
104  if (boxNo <= 1) {
105  indexR = pmtR1[pmtNo-1];
106  indexF = pmtFib1[pmtNo-1];
107  } else {
108  indexR = pmtR2[pmtNo-1];
109  indexF = pmtFib2[pmtNo-1];
110  }
111 
112 #ifdef DebugLog
113  double edep = aStep->GetTotalEnergyDeposit();
114  LogDebug("HFShower") << "HFShowerFibreBundle: Box " << boxNo << " PMT "
115  << pmtNo << " Mapped Indices " << indexR << ", "
116  << indexF << " Edeposit " << edep/MeV << " MeV";
117 #endif
118 
119  double photons = 0;
120  if (indexR >= 0 && indexF > 0) {
121  G4Track *aTrack = aStep->GetTrack();
122  G4ParticleDefinition *particleDef = aTrack->GetDefinition();
123  double stepl = aStep->GetStepLength();
124  double beta = preStepPoint->GetBeta();
125  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
126  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
127  GetTopTransform().TransformAxis(pDir);
128  if (type) {
129  photons = facCone*cherenkov2->computeNPEinPMT(particleDef, beta,
130  localMom.x(), localMom.y(),
131  localMom.z(), stepl);
132  } else {
133  photons = facTube*cherenkov1->computeNPEinPMT(particleDef, beta,
134  localMom.x(), localMom.y(),
135  localMom.z(), stepl);
136  }
137 #ifdef DebugLog
138  LogDebug("HFShower") << "HFShowerFibreBundle::getHits: for particle "
139  << particleDef->GetParticleName() << " Step " << stepl
140  << " Beta " << beta << " Direction " << pDir
141  << " Local " << localMom << " p.e. " << photons;
142 #endif
143 
144  }
145  return photons;
146 }
147 
149 
150  double r = 0.;
151  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
152  r = 0.5*(rTable[indexR]+rTable[indexR+1]);
153 #ifdef DebugLog
154  else
155  LogDebug("HFShower") << "HFShowerFibreBundle::getRadius: R " << indexR
156  << " F " << indexF;
157 #endif
158  if (indexF == 2) r =-r;
159 #ifdef DebugLog
160  LogDebug("HFShower") << "HFShowerFibreBundle: Radius (" << indexR << "/"
161  << indexF << ") " << r;
162 #endif
163  return r;
164 }
165 
166 std::vector<double> HFShowerFibreBundle::getDDDArray(const std::string & str,
167  const DDsvalues_type& sv){
168 
169 #ifdef DebugLog
170  LogDebug("HFShower") << "HFShowerFibreBundle:getDDDArray called for " << str;
171 #endif
172  DDValue value(str);
173  if (DDfetch(&sv,value)) {
174 #ifdef DebugLog
175  LogDebug("HFShower") << value;
176 #endif
177  const std::vector<double> & fvec = value.doubles();
178  int nval = fvec.size();
179  if (nval < 2) {
180  edm::LogError("HFShower") << "HFShowerFibreBundle: # of " << str
181  << " bins " << nval << " < 2 ==> illegal";
182  throw cms::Exception("Unknown", "HFShowerFibreBundle")
183  << "nval < 2 for array " << str << "\n";
184  }
185 
186  return fvec;
187  } else {
188  edm::LogError("HFShower") <<"HFShowerFibreBundle: cannot get array " <<str;
189  throw cms::Exception("Unknown", "HFShowerFibreBundle")
190  << "cannot get array " << str << "\n";
191  }
192 }
#define LogDebug(id)
const double beta
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
double getHits(G4Step *aStep, bool type)
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
int computeNPEinPMT(G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
Definition: HFCherenkov.cc:229
const double MeV
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
std::vector< int > pmtR1
Definition: value.py:1
void initRun(G4ParticleTable *, HcalDDDSimConstants *)
ii
Definition: cuy.py:588
HFShowerFibreBundle(std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
DDsvalues_type mergedSpecifics() const
bool firstChild()
set the current node to the first child ...
std::vector< int > pmtR2
std::vector< int > pmtFib2
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)