CMS 3D CMS Logo

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

#include <MaterialBudgetHcalHistos.h>

Public Member Functions

void fillBeginJob (const DDCompactView &)
 
void fillEndTrack ()
 
void fillPerStep (const G4Step *)
 
void fillStartTrack (const G4Track *)
 
 MaterialBudgetHcalHistos (const edm::ParameterSet &p)
 
virtual ~MaterialBudgetHcalHistos ()
 

Private Member Functions

void book ()
 
void fillHisto (int ii)
 
void fillLayer ()
 
std::vector< double > getDDDArray (const std::string &str, const DDsvalues_type &sv)
 
std::vector< std::string > getNames (DDFilteredView &fv)
 
void hend ()
 
bool isItEC (std::string)
 
bool isItHF (const G4VTouchable *)
 
bool isSensitive (std::string)
 

Private Attributes

int binEta
 
int binPhi
 
double eta
 
double etaHigh
 
double etaLow
 
bool fillHistos
 
std::vector< int > hfLevels
 
std::vector< std::string > hfNames
 
int id
 
double intLen
 
std::vector< double > intLength
 
int layer
 
std::vector< std::string > matList
 
double maxEta
 
TProfile * me100 [maxSet]
 
TProfile2D * me1000 [maxSet]
 
TProfile2D * me1100 [maxSet]
 
TH2F * me1200 [maxSet]
 
TH1F * me1300 [maxSet2]
 
TH2F * me1400 [maxSet2]
 
TProfile * me1500 [maxSet2]
 
TProfile * me200 [maxSet]
 
TProfile * me300 [maxSet]
 
TH1F * me400 [maxSet]
 
TProfile * me500 [maxSet]
 
TProfile * me600 [maxSet]
 
TProfile * me700 [maxSet]
 
TH1F * me800 [maxSet]
 
TProfile2D * me900 [maxSet]
 
int nlayHB
 
int nlayHE
 
int nlayHF
 
int nlayHO
 
double phi
 
bool printSum
 
double radLen
 
std::vector< double > radLength
 
std::vector< std::string > sensitiveEC
 
std::vector< std::string > sensitives
 
double stepLen
 
std::vector< double > stepLength
 
int steps
 

Static Private Attributes

static const int maxSet = 25
 
static const int maxSet2 = 9
 

Detailed Description

Definition at line 19 of file MaterialBudgetHcalHistos.h.

Constructor & Destructor Documentation

MaterialBudgetHcalHistos::MaterialBudgetHcalHistos ( const edm::ParameterSet p)

Definition at line 17 of file MaterialBudgetHcalHistos.cc.

References binEta, binPhi, book(), etaHigh, etaLow, fillHistos, edm::ParameterSet::getUntrackedParameter(), maxEta, pi, and printSum.

17  {
18 
19  binEta = p.getUntrackedParameter<int>("NBinEta", 260);
20  binPhi = p.getUntrackedParameter<int>("NBinPhi", 180);
21  maxEta = p.getUntrackedParameter<double>("MaxEta", 5.2);
22  etaLow = p.getUntrackedParameter<double>("EtaLow", -5.2);
23  etaHigh = p.getUntrackedParameter<double>("EtaHigh", 5.2);
24  fillHistos = p.getUntrackedParameter<bool>("FillHisto", true);
25  printSum = p.getUntrackedParameter<bool>("PrintSummary", false);
26  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: FillHisto : "
27  << fillHistos << " PrintSummary " << printSum
28  << " == Eta plot: NX " << binEta << " Range "
29  << -maxEta << ":" << maxEta <<" Phi plot: NX "
30  << binPhi << " Range " << -pi << ":" << pi
31  << " (Eta limit " << etaLow << ":" << etaHigh
32  <<")";
33  if (fillHistos) book();
34 
35 }
T getUntrackedParameter(std::string const &, T const &) const
const Double_t pi
virtual MaterialBudgetHcalHistos::~MaterialBudgetHcalHistos ( )
inlinevirtual

Member Function Documentation

void MaterialBudgetHcalHistos::book ( )
private

Definition at line 253 of file MaterialBudgetHcalHistos.cc.

