CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HFShowerPMT Class Reference

#include <HFShowerPMT.h>

Public Member Functions

double getHits (const G4Step *aStep)
 
double getRadius ()
 
 HFShowerPMT (const std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
 
void initRun (const HcalDDDSimConstants *)
 
virtual ~HFShowerPMT ()
 

Private Member Functions

std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &)
 

Private Attributes

HFCherenkovcherenkov
 
int indexF
 
int indexR
 
double pePerGeV
 
std::vector< int > pmtFib1
 
std::vector< int > pmtFib2
 
std::vector< int > pmtR1
 
std::vector< int > pmtR2
 
std::vector< double > rTable
 

Detailed Description

Definition at line 20 of file HFShowerPMT.h.

Constructor & Destructor Documentation

HFShowerPMT::HFShowerPMT ( const std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 23 of file HFShowerPMT.cc.

References cherenkov, DDFilteredView::firstChild(), getDDDArray(), edm::ParameterSet::getParameter(), cuy::ii, DDFilteredView::mergedSpecifics(), pePerGeV, pmtFib1, pmtFib2, pmtR1, pmtR2, AlCaHLTBitMon_QueryRunRegistry::string, and relativeConstraints::value.

24  : cherenkov(nullptr) {
25 
26  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShowerPMT");
27  pePerGeV = m_HF.getParameter<double>("PEPerGeVPMT");
28 
29  //Special Geometry parameters
30  std::string attribute = "Volume";
31  std::string value = "HFPMT";
32  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
33  DDFilteredView fv1(cpv,filter1);
34  if (fv1.firstChild()) {
35  DDsvalues_type sv1(fv1.mergedSpecifics());
36  std::vector<double> neta;
37  neta = getDDDArray("indexPMTR",sv1);
38  for (unsigned int ii=0; ii<neta.size(); ii++) {
39  int index = static_cast<int>(neta[ii]);
40  int ir=-1, ifib=-1;
41  if (index >= 0) {
42  ir = index/10; ifib = index%10;
43  }
44  pmtR1.push_back(ir);
45  pmtFib1.push_back(ifib);
46  }
47  neta = getDDDArray("indexPMTL",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  pmtR2.push_back(ir);
55  pmtFib2.push_back(ifib);
56  }
57 #ifdef DebugLog
58  edm::LogInfo("HFShower") << "HFShowerPMT: gets the Index matches for "
59  << neta.size() << " PMTs";
60  for (unsigned int ii=0; ii<neta.size(); ii++) {
61  edm::LogInfo("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = "
62  << pmtR1[ii] << " fibreR[" << ii << "] = "
63  << pmtFib1[ii] << " rIndexL[" << ii << "] = "
64  << pmtR2[ii] << " fibreL[" << ii << "] = "
65  << pmtFib2[ii];
66  }
67 #endif
68  } else {
69  edm::LogWarning("HFShower") << "HFShowerPMT: cannot get filtered "
70  << " view for " << attribute << " matching "
71  << value;
72  }
73 
74  cherenkov = new HFCherenkov(m_HF);
75 }
T getParameter(std::string const &) const
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:41
double pePerGeV
Definition: HFShowerPMT.h:38
std::vector< int > pmtR2
Definition: HFShowerPMT.h:42
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:589
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: HFShowerPMT.cc:158
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:37
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:42
std::vector< int > pmtR1
Definition: HFShowerPMT.h:41
HFShowerPMT::~HFShowerPMT ( )
virtual

Definition at line 77 of file HFShowerPMT.cc.

References cherenkov.

77  {
78  if (cherenkov) delete cherenkov;
79 }
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:37

Member Function Documentation

std::vector< double > HFShowerPMT::getDDDArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 158 of file HFShowerPMT.cc.

References DDfetch(), DDValue::doubles(), Exception, LogDebug, str, and relativeConstraints::value.

Referenced by HFShowerPMT().

159  {
160 
161 #ifdef DebugLog
162  LogDebug("HFShower") << "HFShowerPMT:getDDDArray called for " << str;
163 #endif
164  DDValue value(str);
165  if (DDfetch(&sv,value)) {
166 #ifdef DebugLog
167  LogDebug("HFShower") << value;
168 #endif
169  const std::vector<double> & fvec = value.doubles();
170  int nval = fvec.size();
171  if (nval < 2) {
172  edm::LogError("HFShower") << "HFShowerPMT: # of " << str
173  << " bins " << nval << " < 2 ==> illegal";
174  throw cms::Exception("Unknown", "HFShowerPMT")
175  << "nval < 2 for array " << str << "\n";
176  }
177 
178  return fvec;
179  } else {
180  edm::LogError("HFShower") << "HFShowerPMT: cannot get array " << str;
181  throw cms::Exception("Unknown", "HFShowerPMT")
182  << "cannot get array " << str << "\n";
183  }
184 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
Definition: value.py:1
#define str(s)
double HFShowerPMT::getHits ( const G4Step *  aStep)

Definition at line 94 of file HFShowerPMT.cc.

References pfBoostedDoubleSVAK8TagInfos_cfi::beta, cherenkov, HFCherenkov::computeNPEinPMT(), GeV, indexF, indexR, LogDebug, MeV, pePerGeV, nano_cff::photons, pmtFib1, pmtFib2, pmtR1, and pmtR2.

94  {
95 
96  indexR = indexF = -1;
97 
98  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
99  const G4VTouchable* touch = preStepPoint->GetTouchable();
100  int boxNo = touch->GetReplicaNumber(2);
101  int pmtNo = touch->GetReplicaNumber(1);
102  if (boxNo <= 1) {
103  indexR = pmtR1[pmtNo-1];
104  indexF = pmtFib1[pmtNo-1];
105  } else {
106  indexR = pmtR2[pmtNo-1];
107  indexF = pmtFib2[pmtNo-1];
108  }
109 
110 #ifdef DebugLog
111  double edep = aStep->GetTotalEnergyDeposit();
112  LogDebug("HFShower") << "HFShowerPMT: Box " << boxNo << " PMT "
113  << pmtNo << " Mapped Indices " << indexR << ", "
114  << indexF << " Edeposit " << edep/MeV << " MeV; PE "
115  << edep*pePerGeV/GeV;
116 #endif
117 
118  double photons = 0;
119  if (indexR >= 0 && indexF > 0) {
120  G4Track *aTrack = aStep->GetTrack();
121  G4ParticleDefinition *particleDef = aTrack->GetDefinition();
122  double stepl = aStep->GetStepLength();
123  double beta = preStepPoint->GetBeta();
124  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
125  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
126  GetTopTransform().TransformAxis(pDir);
127  photons = cherenkov->computeNPEinPMT(particleDef, beta, localMom.x(),
128  localMom.y(), localMom.z(), stepl);
129 #ifdef DebugLog
130  LogDebug("HFShower") << "HFShowerPMT::getHits: for particle "
131  << particleDef->GetParticleName() << " Step " << stepl
132  << " Beta " << beta << " Direction " << pDir
133  << " Local " << localMom << " p.e. " << photons;
134 #endif
135 
136  }
137  return photons;
138 }
#define LogDebug(id)
const double GeV
Definition: MathUtil.h:16
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:41
const double MeV
double pePerGeV
Definition: HFShowerPMT.h:38
std::vector< int > pmtR2
Definition: HFShowerPMT.h:42
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:37
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:42
std::vector< int > pmtR1
Definition: HFShowerPMT.h:41
int computeNPEinPMT(const G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
Definition: HFCherenkov.cc:229
double HFShowerPMT::getRadius ( )

Definition at line 140 of file HFShowerPMT.cc.

References indexF, indexR, LogDebug, alignCSCRings::r, and rTable.

140  {
141 
142  double r = 0.;
143  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
144  r = 0.5*(rTable[indexR]+rTable[indexR+1]);
145 #ifdef DebugLog
146  else
147  LogDebug("HFShower") << "HFShowerPMT::getRadius: R " << indexR
148  << " F " << indexF;
149 #endif
150  if (indexF == 2) r =-r;
151 #ifdef DebugLog
152  LogDebug("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF
153  << ") " << r;
154 #endif
155  return r;
156 }
#define LogDebug(id)
std::vector< double > rTable
Definition: HFShowerPMT.h:40
void HFShowerPMT::initRun ( const HcalDDDSimConstants hcons)

Definition at line 81 of file HFShowerPMT.cc.

References HcalDDDSimConstants::getRTableHF(), and rTable.

81  {
82 
83  // Special Geometry parameters
84  rTable = hcons->getRTableHF();
85  std::stringstream sss;
86  for (unsigned int ig=0; ig<rTable.size(); ++ig) {
87  if(ig/10*10 == ig) { sss << "\n"; }
88  sss << " " << rTable[ig]/cm;
89  }
90  edm::LogInfo("HFShowerPMT") << "HFShowerPMT: " << rTable.size()
91  << " rTable(cm):" << sss.str();
92 }
std::vector< double > rTable
Definition: HFShowerPMT.h:40
const std::vector< double > & getRTableHF() const

Member Data Documentation

HFCherenkov* HFShowerPMT::cherenkov
private

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), HFShowerPMT(), and ~HFShowerPMT().

int HFShowerPMT::indexF
private

Definition at line 39 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

int HFShowerPMT::indexR
private

Definition at line 39 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

double HFShowerPMT::pePerGeV
private

Definition at line 38 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtFib1
private

Definition at line 41 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtFib2
private

Definition at line 42 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtR1
private

Definition at line 41 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtR2
private

Definition at line 42 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<double> HFShowerPMT::rTable
private

Definition at line 40 of file HFShowerPMT.h.

Referenced by getRadius(), and initRun().