CMS 3D CMS Logo

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

#include <HFShowerFibreBundle.h>

Public Member Functions

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

Private Member Functions

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

Private Attributes

HFCherenkovcherenkov1
 
HFCherenkovcherenkov2
 
double facCone
 
double facTube
 
int indexF
 
int indexR
 
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 HFShowerFibreBundle.h.

Constructor & Destructor Documentation

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

Definition at line 23 of file HFShowerFibreBundle.cc.

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

25  {
26 
27  edm::ParameterSet m_HF1 = p.getParameter<edm::ParameterSet>("HFShowerStraightBundle");
28  facTube = m_HF1.getParameter<double>("FactorBundle");
29  cherenkov1 = new HFCherenkov(m_HF1);
30  edm::ParameterSet m_HF2 = p.getParameter<edm::ParameterSet>("HFShowerConicalBundle");
31  facCone = m_HF2.getParameter<double>("FactorBundle");
32  cherenkov2 = new HFCherenkov(m_HF2);
33  edm::LogInfo("HFShower") << "HFShowerFibreBundle intialized with factors: "
34  << facTube << " for the straight portion and "
35  << facCone << " for the curved portion";
36 
37  //Special Geometry parameters
38  std::string attribute = "Volume";
39  std::string value = "HFPMT";
40  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
41  DDFilteredView fv1(cpv,filter1);
42  if (fv1.firstChild()) {
43  DDsvalues_type sv1(fv1.mergedSpecifics());
44  std::vector<double> neta;
45  neta = getDDDArray("indexPMTR",sv1);
46  for (unsigned int ii=0; ii<neta.size(); ii++) {
47  int index = static_cast<int>(neta[ii]);
48  int ir=-1, ifib=-1;
49  if (index >= 0) {
50  ir = index/10; ifib = index%10;
51  }
52  pmtR1.push_back(ir);
53  pmtFib1.push_back(ifib);
54  }
55  neta = getDDDArray("indexPMTL",sv1);
56  for (unsigned int ii=0; ii<neta.size(); ii++) {
57  int index = static_cast<int>(neta[ii]);
58  int ir=-1, ifib=-1;
59  if (index >= 0) {
60  ir = index/10; ifib = index%10;
61  }
62  pmtR2.push_back(ir);
63  pmtFib2.push_back(ifib);
64  }
65  edm::LogInfo("HFShower") << "HFShowerFibreBundle: gets the Index matches "
66  << "for " << neta.size() << " PMTs";
67  for (unsigned int ii=0; ii<neta.size(); ii++)
68  edm::LogInfo("HFShower") << "HFShowerFibreBundle: rIndexR[" << ii
69  << "] = " << pmtR1[ii] << " fibreR[" << ii
70  << "] = " << pmtFib1[ii] << " rIndexL[" << ii
71  << "] = " << pmtR2[ii] << " fibreL[" << ii
72  << "] = " << pmtFib2[ii];
73  } else {
74  edm::LogWarning("HFShower") << "HFShowerFibreBundle: cannot get filtered "
75  << " view for " << attribute << " matching "
76  << value;
77  }
78 
79 }
T getParameter(std::string const &) const
std::vector< int > pmtFib1
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
ii
Definition: cuy.py:589
std::vector< int > pmtR2
std::vector< int > pmtFib2
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
HFShowerFibreBundle::~HFShowerFibreBundle ( )
virtual

Definition at line 81 of file HFShowerFibreBundle.cc.

References cherenkov1, and cherenkov2.

81  {
82  delete cherenkov1;
83  delete cherenkov2;
84 }

Member Function Documentation

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

Definition at line 169 of file HFShowerFibreBundle.cc.

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

Referenced by HFShowerFibreBundle().

170  {
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)
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 HFShowerFibreBundle::getHits ( const G4Step *  aStep,
bool  type 
)

Definition at line 99 of file HFShowerFibreBundle.cc.

References pfBoostedDoubleSVAK8TagInfos_cfi::beta, cherenkov1, cherenkov2, HFCherenkov::computeNPEinPMT(), facCone, facTube, indexF, indexR, LogDebug, MeV, nano_cff::photons, pmtFib1, pmtFib2, pmtR1, and pmtR2.

99  {
100 
101  indexR = indexF = -1;
102 
103  const 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  const G4Track *aTrack = aStep->GetTrack();
125  const 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 }
#define LogDebug(id)
type
Definition: HCALResponse.h:21
std::vector< int > pmtFib1
const double MeV
std::vector< int > pmtR1
std::vector< int > pmtR2
std::vector< int > pmtFib2
int computeNPEinPMT(const G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
Definition: HFCherenkov.cc:229
double HFShowerFibreBundle::getRadius ( )

Definition at line 151 of file HFShowerFibreBundle.cc.

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

151  {
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 }
#define LogDebug(id)
std::vector< double > rTable
void HFShowerFibreBundle::initRun ( const HcalDDDSimConstants hcons)

Definition at line 86 of file HFShowerFibreBundle.cc.

References HcalDDDSimConstants::getRTableHF(), and rTable.

86  {
87 
88  // Special Geometry parameters
89  rTable = hcons->getRTableHF();
90  std::stringstream sss;
91  for (unsigned int ig=0; ig<rTable.size(); ig++) {
92  if(ig/10*10 == ig) { sss << "\n"; }
93  sss << " " << rTable[ig]/cm;
94  }
95  edm::LogInfo("HFShower") << "HFShowerFibreBundle: " << rTable.size()
96  << " rTable(cm):" << sss.str();
97 }
std::vector< double > rTable
const std::vector< double > & getRTableHF() const

Member Data Documentation

HFCherenkov* HFShowerFibreBundle::cherenkov1
private

Definition at line 37 of file HFShowerFibreBundle.h.

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

HFCherenkov * HFShowerFibreBundle::cherenkov2
private

Definition at line 37 of file HFShowerFibreBundle.h.

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

double HFShowerFibreBundle::facCone
private

Definition at line 38 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

double HFShowerFibreBundle::facTube
private

Definition at line 38 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

int HFShowerFibreBundle::indexF
private

Definition at line 39 of file HFShowerFibreBundle.h.

Referenced by getHits(), and getRadius().

int HFShowerFibreBundle::indexR
private

Definition at line 39 of file HFShowerFibreBundle.h.

Referenced by getHits(), and getRadius().

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

Definition at line 41 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 42 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 41 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 42 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 40 of file HFShowerFibreBundle.h.

Referenced by getRadius(), and initRun().