References binEta, binPhi, mps_fire::i, edm::Service< T >::isAvailable(), TFileService::make(), maxEta, MuonErrorMatrixAnalyzer_cfi::maxPhi, maxSet, maxSet2, me100, me1000, me1100, me1200, me1300, me1400, me1500, me200, me300, me400, me500, me600, me700, me800, me900, pi, AlCaHLTBitMon_QueryRunRegistry::string, and compare::tfile.

Referenced by MaterialBudgetHcalHistos(), and ~MaterialBudgetHcalHistos().

253  {
254 
255  // Book histograms
257 
258  if ( !tfile.isAvailable() )
259  throw cms::Exception("BadConfig") << "TFileService unavailable: "
260  << "please add it to config file";
261 
262  double maxPhi=pi;
263  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Booking user "
264  << "histos === with " << binEta << " bins "
265  << "in eta from " << -maxEta << " to "
266  << maxEta << " and " << binPhi << " bins "
267  << "in phi from " << -maxPhi << " to "
268  << maxPhi;
269 
270  std::string iter;
271  // total X0
272  for (int i=0; i<maxSet; i++) {
273  iter = std::to_string(i);
274  me100[i] = tfile->make<TProfile>(std::to_string(i + 100).c_str(),
275  ("MB(X0) prof Eta in region " + iter).c_str(), binEta, -maxEta, maxEta);
276  me200[i] = tfile->make<TProfile>(std::to_string(i + 200).c_str(),
277  ("MB(L0) prof Eta in region " + iter).c_str(), binEta, -maxEta, maxEta);
278  me300[i] = tfile->make<TProfile>(std::to_string(i + 300).c_str(),
279  ("MB(Step) prof Eta in region " + iter).c_str(), binEta, -maxEta, maxEta);
280  me400[i] = tfile->make<TH1F>(std::to_string(i + 400).c_str(),
281  ("Eta in region " + iter).c_str(), binEta, -maxEta, maxEta);
282  me500[i] = tfile->make<TProfile>(std::to_string(i + 500).c_str(),
283  ("MB(X0) prof Ph in region " + iter).c_str(), binPhi, -maxPhi, maxPhi);
284  me600[i] = tfile->make<TProfile>(std::to_string(i + 600).c_str(),
285  ("MB(L0) prof Ph in region " + iter).c_str(), binPhi, -maxPhi, maxPhi);
286  me700[i] = tfile->make<TProfile>(std::to_string(i + 700).c_str(),
287  ("MB(Step) prof Ph in region " + iter).c_str(), binPhi, -maxPhi, maxPhi);
288  me800[i] = tfile->make<TH1F>(std::to_string(i + 800).c_str(),
289  ("Phi in region " + iter).c_str(), binPhi, -maxPhi, maxPhi);
290  me900[i] = tfile->make<TProfile2D>(std::to_string(i + 900).c_str(),
291  ("MB(X0) prof Eta Phi in region " + iter).c_str(), binEta/2, -maxEta, maxEta,
292  binPhi/2, -maxPhi, maxPhi);
293  me1000[i]= tfile->make<TProfile2D>(std::to_string(i + 1000).c_str(),
294  ("MB(L0) prof Eta Phi in region " + iter).c_str(), binEta/2, -maxEta, maxEta,
295  binPhi/2, -maxPhi, maxPhi);
296  me1100[i]= tfile->make<TProfile2D>(std::to_string(i + 1100).c_str(),
297  ("MB(Step) prof Eta Phi in region " + iter).c_str(), binEta/2, -maxEta, maxEta,
298  binPhi/2, -maxPhi, maxPhi);
299  me1200[i]= tfile->make<TH2F>(std::to_string(i + 1200).c_str(),
300  ("Eta vs Phi in region " + iter).c_str(), binEta/2, -maxEta, maxEta,
301  binPhi/2, -maxPhi, maxPhi);
302  }
303  for (int i=0; i<maxSet2; i++) {
304  iter = std::to_string(i);
305  me1300[i]= tfile->make<TH1F>(std::to_string(i + 1300).c_str(),
306  ("Events with layers Hit (0 all, 1 HB, ..) for " + iter).c_str(), binEta, -maxEta, maxEta);
307  me1400[i]= tfile->make<TH2F>(std::to_string(i + 1400).c_str(),
308  ("Eta vs Phi for layers hit in " + iter).c_str(), binEta/2, -maxEta, maxEta,
309  binPhi/2, -maxPhi, maxPhi);
310  me1500[i]= tfile->make<TProfile>(std::to_string(i + 1500).c_str(),
311  ("Number of layers crossed (0 all, 1 HB, ..) for " + iter).c_str(), binEta, -maxEta, maxEta);
312  }
313 
314  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Booking user "
315  << "histos done ===";
316 
317 }
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
const Double_t pi
bool isAvailable() const
Definition: Service.h:46
void MaterialBudgetHcalHistos::fillBeginJob ( const DDCompactView cpv)

