CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HFShowerPMT Class Reference

#include <HFShowerPMT.h>

Public Member Functions

double getHits (G4Step *aStep)
 
double getRadius ()
 
 HFShowerPMT (std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
 
void initRun (G4ParticleTable *, 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 ( std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 22 of file HFShowerPMT.cc.

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

23  : cherenkov(0) {
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  DDSpecificsFilter filter1;
32  DDValue ddv1(attribute,value,0);
33  filter1.setCriteria(ddv1,DDCompOp::equals);
34  DDFilteredView fv1(cpv);
35  fv1.addFilter(filter1);
36  if (fv1.firstChild()) {
37  DDsvalues_type sv1(fv1.mergedSpecifics());
38  std::vector<double> neta;
39  neta = getDDDArray("indexPMTR",sv1);
40  for (unsigned int ii=0; ii<neta.size(); ii++) {
41  int index = static_cast<int>(neta[ii]);
42  int ir=-1, ifib=-1;
43  if (index >= 0) {
44  ir = index/10; ifib = index%10;
45  }
46  pmtR1.push_back(ir);
47  pmtFib1.push_back(ifib);
48  }
49  neta = getDDDArray("indexPMTL",sv1);
50  for (unsigned int ii=0; ii<neta.size(); ii++) {
51  int index = static_cast<int>(neta[ii]);
52  int ir=-1, ifib=-1;
53  if (index >= 0) {
54  ir = index/10; ifib = index%10;
55  }
56  pmtR2.push_back(ir);
57  pmtFib2.push_back(ifib);
58  }
59  edm::LogInfo("HFShower") << "HFShowerPMT: gets the Index matches for "
60  << neta.size() << " PMTs";
61  for (unsigned int ii=0; ii<neta.size(); ii++)
62  edm::LogInfo("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = "
63  << pmtR1[ii] << " fibreR[" << ii << "] = "
64  << pmtFib1[ii] << " rIndexL[" << ii << "] = "
65  << pmtR2[ii] << " fibreL[" << ii << "] = "
66  << pmtFib2[ii];
67  } else {
68  edm::LogWarning("HFShower") << "HFShowerPMT: cannot get filtered "
69  << " view for " << attribute << " matching "
70  << value;
71  }
72 
73  cherenkov = new HFCherenkov(m_HF);
74 }
T getParameter(std::string const &) const
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:41
int ii
Definition: cuy.py:588
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:19
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: HFShowerPMT.cc:155
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:37
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:42
std::vector< int > pmtR1
Definition: HFShowerPMT.h:41
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:245
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:32
HFShowerPMT::~HFShowerPMT ( )
virtual

Definition at line 76 of file HFShowerPMT.cc.

References cherenkov.

76  {
77  if (cherenkov) delete cherenkov;
78 }
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 155 of file HFShowerPMT.cc.

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

Referenced by HFShowerPMT().

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

Definition at line 91 of file HFShowerPMT.cc.

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

Referenced by HCalSD::getHitPMT().

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

Definition at line 137 of file HFShowerPMT.cc.

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

Referenced by HCalSD::getHitPMT().

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

Definition at line 80 of file HFShowerPMT.cc.

References HcalDDDSimConstants::getRTableHF(), and rTable.

Referenced by HCalSD::initRun().

80  {
81 
82  // Special Geometry parameters
83  rTable = hcons->getRTableHF();
84  edm::LogInfo("HFShower") << "HFShowerPMT: " << rTable.size()
85  << " rTable (cm)";
86  for (unsigned int ig=0; ig<rTable.size(); ig++)
87  edm::LogInfo("HFShower") << "HFShowerPMT: rTable[" << ig << "] = "
88  << rTable[ig]/cm << " cm";
89 }
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().