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)
 
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 19 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, DDSpecificsFilter::equals, edm::hlt::Exception, DDFilteredView::firstChild(), getDDDArray(), edm::ParameterSet::getParameter(), cuy::ii, cmsHarvester::index, DDFilteredView::mergedSpecifics(), mergeVDriftHistosByStation::name, pePerGeV, pmtFib1, pmtFib2, pmtR1, pmtR2, rTable, DDSpecificsFilter::setCriteria(), 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  G4String attribute = "ReadOutName";
29  G4String value = name;
30  DDSpecificsFilter filter0;
31  DDValue ddv0(attribute,value,0);
33  DDFilteredView fv0(cpv);
34  fv0.addFilter(filter0);
35  if (fv0.firstChild()) {
36  DDsvalues_type sv0(fv0.mergedSpecifics());
37 
38  //Special Geometry parameters
39  rTable = getDDDArray("rTable",sv0);
40  edm::LogInfo("HFShower") << "HFShowerPMT: " << rTable.size()
41  << " rTable (cm)";
42  for (unsigned int ig=0; ig<rTable.size(); ig++)
43  edm::LogInfo("HFShower") << "HFShowerPMT: rTable[" << ig << "] = "
44  << rTable[ig]/cm << " cm";
45  } else {
46  edm::LogError("HFShower") << "HFShowerPMT: cannot get filtered "
47  << " view for " << attribute << " matching "
48  << value;
49  throw cms::Exception("Unknown", "HFShowerPMT")
50  << "cannot match " << attribute << " to " << name <<"\n";
51  }
52 
53  attribute = "Volume";
54  value = "HFPMT";
55  DDSpecificsFilter filter1;
56  DDValue ddv1(attribute,value,0);
58  DDFilteredView fv1(cpv);
59  fv1.addFilter(filter1);
60  if (fv1.firstChild()) {
61  DDsvalues_type sv1(fv1.mergedSpecifics());
62  std::vector<double> neta;
63  neta = getDDDArray("indexPMTR",sv1);
64  for (unsigned int ii=0; ii<neta.size(); ii++) {
65  int index = static_cast<int>(neta[ii]);
66  int ir=-1, ifib=-1;
67  if (index >= 0) {
68  ir = index/10; ifib = index%10;
69  }
70  pmtR1.push_back(ir);
71  pmtFib1.push_back(ifib);
72  }
73  neta = getDDDArray("indexPMTL",sv1);
74  for (unsigned int ii=0; ii<neta.size(); ii++) {
75  int index = static_cast<int>(neta[ii]);
76  int ir=-1, ifib=-1;
77  if (index >= 0) {
78  ir = index/10; ifib = index%10;
79  }
80  pmtR2.push_back(ir);
81  pmtFib2.push_back(ifib);
82  }
83  edm::LogInfo("HFShower") << "HFShowerPMT: gets the Index matches for "
84  << neta.size() << " PMTs";
85  for (unsigned int ii=0; ii<neta.size(); ii++)
86  edm::LogInfo("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = "
87  << pmtR1[ii] << " fibreR[" << ii << "] = "
88  << pmtFib1[ii] << " rIndexL[" << ii << "] = "
89  << pmtR2[ii] << " fibreL[" << ii << "] = "
90  << pmtFib2[ii];
91  } else {
92  edm::LogWarning("HFShower") << "HFShowerPMT: cannot get filtered "
93  << " view for " << attribute << " matching "
94  << value;
95  }
96 
97  cherenkov = new HFCherenkov(m_HF);
98 }
T getParameter(std::string const &) const
std::vector< double > rTable
Definition: HFShowerPMT.h:38
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:39
int ii
Definition: cuy.py:588
double pePerGeV
Definition: HFShowerPMT.h:36
std::vector< int > pmtR2
Definition: HFShowerPMT.h:40
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:168
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:35
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:40
std::vector< int > pmtR1
Definition: HFShowerPMT.h:39
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
HFShowerPMT::~HFShowerPMT ( )
virtual

Definition at line 100 of file HFShowerPMT.cc.

References cherenkov.

100  {
101  if (cherenkov) delete cherenkov;
102 }
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:35

Member Function Documentation

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

Definition at line 168 of file HFShowerPMT.cc.

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

Referenced by HFShowerPMT().