Definition at line 37 of file MaterialBudgetHcalHistos.cc.

References fillHistos, DDFilteredView::firstChild(), getDDDArray(), getNames(), hfLevels, hfNames, mps_fire::i, gen::k, hcalDigis_cfi::level, DDFilteredView::mergedSpecifics(), sensitiveEC, sensitives, AlCaHLTBitMon_QueryRunRegistry::string, and groupFilesInBlocks::temp.

Referenced by MaterialBudgetHcal::update(), and ~MaterialBudgetHcalHistos().

37  {
38 
39  if (fillHistos) {
40  std::string attribute = "ReadOutName";
41  std::string value = "HcalHits";
42  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
43  DDFilteredView fv1(cpv,filter1);
44  sensitives = getNames(fv1);
45  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Names to be "
46  << "tested for " << attribute << " = "
47  << value << " has " << sensitives.size()
48  << " elements";
49  for (unsigned int i=0; i<sensitives.size(); i++)
50  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: sensitives["
51  << i << "] = " << sensitives[i];
52 
53  attribute = "Volume";
54  value = "HF";
55  DDSpecificsMatchesValueFilter filter2{DDValue(attribute,value,0)};
56  DDFilteredView fv2(cpv,filter2);
57  hfNames = getNames(fv2);
58  fv2.firstChild();
59  DDsvalues_type sv(fv2.mergedSpecifics());
60  std::vector<double> temp = getDDDArray("Levels",sv);
61  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Names to be "
62  << "tested for " << attribute << " = "
63  << value << " has " << hfNames.size()
64  << " elements";
65  for (unsigned int i=0; i<hfNames.size(); i++) {
66  int level = static_cast<int>(temp[i]);
67  hfLevels.push_back(level);
68  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: HF[" << i
69  << "] = " << hfNames[i] << " at level "
70  << hfLevels[i];
71  }
72 
73  std::string ecalRO[2] = {"EcalHitsEB", "EcalHitsEE"};
74  attribute = "ReadOutName";
75  for (int k=0; k<2; k++) {
76  value = ecalRO[k];
77  DDSpecificsMatchesValueFilter filter3{DDValue(attribute,value,0)};
78  DDFilteredView fv3(cpv,filter3);
79  std::vector<std::string> senstmp = getNames(fv3);
80  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Names to be"
81  << " tested for " << attribute << " = "
82  << value << " has " << senstmp.size()
83  << " elements";
84  for (unsigned int i=0; i<senstmp.size(); i++)
85  sensitiveEC.push_back(senstmp[i]);
86  }
87  for (unsigned int i=0; i<sensitiveEC.size(); i++)
88  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos:sensitiveEC["
89  << i << "] = " << sensitiveEC[i];
90  }
91 }
std::vector< std::string > getNames(DDFilteredView &fv)
std::vector< std::string > hfNames
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< std::string > sensitiveEC
Definition: value.py:1
int k[5][pyjets_maxn]
std::vector< std::string > sensitives
std::vector< double > getDDDArray(const std::string &str, const DDsvalues_type &sv)
void MaterialBudgetHcalHistos::fillEndTrack ( )

Definition at line 236 of file MaterialBudgetHcalHistos.cc.

References fillHisto(), fillHistos, fillLayer(), cuy::ii, intLength, matList, maxSet, nlayHB, nlayHE, nlayHF, nlayHO, printSum, radLength, and stepLength.

Referenced by MaterialBudgetHcal::update(), and ~MaterialBudgetHcalHistos().

236  {
237  edm::LogInfo("MaterialBudget") << "Number of layers hit in HB:" << nlayHB
238  << " HE:" << nlayHE << " HO:" << nlayHO
239  << " HF:" << nlayHF;
240  if (fillHistos) {
241  fillHisto(maxSet-1);
242  fillLayer();
243  }
244  if (printSum) {
245  for (unsigned int ii=0; ii<matList.size(); ii++) {
246  edm::LogInfo("MaterialBudget") << matList[ii] << "\t" << stepLength[ii]
247  << "\t" << radLength[ii] << "\t"
248  << intLength[ii];
249  }
250  }
251 }
std::vector< double > stepLength
std::vector< std::string > matList
std::vector< double > radLength
ii
Definition: cuy.py:588
std::vector< double > intLength
void MaterialBudgetHcalHistos::fillHisto ( int  ii)
private

Definition at line 319 of file MaterialBudgetHcalHistos.cc.

References eta, etaHigh, etaLow, cuy::ii, intLen, LogDebug, maxSet, me100, me1000, me1100, me1200, me200, me300, me400, me500, me600, me700, me800, me900, phi, radLen, and stepLen.

Referenced by fillEndTrack(), fillPerStep(), and ~MaterialBudgetHcalHistos().

319  {
320 
321  LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos:FillHisto called "
322  << "with index " << ii << " integrated step "
323  << stepLen << " X0 " << radLen << " Lamda "
324  << intLen;
325 
326  if (ii >=0 && ii < maxSet) {
327  me100[ii]->Fill(eta, radLen);
328  me200[ii]->Fill(eta, intLen);
329  me300[ii]->Fill(eta, stepLen);
330  me400[ii]->Fill(eta);
331 
332  if (eta >= etaLow && eta <= etaHigh) {
333  me500[ii]->Fill(phi, radLen);
334  me600[ii]->Fill(phi, intLen);
335  me700[ii]->Fill(phi, stepLen);
336  me800[ii]->Fill(phi);
337  }
338 
339  me900[ii]->Fill(eta, phi, radLen);
340  me1000[ii]->Fill(eta, phi, intLen);
341  me1100[ii]->Fill(eta, phi, stepLen);
342  me1200[ii]->Fill(eta, phi);
343 
344  }
345 }
#define LogDebug(id)
ii
Definition: cuy.py:588
void MaterialBudgetHcalHistos::fillLayer ( )
private

Definition at line 347 of file MaterialBudgetHcalHistos.cc.

References funct::abs(), eta, me1300, me1400, me1500, nlayHB, nlayHE, nlayHF, nlayHO, and phi.

Referenced by fillEndTrack(), and ~MaterialBudgetHcalHistos().

347  {
348 
349  me1300[0]->Fill(eta);
350  me1400[0]->Fill(eta,phi);
351  if (nlayHB > 0) {
352  me1300[1]->Fill(eta);
353  me1400[1]->Fill(eta,phi);
354  }
355  if (nlayHB >= 16) {
356  me1300[2]->Fill(eta);
357  me1400[2]->Fill(eta,phi);
358  }
359  if (nlayHE > 0) {
360  me1300[3]->Fill(eta);
361  me1400[3]->Fill(eta,phi);
362  }
363  if (nlayHE >= 16) {
364  me1300[4]->Fill(eta);
365  me1400[4]->Fill(eta,phi);
366  }
367  if (nlayHO > 0) {
368  me1300[5]->Fill(eta);
369  me1400[5]->Fill(eta,phi);
370  }
371  if (nlayHO >= 2) {
372  me1300[6]->Fill(eta);
373  me1400[6]->Fill(eta,phi);
374  }
375  if (nlayHF > 0) {
376  me1300[7]->Fill(eta);
377  me1400[7]->Fill(eta,phi);
378  }
379  if (nlayHB > 0 || nlayHE > 0 || (nlayHF > 0 && std::abs(eta) > 3.0)) {
380  me1300[8]->Fill(eta);
381  me1400[8]->Fill(eta,phi);
382  }
383  me1500[0]->Fill(eta,(double)(nlayHB+nlayHO+nlayHE+nlayHF));
384  me1500[1]->Fill(eta,(double)(nlayHB));
385  me1500[2]->Fill(eta,(double)(nlayHE));
386  me1500[4]->Fill(eta,(double)(nlayHF));
387 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void MaterialBudgetHcalHistos::fillPerStep ( const G4Step *  aStep)

Definition at line 124 of file MaterialBudgetHcalHistos.cc.

References funct::abs(), eta, fillHisto(), fillHistos, runEdmFileComparison::found, g, id, cuy::ii, intLen, intLength, isItEC(), isItHF(), isSensitive(), layer, LogDebug, matList, dataset::name, nlayHB, nlayHE, nlayHF, nlayHO, phi, printSum, radLen, radLength, stepLen, stepLength, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by MaterialBudgetHcal::update(), and ~MaterialBudgetHcalHistos().

124  {
125 
126  G4Material * material = aStep->GetPreStepPoint()->GetMaterial();
127  double step = aStep->GetStepLength();
128  double radl = material->GetRadlen();
129  double intl = material->GetNuclearInterLength();
130  double density = material->GetDensity() / (g/cm3);
131 
132  int idOld = id;
133  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
134  std::string name = touch->GetVolume(0)->GetName();
135  std::string matName = material->GetName();
136  if (printSum) {
137  bool found = false;
138  for (unsigned int ii=0; ii<matList.size(); ii++) {
139  if (matList[ii] == matName) {
140  stepLength[ii] += step;
141  radLength[ii] += (step/radl);
142  intLength[ii] += (step/intl);
143  found = true;
144  break;
145  }
146  }
147  if (!found) {
148  matList.push_back(matName);
149  stepLength.push_back(step);
150  radLength.push_back(step/radl);
151  intLength.push_back(step/intl);
152  }
153  edm::LogInfo("MaterialBudget") << name << " " << step << " " << matName
154  << " " << stepLen << " " << step/radl << " "
155  << radLen << " " <<step/intl << " " <<intLen;
156  } else {
157  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Step at "
158  << name << " Length " << step << " in "
159  << matName << " of density " << density
160  << " g/cc; Radiation Length " <<radl <<" mm;"
161  << " Interaction Length " << intl << " mm\n"
162  << " Position "
163  << aStep->GetPreStepPoint()->GetPosition()
164  << " Cylindrical R "
165  <<aStep->GetPreStepPoint()->GetPosition().perp()
166  << " Length (so far) " << stepLen << " L/X0 "
167  << step/radl << "/" << radLen << " L/Lambda "
168  << step/intl << "/" << intLen;
169  }
170 
171  int det=0, lay=0;
172  if (fillHistos) {
173  if (isItEC(name)) {
174  det = 1;
175  lay = 1;
176  } else {
177  if (isSensitive(name)) {
178  if (isItHF(touch)) {
179  det = 5;
180  lay = 21;
181  if (lay != layer) nlayHF++;
182  } else {
183  det = (touch->GetReplicaNumber(1))/1000;
184  lay = (touch->GetReplicaNumber(0)/10)%100 + 3;
185  if (det == 4) {
186  double abeta = std::abs(eta);
187  if (abeta < 1.479) lay = layer + 1;
188  else lay--;
189  if (lay < 3) lay = 3;
190  if (lay == layer) lay++;
191  if (lay > 20) lay = 20;
192  if (lay != layer) nlayHE++;
193  } else if (lay != layer) {
194  if (lay >= 20) nlayHO++;
195  else nlayHB++;
196  }
197  }
198  LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos: Det " << det
199  << " Layer " << lay << " Eta " << eta
200  << " Phi " << phi/deg;
201  } else if (layer == 1) {
202  det = -1;
203  lay = 2;
204  }
205  }
206  if (det != 0) {
207  if (lay != layer) {
208  id = lay;
209  layer = lay;
210  }
211  }
212 
213  if (id > idOld) {
214  // edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Step at " << name;
215  fillHisto(id-1);
216  }
217  }
218 
219  stepLen += step;
220  radLen += step/radl;
221  intLen += step/intl;
222  if (fillHistos) {
223  if (layer == 21 && det == 5) {
224  if (!isItHF(aStep->GetPostStepPoint()->GetTouchable())) {
225  LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos: After HF in "
226  << aStep->GetPostStepPoint()->GetTouchable()->GetVolume(0)->GetName();
227  fillHisto(id);
228  id++;
229  layer = 0;
230  }
231  }
232  }
233 }
#define LogDebug(id)
std::vector< double > stepLength
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
std::vector< std::string > matList
std::vector< double > radLength
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ii
Definition: cuy.py:588
step
std::vector< double > intLength
bool isItHF(const G4VTouchable *)
void MaterialBudgetHcalHistos::fillStartTrack ( const G4Track *  aTrack)

Definition at line 93 of file MaterialBudgetHcalHistos.cc.

References dir, eta, GeV, createfilelist::int, intLen, intLength, layer, matList, nlayHB, nlayHE, nlayHF, nlayHO, phi, printSum, radLen, radLength, stepLen, stepLength, and steps.

Referenced by MaterialBudgetHcal::update(), and ~MaterialBudgetHcalHistos().

93  {
94 
95  id = layer = steps = 0;
96  radLen = intLen = stepLen = 0;
97  nlayHB = nlayHE = nlayHF = nlayHO = 0;
98 
99  const G4ThreeVector& dir = aTrack->GetMomentum() ;
100  if (dir.theta() != 0 ) {
101  eta = dir.eta();
102  } else {
103  eta = -99;
104  }
105  phi = dir.phi();
106  double theEnergy = aTrack->GetTotalEnergy();
107  int theID = (int)(aTrack->GetDefinition()->GetPDGEncoding());
108 
109  if (printSum) {
110  matList.clear();
111  stepLength.clear();
112  radLength.clear();
113  intLength.clear();
114  }
115 
116  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Track "
117  << aTrack->GetTrackID() << " Code " << theID
118  << " Energy " << theEnergy/GeV << " GeV; Eta "
119  << eta << " Phi " << phi/deg << " PT "
120  << dir.perp()/GeV << " GeV *****";
121 }
const double GeV
Definition: MathUtil.h:16
std::vector< double > stepLength
std::vector< std::string > matList
std::vector< double > radLength
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< double > intLength
std::vector< double > MaterialBudgetHcalHistos::getDDDArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 410 of file MaterialBudgetHcalHistos.cc.

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

Referenced by fillBeginJob(), and ~MaterialBudgetHcalHistos().

411  {
412 
413  LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos:getDDDArray called "
414  << "for " << str;
415  DDValue value(str);
416  if (DDfetch(&sv,value)) {
417  LogDebug("MaterialBudget") << value;
418  const std::vector<double> & fvec = value.doubles();
419  int nval = fvec.size();
420  if (nval < 1) {
421  edm::LogError("MaterialBudget") << "MaterialBudgetHcalHistos : # of "
422  << str << " bins " << nval
423  << " < 1 ==> illegal";
424  throw cms::Exception("Unknown", "MaterialBudgetHcalHistos")
425  << "nval < 1 for array " << str << "\n";
426  }
427 
428  return fvec;
429  } else {
430  edm::LogError("MaterialBudget") << "MaterialBudgetHcalHistos : cannot get "
431  << "array " << str;
432  throw cms::Exception("Unknown", "MaterialBudgetHcalHistos")
433  << "cannot get array " << str << "\n";
434  }
435 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:80
Definition: value.py:1
std::vector< std::string > MaterialBudgetHcalHistos::getNames ( DDFilteredView fv)
private

Definition at line 394 of file MaterialBudgetHcalHistos.cc.

References DDFilteredView::firstChild(), mps_fire::i, cmsBatch::log, DDFilteredView::logicalPart(), DDName::name(), DDBase< N, C >::name(), DDFilteredView::next(), convertSQLiteXML::ok, AlCaHLTBitMon_QueryRunRegistry::string, and tmp.

Referenced by fillBeginJob(), and ~MaterialBudgetHcalHistos().

394  {
395 
396  std::vector<std::string> tmp;
397  bool dodet = fv.firstChild();
398  while (dodet) {
399  const DDLogicalPart & log = fv.logicalPart();
400  std::string namx = log.name().name();
401  bool ok = true;
402  for (unsigned int i=0; i<tmp.size(); i++)
403  if (namx == tmp[i]) ok = false;
404  if (ok) tmp.push_back(namx);
405  dodet = fv.next();
406  }
407  return tmp;
408 }
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const N & name() const
Definition: DDBase.h:78
bool next()
set current node to the next node in the filtered tree
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:92
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
bool firstChild()
set the current node to the first child ...
const std::string & name() const
Returns the name.
Definition: DDName.cc:90
void MaterialBudgetHcalHistos::hend ( )
private

Definition at line 389 of file MaterialBudgetHcalHistos.cc.

Referenced by ~MaterialBudgetHcalHistos().

389  {
390  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Save user "
391  << "histos ===";
392 }
bool MaterialBudgetHcalHistos::isItEC ( std::string  name)
private

Definition at line 459 of file MaterialBudgetHcalHistos.cc.

References sensitiveEC.

Referenced by fillPerStep(), and ~MaterialBudgetHcalHistos().

459  {
460 
461  std::vector<std::string>::const_iterator it = sensitiveEC.begin();
462  std::vector<std::string>::const_iterator itEnd = sensitiveEC.end();
463  for (; it != itEnd; ++it)
464  if (name == *it) return true;
465  return false;
466 }
std::vector< std::string > sensitiveEC
bool MaterialBudgetHcalHistos::isItHF ( const G4VTouchable *  touch)
private

Definition at line 446 of file MaterialBudgetHcalHistos.cc.

References hfLevels, hfNames, jetCorrFactors_cfi::levels, dataset::name, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by fillPerStep(), and ~MaterialBudgetHcalHistos().

446  {
447 
448  // std::vector<std::string>::const_iterator it = hfNames.begin();
449  int levels = ((touch->GetHistoryDepth())+1);
450  for (unsigned int it=0; it < hfNames.size(); it++) {
451  if (levels >= hfLevels[it]) {
452  std::string name = touch->GetVolume(levels-hfLevels[it])->GetName();
453  if (name == hfNames[it]) return true;
454  }
455  }
456  return false;
457 }
std::vector< std::string > hfNames
levels
correction levels
bool MaterialBudgetHcalHistos::isSensitive ( std::string  name)
private

Definition at line 437 of file MaterialBudgetHcalHistos.cc.

References sensitives.

Referenced by fillPerStep(), and ~MaterialBudgetHcalHistos().

437  {
438 
439  std::vector<std::string>::const_iterator it = sensitives.begin();
440  std::vector<std::string>::const_iterator itEnd = sensitives.end();
441  for (; it != itEnd; ++it)
442  if (name == *it) return true;
443  return false;
444 }
std::vector< std::string > sensitives

Member Data Documentation

int MaterialBudgetHcalHistos::binEta
private

Definition at line 50 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and MaterialBudgetHcalHistos().

int MaterialBudgetHcalHistos::binPhi
private

Definition at line 50 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and MaterialBudgetHcalHistos().

double MaterialBudgetHcalHistos::eta
private
double MaterialBudgetHcalHistos::etaHigh
private

Definition at line 51 of file MaterialBudgetHcalHistos.h.

Referenced by fillHisto(), and MaterialBudgetHcalHistos().

double MaterialBudgetHcalHistos::etaLow
private

Definition at line 51 of file MaterialBudgetHcalHistos.h.

Referenced by fillHisto(), and MaterialBudgetHcalHistos().

bool MaterialBudgetHcalHistos::fillHistos
private
std::vector<int> MaterialBudgetHcalHistos::hfLevels
private

Definition at line 48 of file MaterialBudgetHcalHistos.h.

Referenced by fillBeginJob(), and isItHF().

std::vector<std::string> MaterialBudgetHcalHistos::hfNames
private

Definition at line 47 of file MaterialBudgetHcalHistos.h.

Referenced by fillBeginJob(), and isItHF().

int MaterialBudgetHcalHistos::id
private

Definition at line 60 of file MaterialBudgetHcalHistos.h.

Referenced by fillPerStep().

double MaterialBudgetHcalHistos::intLen
private

Definition at line 61 of file MaterialBudgetHcalHistos.h.

Referenced by fillHisto(), fillPerStep(), and fillStartTrack().

std::vector<double> MaterialBudgetHcalHistos::intLength
private

Definition at line 53 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillPerStep(), and fillStartTrack().

int MaterialBudgetHcalHistos::layer
private
std::vector<std::string> MaterialBudgetHcalHistos::matList
private

Definition at line 52 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillPerStep(), and fillStartTrack().

double MaterialBudgetHcalHistos::maxEta
private

Definition at line 51 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and MaterialBudgetHcalHistos().

const int MaterialBudgetHcalHistos::maxSet = 25
staticprivate

Definition at line 46 of file MaterialBudgetHcalHistos.h.

Referenced by book(), fillEndTrack(), and fillHisto().

const int MaterialBudgetHcalHistos::maxSet2 = 9
staticprivate

Definition at line 46 of file MaterialBudgetHcalHistos.h.

Referenced by book().

TProfile* MaterialBudgetHcalHistos::me100[maxSet]
private

Definition at line 56 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile2D * MaterialBudgetHcalHistos::me1000[maxSet]
private

Definition at line 59 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile2D * MaterialBudgetHcalHistos::me1100[maxSet]
private

Definition at line 59 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TH2F* MaterialBudgetHcalHistos::me1200[maxSet]
private

Definition at line 55 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TH1F * MaterialBudgetHcalHistos::me1300[maxSet2]
private

Definition at line 54 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillLayer().

TH2F * MaterialBudgetHcalHistos::me1400[maxSet2]
private

Definition at line 55 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillLayer().

TProfile* MaterialBudgetHcalHistos::me1500[maxSet2]
private

Definition at line 58 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillLayer().

TProfile * MaterialBudgetHcalHistos::me200[maxSet]
private

Definition at line 56 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile * MaterialBudgetHcalHistos::me300[maxSet]
private

Definition at line 56 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TH1F* MaterialBudgetHcalHistos::me400[maxSet]
private

Definition at line 54 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile* MaterialBudgetHcalHistos::me500[maxSet]
private

Definition at line 57 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile * MaterialBudgetHcalHistos::me600[maxSet]
private

Definition at line 57 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile * MaterialBudgetHcalHistos::me700[maxSet]
private

Definition at line 57 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TH1F * MaterialBudgetHcalHistos::me800[maxSet]
private

Definition at line 54 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

TProfile2D* MaterialBudgetHcalHistos::me900[maxSet]
private

Definition at line 59 of file MaterialBudgetHcalHistos.h.

Referenced by book(), and fillHisto().

int MaterialBudgetHcalHistos::nlayHB
private

Definition at line 63 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillLayer(), fillPerStep(), and fillStartTrack().

int MaterialBudgetHcalHistos::nlayHE
private

Definition at line 63 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillLayer(), fillPerStep(), and fillStartTrack().

int MaterialBudgetHcalHistos::nlayHF
private

Definition at line 63 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillLayer(), fillPerStep(), and fillStartTrack().

int MaterialBudgetHcalHistos::nlayHO
private

Definition at line 63 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillLayer(), fillPerStep(), and fillStartTrack().

double MaterialBudgetHcalHistos::phi
private
bool MaterialBudgetHcalHistos::printSum
private
double MaterialBudgetHcalHistos::radLen
private

Definition at line 61 of file MaterialBudgetHcalHistos.h.

Referenced by fillHisto(), fillPerStep(), and fillStartTrack().

std::vector<double> MaterialBudgetHcalHistos::radLength
private

Definition at line 53 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillPerStep(), and fillStartTrack().

std::vector<std::string> MaterialBudgetHcalHistos::sensitiveEC
private

Definition at line 47 of file MaterialBudgetHcalHistos.h.

Referenced by fillBeginJob(), and isItEC().

std::vector<std::string> MaterialBudgetHcalHistos::sensitives
private

Definition at line 47 of file MaterialBudgetHcalHistos.h.

Referenced by fillBeginJob(), and isSensitive().

double MaterialBudgetHcalHistos::stepLen
private

Definition at line 61 of file MaterialBudgetHcalHistos.h.

Referenced by fillHisto(), fillPerStep(), and fillStartTrack().

std::vector<double> MaterialBudgetHcalHistos::stepLength
private

Definition at line 53 of file MaterialBudgetHcalHistos.h.

Referenced by fillEndTrack(), fillPerStep(), and fillStartTrack().

int MaterialBudgetHcalHistos::steps
private

Definition at line 60 of file MaterialBudgetHcalHistos.h.

Referenced by fillStartTrack().