test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  DDSpecificsFilter filter1;
40  DDValue ddv1(attribute,value,0);
41  filter1.setCriteria(ddv1,DDCompOp::equals);
42  DDFilteredView fv1(cpv);
43  fv1.addFilter(filter1);
44  if (fv1.firstChild()) {
46  std::vector<double> neta;
47  neta = getDDDArray("indexPMTR",sv1);
48  for (unsigned int ii=0; ii<neta.size(); ii++) {
49  int index = static_cast<int>(neta[ii]);
50  int ir=-1, ifib=-1;
51  if (index >= 0) {
52  ir = index/10; ifib = index%10;
53  }
54  pmtR1.push_back(ir);
55  pmtFib1.push_back(ifib);
56  }
57  neta = getDDDArray("indexPMTL",sv1);
58  for (unsigned int ii=0; ii<neta.size(); ii++) {
59  int index = static_cast<int>(neta[ii]);
60  int ir=-1, ifib=-1;
61  if (index >= 0) {
62  ir = index/10; ifib = index%10;
63  }
64  pmtR2.push_back(ir);
65  pmtFib2.push_back(ifib);
66  }
67  edm::LogInfo("HFShower") << "HFShowerFibreBundle: gets the Index matches "
68  << "for " << neta.size() << " PMTs";
69  for (unsigned int ii=0; ii<neta.size(); ii++)
70  edm::LogInfo("HFShower") << "HFShowerFibreBundle: rIndexR[" << ii
71  << "] = " << pmtR1[ii] << " fibreR[" << ii
72  << "] = " << pmtFib1[ii] << " rIndexL[" << ii
73  << "] = " << pmtR2[ii] << " fibreL[" << ii
74  << "] = " << pmtFib2[ii];
75  } else {
76  edm::LogWarning("HFShower") << "HFShowerFibreBundle: cannot get filtered "
77  << " view for " << attribute << " matching "
78  << value;
79  }
80 
81 }
82 
84  delete cherenkov1;
85  delete cherenkov2;
86 }
87 
88 void HFShowerFibreBundle::initRun(G4ParticleTable *, HcalDDDSimConstants* hcons) {
89 
90  // Special Geometry parameters
91  rTable = hcons->getRTableHF();
92  edm::LogInfo("HFShower") << "HFShowerFibreBundle: " << rTable.size()
93  << " rTable (cm)";
94  for (unsigned int ig=0; ig<rTable.size(); ig++)
95  edm::LogInfo("HFShower") << "HFShowerFibreBundle: rTable[" << ig << "] = "
96  << rTable[ig]/cm << " cm";
97 }
98 
99 double HFShowerFibreBundle::getHits(G4Step * aStep, bool type) {
100 
101  indexR = indexF = -1;
102 
103  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
104  const G4VTouchable* touch = preStepPoint->GetTouchable();
105  int boxNo = touch->GetReplicaNumber(1);
106  int pmtNo = touch->GetReplicaNumber(0);
107  if (boxNo <= 1) {
108  indexR = pmtR1[pmtNo-1];
109  indexF = pmtFib1[pmtNo-1];
110  } else {
111  indexR = pmtR2[pmtNo-1];
112  indexF = pmtFib2[pmtNo-1];
113  }
114 
115 #ifdef DebugLog
116  double edep = aStep->GetTotalEnergyDeposit();
117  LogDebug("HFShower") << "HFShowerFibreBundle: Box " << boxNo << " PMT "
118  << pmtNo << " Mapped Indices " << indexR << ", "
119  << indexF << " Edeposit " << edep/MeV << " MeV";
120 #endif
121 
122  double photons = 0;
123  if (indexR >= 0 && indexF > 0) {
124  G4Track *aTrack = aStep->GetTrack();
125  G4ParticleDefinition *particleDef = aTrack->GetDefinition();
126  double stepl = aStep->GetStepLength();
127  double beta = preStepPoint->GetBeta();
128  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
129  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
130  GetTopTransform().TransformAxis(pDir);
131  if (type) {
132  photons = facCone*cherenkov2->computeNPEinPMT(particleDef, beta,
133  localMom.x(), localMom.y(),
134  localMom.z(), stepl);
135  } else {
136  photons = facTube*cherenkov1->computeNPEinPMT(particleDef, beta,
137  localMom.x(), localMom.y(),
138  localMom.z(), stepl);
139  }
140 #ifdef DebugLog
141  LogDebug("HFShower") << "HFShowerFibreBundle::getHits: for particle "
142  << particleDef->GetParticleName() << " Step " << stepl
143  << " Beta " << beta << " Direction " << pDir
144  << " Local " << localMom << " p.e. " << photons;
145 #endif
146 
147  }
148  return photons;
149 }
150 
152 
153  double r = 0.;
154  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
155  r = 0.5*(rTable[indexR]+rTable[indexR+1]);
156 #ifdef DebugLog
157  else
158  LogDebug("HFShower") << "HFShowerFibreBundle::getRadius: R " << indexR
159  << " F " << indexF;
160 #endif
161  if (indexF == 2) r =-r;
162 #ifdef DebugLog
163  LogDebug("HFShower") << "HFShowerFibreBundle: Radius (" << indexR << "/"
164  << indexF << ") " << r;
165 #endif
166  return r;
167 }
168 
169 std::vector<double> HFShowerFibreBundle::getDDDArray(const std::string & str,
170  const DDsvalues_type& sv){
171 
172 #ifdef DebugLog
173  LogDebug("HFShower") << "HFShowerFibreBundle:getDDDArray called for " << str;
174 #endif
175  DDValue value(str);
176  if (DDfetch(&sv,value)) {
177 #ifdef DebugLog
178  LogDebug("HFShower") << value;
179 #endif
180  const std::vector<double> & fvec = value.doubles();
181  int nval = fvec.size();
182  if (nval < 2) {
183  edm::LogError("HFShower") << "HFShowerFibreBundle: # of " << str
184  << " bins " << nval << " < 2 ==> illegal";
185  throw cms::Exception("Unknown", "HFShowerFibreBundle")
186  << "nval < 2 for array " << str << "\n";
187  }
188 
189  return fvec;
190  } else {
191  edm::LogError("HFShower") <<"HFShowerFibreBundle: cannot get array " <<str;
192  throw cms::Exception("Unknown", "HFShowerFibreBundle")
193  << "cannot get array " << str << "\n";
194  }
195 }
#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:139
std::vector< int > pmtFib1
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
std::vector< double > rTable
double getHits(G4Step *aStep, bool type)
int ii
Definition: cuy.py:588
type of data representation of DDCompactView
Definition: DDCompactView.h:77
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:80
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:19
std::vector< int > pmtR1
void initRun(G4ParticleTable *, HcalDDDSimConstants *)
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
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:245
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:32