169  {
170 
171 #ifdef DebugLog
172  LogDebug("HFShower") << "HFShowerPMT:getDDDArray called for " << str;
173 #endif
174  DDValue value(str);
175  if (DDfetch(&sv,value)) {
176 #ifdef DebugLog
177  LogDebug("HFShower") << value;
178 #endif
179  const std::vector<double> & fvec = value.doubles();
180  int nval = fvec.size();
181  if (nval < 2) {
182  edm::LogError("HFShower") << "HFShowerPMT: # of " << str
183  << " bins " << nval << " < 2 ==> illegal";
184  throw cms::Exception("Unknown", "HFShowerPMT")
185  << "nval < 2 for array " << str << "\n";
186  }
187 
188  return fvec;
189  } else {
190  edm::LogError("HFShower") << "HFShowerPMT: cannot get array " << str;
191  throw cms::Exception("Unknown", "HFShowerPMT")
192  << "cannot get array " << str << "\n";
193  }
194 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
double HFShowerPMT::getHits ( G4Step *  aStep)

Definition at line 104 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().

104  {
105 
106  indexR = indexF = -1;
107 
108  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
109  const G4VTouchable* touch = preStepPoint->GetTouchable();
110  int boxNo = touch->GetReplicaNumber(2);
111  int pmtNo = touch->GetReplicaNumber(1);
112  if (boxNo <= 1) {
113  indexR = pmtR1[pmtNo-1];
114  indexF = pmtFib1[pmtNo-1];
115  } else {
116  indexR = pmtR2[pmtNo-1];
117  indexF = pmtFib2[pmtNo-1];
118  }
119 
120 #ifdef DebugLog
121  double edep = aStep->GetTotalEnergyDeposit();
122  LogDebug("HFShower") << "HFShowerPMT: Box " << boxNo << " PMT "
123  << pmtNo << " Mapped Indices " << indexR << ", "
124  << indexF << " Edeposit " << edep/MeV << " MeV; PE "
125  << edep*pePerGeV/GeV;
126 #endif
127 
128  double photons = 0;
129  if (indexR >= 0 && indexF > 0) {
130  G4Track *aTrack = aStep->GetTrack();
131  G4ParticleDefinition *particleDef = aTrack->GetDefinition();
132  double stepl = aStep->GetStepLength();
133  double beta = preStepPoint->GetBeta();
134  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
135  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
136  GetTopTransform().TransformAxis(pDir);
137  photons = cherenkov->computeNPEinPMT(particleDef, beta, localMom.x(),
138  localMom.y(), localMom.z(), stepl);
139 #ifdef DebugLog
140  LogDebug("HFShower") << "HFShowerPMT::getHits: for particle "
141  << particleDef->GetParticleName() << " Step " << stepl
142  << " Beta " << beta << " Direction " << pDir
143  << " Local " << localMom << " p.e. " << photons;
144 #endif
145 
146  }
147  return photons;
148 }
#define LogDebug(id)
const double beta
const double GeV
Definition: MathUtil.h:16
std::vector< int > pmtFib1
Definition: HFShowerPMT.h:39
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:36
std::vector< int > pmtR2
Definition: HFShowerPMT.h:40
HFCherenkov * cherenkov
Definition: HFShowerPMT.h:35
std::vector< int > pmtFib2
Definition: HFShowerPMT.h:40
std::vector< int > pmtR1
Definition: HFShowerPMT.h:39
double HFShowerPMT::getRadius ( )

Definition at line 150 of file HFShowerPMT.cc.

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

Referenced by HCalSD::getHitPMT().

150  {
151 
152  double r = 0.;
153  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
154  r = 0.5*(rTable[indexR]+rTable[indexR+1]);
155 #ifdef DebugLog
156  else
157  LogDebug("HFShower") << "HFShowerPMT::getRadius: R " << indexR
158  << " F " << indexF;
159 #endif
160  if (indexF == 2) r =-r;
161 #ifdef DebugLog
162  LogDebug("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF
163  << ") " << r;
164 #endif
165  return r;
166 }
#define LogDebug(id)
std::vector< double > rTable
Definition: HFShowerPMT.h:38

Member Data Documentation

HFCherenkov* HFShowerPMT::cherenkov
private

Definition at line 35 of file HFShowerPMT.h.

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

int HFShowerPMT::indexF
private

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

int HFShowerPMT::indexR
private

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

double HFShowerPMT::pePerGeV
private

Definition at line 36 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

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

Definition at line 39 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

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

Definition at line 40 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

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

Definition at line 39 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

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

Definition at line 40 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

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

Definition at line 38 of file HFShowerPMT.h.

Referenced by getRadius(), and HFShowerPMT().