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
HFShowerFibreBundle Class Reference

#include <HFShowerFibreBundle.h>

Public Member Functions

double getHits (G4Step *aStep, bool type)
 
double getRadius ()
 
 HFShowerFibreBundle (std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
 
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 19 of file HFShowerFibreBundle.h.

Constructor & Destructor Documentation

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

Definition at line 22 of file HFShowerFibreBundle.cc.

References DDFilteredView::addFilter(), cherenkov1, cherenkov2, DDSpecificsFilter::equals, edm::hlt::Exception, facCone, facTube, DDFilteredView::firstChild(), getDDDArray(), edm::ParameterSet::getParameter(), getHLTprescales::index, DDFilteredView::mergedSpecifics(), mergeVDriftHistosByStation::name, pmtFib1, pmtFib2, pmtR1, pmtR2, rTable, DDSpecificsFilter::setCriteria(), and relativeConstraints::value.

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

Definition at line 107 of file HFShowerFibreBundle.cc.

References cherenkov1, and cherenkov2.

107  {
108  delete cherenkov1;
109  delete cherenkov2;
110 }

Member Function Documentation

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

Definition at line 182 of file HFShowerFibreBundle.cc.

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

Referenced by HFShowerFibreBundle().

183  {
184 
185 #ifdef DebugLog
186  LogDebug("HFShower") << "HFShowerFibreBundle:getDDDArray called for " << str;
187 #endif
188  DDValue value(str);
189  if (DDfetch(&sv,value)) {
190 #ifdef DebugLog
191  LogDebug("HFShower") << value;
192 #endif
193  const std::vector<double> & fvec = value.doubles();
194  int nval = fvec.size();
195  if (nval < 2) {
196  edm::LogError("HFShower") << "HFShowerFibreBundle: # of " << str
197  << " bins " << nval << " < 2 ==> illegal";
198  throw cms::Exception("Unknown", "HFShowerFibreBundle")
199  << "nval < 2 for array " << str << "\n";
200  }
201 
202  return fvec;
203  } else {
204  edm::LogError("HFShower") <<"HFShowerFibreBundle: cannot get array " <<str;
205  throw cms::Exception("Unknown", "HFShowerFibreBundle")
206  << "cannot get array " << str << "\n";
207  }
208 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
double HFShowerFibreBundle::getHits ( G4Step *  aStep,
bool  type 
)

Definition at line 112 of file HFShowerFibreBundle.cc.

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

Referenced by HCalSD::getHitFibreBundle().

112  {
113 
114  indexR = indexF = -1;
115 
116  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
117  const G4VTouchable* touch = preStepPoint->GetTouchable();
118  int boxNo = touch->GetReplicaNumber(1);
119  int pmtNo = touch->GetReplicaNumber(0);
120  if (boxNo <= 1) {
121  indexR = pmtR1[pmtNo-1];
122  indexF = pmtFib1[pmtNo-1];
123  } else {
124  indexR = pmtR2[pmtNo-1];
125  indexF = pmtFib2[pmtNo-1];
126  }
127 
128 #ifdef DebugLog
129  double edep = aStep->GetTotalEnergyDeposit();
130  LogDebug("HFShower") << "HFShowerFibreBundle: Box " << boxNo << " PMT "
131  << pmtNo << " Mapped Indices " << indexR << ", "
132  << indexF << " Edeposit " << edep/MeV << " MeV";
133 #endif
134 
135  double photons = 0;
136  if (indexR >= 0 && indexF > 0) {
137  G4Track *aTrack = aStep->GetTrack();
138  G4ParticleDefinition *particleDef = aTrack->GetDefinition();
139  double stepl = aStep->GetStepLength();
140  double beta = preStepPoint->GetBeta();
141  G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
142  G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
143  GetTopTransform().TransformAxis(pDir);
144  if (type) {
145  photons = facCone*cherenkov2->computeNPEinPMT(particleDef, beta,
146  localMom.x(), localMom.y(),
147  localMom.z(), stepl);
148  } else {
149  photons = facTube*cherenkov1->computeNPEinPMT(particleDef, beta,
150  localMom.x(), localMom.y(),
151  localMom.z(), stepl);
152  }
153 #ifdef DebugLog
154  LogDebug("HFShower") << "HFShowerFibreBundle::getHits: for particle "
155  << particleDef->GetParticleName() << " Step " << stepl
156  << " Beta " << beta << " Direction " << pDir
157  << " Local " << localMom << " p.e. " << photons;
158 #endif
159 
160  }
161  return photons;
162 }
#define LogDebug(id)
const double beta
type
Definition: HCALResponse.h:22
std::vector< int > pmtFib1
int computeNPEinPMT(G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
Definition: HFCherenkov.cc:183
std::vector< int > pmtR1
std::vector< int > pmtR2
std::vector< int > pmtFib2
double HFShowerFibreBundle::getRadius ( )

Definition at line 164 of file HFShowerFibreBundle.cc.

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

Referenced by HCalSD::getHitFibreBundle().

164  {
165 
166  double r = 0.;
167  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
168  r = 0.5*(rTable[indexR]+rTable[indexR+1]);
169 #ifdef DebugLog
170  else
171  LogDebug("HFShower") << "HFShowerFibreBundle::getRadius: R " << indexR
172  << " F " << indexF;
173 #endif
174  if (indexF == 2) r =-r;
175 #ifdef DebugLog
176  LogDebug("HFShower") << "HFShowerFibreBundle: Radius (" << indexR << "/"
177  << indexF << ") " << r;
178 #endif
179  return r;
180 }
#define LogDebug(id)
std::vector< double > rTable

Member Data Documentation

HFCherenkov* HFShowerFibreBundle::cherenkov1
private

Definition at line 35 of file HFShowerFibreBundle.h.

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

HFCherenkov * HFShowerFibreBundle::cherenkov2
private

Definition at line 35 of file HFShowerFibreBundle.h.

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

double HFShowerFibreBundle::facCone
private

Definition at line 36 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

double HFShowerFibreBundle::facTube
private

Definition at line 36 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

int HFShowerFibreBundle::indexF
private

Definition at line 37 of file HFShowerFibreBundle.h.

Referenced by getHits(), and getRadius().

int HFShowerFibreBundle::indexR
private

Definition at line 37 of file HFShowerFibreBundle.h.

Referenced by getHits(), and getRadius().

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

Definition at line 39 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 40 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 39 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 40 of file HFShowerFibreBundle.h.

Referenced by getHits(), and HFShowerFibreBundle().

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

Definition at line 38 of file HFShowerFibreBundle.h.

Referenced by getRadius(), and HFShowerFibreBundle().