#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 |
static const int | maxSet2 = 9 |
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 |
Definition at line 19 of file MaterialBudgetHcalHistos.h.
MaterialBudgetHcalHistos::MaterialBudgetHcalHistos | ( | const edm::ParameterSet & | p | ) |
Definition at line 15 of file MaterialBudgetHcalHistos.cc.
References binEta, binPhi, book(), etaHigh, etaLow, fillHistos, edm::ParameterSet::getUntrackedParameter(), maxEta, pi, and printSum.
{ binEta = p.getUntrackedParameter<int>("NBinEta", 260); binPhi = p.getUntrackedParameter<int>("NBinPhi", 180); maxEta = p.getUntrackedParameter<double>("MaxEta", 5.2); etaLow = p.getUntrackedParameter<double>("EtaLow", -5.2); etaHigh = p.getUntrackedParameter<double>("EtaHigh", 5.2); fillHistos = p.getUntrackedParameter<bool>("FillHisto", true); printSum = p.getUntrackedParameter<bool>("PrintSummary", false); edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: FillHisto : " << fillHistos << " PrintSummary " << printSum << " == Eta plot: NX " << binEta << " Range " << -maxEta << ":" << maxEta <<" Phi plot: NX " << binPhi << " Range " << -pi << ":" << pi << " (Eta limit " << etaLow << ":" << etaHigh <<")"; if (fillHistos) book(); }
virtual MaterialBudgetHcalHistos::~MaterialBudgetHcalHistos | ( | ) | [inline, virtual] |
void MaterialBudgetHcalHistos::book | ( | ) | [private] |
Definition at line 260 of file MaterialBudgetHcalHistos.cc.
References binEta, binPhi, i, edm::Service< T >::isAvailable(), maxEta, maxSet, maxSet2, me100, me1000, me1100, me1200, me1300, me1400, me1500, me200, me300, me400, me500, me600, me700, me800, me900, mergeVDriftHistosByStation::name, pi, and indexGen::title.
Referenced by MaterialBudgetHcalHistos().
{ // Book histograms edm::Service<TFileService> tfile; if ( !tfile.isAvailable() ) throw cms::Exception("BadConfig") << "TFileService unavailable: " << "please add it to config file"; double maxPhi=pi; edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Booking user " << "histos === with " << binEta << " bins " << "in eta from " << -maxEta << " to " << maxEta << " and " << binPhi << " bins " << "in phi from " << -maxPhi << " to " << maxPhi; char name[10], title[40]; // total X0 for (int i=0; i<maxSet; i++) { sprintf(name, "%d", i+100); sprintf(title, "MB(X0) prof Eta in region %d", i); me100[i] = tfile->make<TProfile>(name, title, binEta, -maxEta, maxEta); sprintf(name, "%d", i+200); sprintf(title, "MB(L0) prof Eta in region %d", i); me200[i] = tfile->make<TProfile>(name, title, binEta, -maxEta, maxEta); sprintf(name, "%d", i+300); sprintf(title, "MB(Step) prof Eta in region %d", i); me300[i] = tfile->make<TProfile>(name, title, binEta, -maxEta, maxEta); sprintf(name, "%d", i+400); sprintf(title, "Eta in region %d", i); me400[i] = tfile->make<TH1F>(name, title, binEta, -maxEta, maxEta); sprintf(name, "%d", i+500); sprintf(title, "MB(X0) prof Ph in region %d", i); me500[i] = tfile->make<TProfile>(name, title, binPhi, -maxPhi, maxPhi); sprintf(name, "%d", i+600); sprintf(title, "MB(L0) prof Ph in region %d", i); me600[i] = tfile->make<TProfile>(name, title, binPhi, -maxPhi, maxPhi); sprintf(name, "%d", i+700); sprintf(title, "MB(Step) prof Ph in region %d", i); me700[i] = tfile->make<TProfile>(name, title, binPhi, -maxPhi, maxPhi); sprintf(name, "%d", i+800); sprintf(title, "Phi in region %d", i); me800[i] = tfile->make<TH1F>(name, title, binPhi, -maxPhi, maxPhi); sprintf(name, "%d", i+900); sprintf(title, "MB(X0) prof Eta Phi in region %d", i); me900[i] = tfile->make<TProfile2D>(name, title, binEta/2, -maxEta, maxEta, binPhi/2, -maxPhi, maxPhi); sprintf(name, "%d", i+1000); sprintf(title, "MB(L0) prof Eta Phi in region %d", i); me1000[i]= tfile->make<TProfile2D>(name, title, binEta/2, -maxEta, maxEta, binPhi/2, -maxPhi, maxPhi); sprintf(name, "%d", i+1100); sprintf(title, "MB(Step) prof Eta Phi in region %d", i); me1100[i]= tfile->make<TProfile2D>(name, title, binEta/2, -maxEta, maxEta, binPhi/2, -maxPhi, maxPhi); sprintf(name, "%d", i+1200); sprintf(title, "Eta vs Phi in region %d", i); me1200[i]= tfile->make<TH2F>(name, title, binEta/2, -maxEta, maxEta, binPhi/2, -maxPhi, maxPhi); } for (int i=0; i<maxSet2; i++) { sprintf(name, "%d", i+1300); sprintf(title, "Events with layers Hit (0 all, 1 HB, ..) for %d", i); me1300[i]= tfile->make<TH1F>(name, title, binEta, -maxEta, maxEta); sprintf(name, "%d", i+1400); sprintf(title, "Eta vs Phi for layers hit in %d", i); me1400[i]= tfile->make<TH2F>(name, title, binEta/2, -maxEta, maxEta, binPhi/2, -maxPhi, maxPhi); sprintf(name, "%d", i+1500); sprintf(title, "Number of layers crossed (0 all, 1 HB, ..) for %d", i); me1500[i]= tfile->make<TProfile>(name, title, binEta, -maxEta, maxEta); } edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Booking user " << "histos done ==="; }
void MaterialBudgetHcalHistos::fillBeginJob | ( | const DDCompactView & | cpv | ) |
Definition at line 35 of file MaterialBudgetHcalHistos.cc.
References DDFilteredView::addFilter(), DDSpecificsFilter::equals, fillHistos, DDFilteredView::firstChild(), getDDDArray(), getNames(), hfLevels, hfNames, i, gen::k, testEve_cfg::level, DDFilteredView::mergedSpecifics(), sensitiveEC, sensitives, DDSpecificsFilter::setCriteria(), cond::rpcobtemp::temp, and relativeConstraints::value.
Referenced by MaterialBudgetHcal::update().
{ if (fillHistos) { std::string attribute = "ReadOutName"; std::string value = "HcalHits"; DDSpecificsFilter filter1; DDValue ddv1(attribute,value,0); filter1.setCriteria(ddv1,DDSpecificsFilter::equals); DDFilteredView fv1(cpv); fv1.addFilter(filter1); sensitives = getNames(fv1); edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Names to be " << "tested for " << attribute << " = " << value << " has " << sensitives.size() << " elements"; for (unsigned int i=0; i<sensitives.size(); i++) edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: sensitives[" << i << "] = " << sensitives[i]; attribute = "Volume"; value = "HF"; DDSpecificsFilter filter2; DDValue ddv2(attribute,value,0); filter2.setCriteria(ddv2,DDSpecificsFilter::equals); DDFilteredView fv2(cpv); fv2.addFilter(filter2); hfNames = getNames(fv2); fv2.firstChild(); DDsvalues_type sv(fv2.mergedSpecifics()); std::vector<double> temp = getDDDArray("Levels",sv); edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Names to be " << "tested for " << attribute << " = " << value << " has " << hfNames.size() << " elements"; for (unsigned int i=0; i<hfNames.size(); i++) { int level = static_cast<int>(temp[i]); hfLevels.push_back(level); edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: HF[" << i << "] = " << hfNames[i] << " at level " << hfLevels[i]; } std::string ecalRO[2] = {"EcalHitsEB", "EcalHitsEE"}; attribute = "ReadOutName"; for (int k=0; k<2; k++) { value = ecalRO[k]; DDSpecificsFilter filter3; DDValue ddv3(attribute,value,0); filter3.setCriteria(ddv3,DDSpecificsFilter::equals); DDFilteredView fv3(cpv); fv3.addFilter(filter3); std::vector<std::string> senstmp = getNames(fv3); edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Names to be" << " tested for " << attribute << " = " << value << " has " << senstmp.size() << " elements"; for (unsigned int i=0; i<senstmp.size(); i++) sensitiveEC.push_back(senstmp[i]); } for (unsigned int i=0; i<sensitiveEC.size(); i++) edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos:sensitiveEC[" << i << "] = " << sensitiveEC[i]; } }
void MaterialBudgetHcalHistos::fillEndTrack | ( | ) |
Definition at line 243 of file MaterialBudgetHcalHistos.cc.
References fillHisto(), fillHistos, fillLayer(), intLength, matList, maxSet, nlayHB, nlayHE, nlayHF, nlayHO, printSum, radLength, and stepLength.
Referenced by MaterialBudgetHcal::update().
{ edm::LogInfo("MaterialBudget") << "Number of layers hit in HB:" << nlayHB << " HE:" << nlayHE << " HO:" << nlayHO << " HF:" << nlayHF; if (fillHistos) { fillHisto(maxSet-1); fillLayer(); } if (printSum) { for (unsigned int ii=0; ii<matList.size(); ii++) { edm::LogInfo("MaterialBudget") << matList[ii] << "\t" << stepLength[ii] << "\t" << radLength[ii] << "\t" << intLength[ii]; } } }
void MaterialBudgetHcalHistos::fillHisto | ( | int | ii | ) | [private] |
Definition at line 339 of file MaterialBudgetHcalHistos.cc.
References eta, etaHigh, etaLow, intLen, LogDebug, maxSet, me100, me1000, me1100, me1200, me200, me300, me400, me500, me600, me700, me800, me900, phi, radLen, and stepLen.
Referenced by fillEndTrack(), and fillPerStep().
{ LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos:FillHisto called " << "with index " << ii << " integrated step " << stepLen << " X0 " << radLen << " Lamda " << intLen; if (ii >=0 && ii < maxSet) { me100[ii]->Fill(eta, radLen); me200[ii]->Fill(eta, intLen); me300[ii]->Fill(eta, stepLen); me400[ii]->Fill(eta); if (eta >= etaLow && eta <= etaHigh) { me500[ii]->Fill(phi, radLen); me600[ii]->Fill(phi, intLen); me700[ii]->Fill(phi, stepLen); me800[ii]->Fill(phi); } me900[ii]->Fill(eta, phi, radLen); me1000[ii]->Fill(eta, phi, intLen); me1100[ii]->Fill(eta, phi, stepLen); me1200[ii]->Fill(eta, phi); } }
void MaterialBudgetHcalHistos::fillLayer | ( | ) | [private] |
Definition at line 367 of file MaterialBudgetHcalHistos.cc.
References abs, eta, me1300, me1400, me1500, nlayHB, nlayHE, nlayHF, nlayHO, and phi.
Referenced by fillEndTrack().
{ me1300[0]->Fill(eta); me1400[0]->Fill(eta,phi); if (nlayHB > 0) { me1300[1]->Fill(eta); me1400[1]->Fill(eta,phi); } if (nlayHB >= 16) { me1300[2]->Fill(eta); me1400[2]->Fill(eta,phi); } if (nlayHE > 0) { me1300[3]->Fill(eta); me1400[3]->Fill(eta,phi); } if (nlayHE >= 16) { me1300[4]->Fill(eta); me1400[4]->Fill(eta,phi); } if (nlayHO > 0) { me1300[5]->Fill(eta); me1400[5]->Fill(eta,phi); } if (nlayHO >= 2) { me1300[6]->Fill(eta); me1400[6]->Fill(eta,phi); } if (nlayHF > 0) { me1300[7]->Fill(eta); me1400[7]->Fill(eta,phi); } if (nlayHB > 0 || nlayHE > 0 || (nlayHF > 0 && std::abs(eta) > 3.0)) { me1300[8]->Fill(eta); me1400[8]->Fill(eta,phi); } me1500[0]->Fill(eta,(double)(nlayHB+nlayHO+nlayHE+nlayHF)); me1500[1]->Fill(eta,(double)(nlayHB)); me1500[2]->Fill(eta,(double)(nlayHE)); me1500[4]->Fill(eta,(double)(nlayHF)); }
void MaterialBudgetHcalHistos::fillPerStep | ( | const G4Step * | aStep | ) |
Definition at line 131 of file MaterialBudgetHcalHistos.cc.
References abs, eta, fillHisto(), fillHistos, newFWLiteAna::found, g, id, intLen, intLength, isItEC(), isItHF(), isSensitive(), layer, LogDebug, matList, mergeVDriftHistosByStation::name, nlayHB, nlayHE, nlayHF, nlayHO, phi, printSum, radLen, radLength, launcher::step, stepLen, and stepLength.
Referenced by MaterialBudgetHcal::update().
{ G4Material * material = aStep->GetPreStepPoint()->GetMaterial(); double step = aStep->GetStepLength(); double radl = material->GetRadlen(); double intl = material->GetNuclearInterLength(); double density = material->GetDensity() / (g/cm3); int idOld = id; const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable(); std::string name = touch->GetVolume(0)->GetName(); std::string matName = material->GetName(); if (printSum) { bool found = false; for (unsigned int ii=0; ii<matList.size(); ii++) { if (matList[ii] == matName) { stepLength[ii] += step; radLength[ii] += (step/radl); intLength[ii] += (step/intl); found = true; break; } } if (!found) { matList.push_back(matName); stepLength.push_back(step); radLength.push_back(step/radl); intLength.push_back(step/intl); } edm::LogInfo("MaterialBudget") << name << " " << step << " " << matName << " " << stepLen << " " << step/radl << " " << radLen << " " <<step/intl << " " <<intLen; } else { edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Step at " << name << " Length " << step << " in " << matName << " of density " << density << " g/cc; Radiation Length " <<radl <<" mm;" << " Interaction Length " << intl << " mm\n" << " Position " << aStep->GetPreStepPoint()->GetPosition() << " Cylindrical R " <<aStep->GetPreStepPoint()->GetPosition().perp() << " Length (so far) " << stepLen << " L/X0 " << step/radl << "/" << radLen << " L/Lambda " << step/intl << "/" << intLen; } int det=0, lay=0; if (fillHistos) { if (isItEC(name)) { det = 1; lay = 1; } else { if (isSensitive(name)) { if (isItHF(touch)) { det = 5; lay = 21; if (lay != layer) nlayHF++; } else { det = (touch->GetReplicaNumber(1))/1000; lay = (touch->GetReplicaNumber(0)/10)%100 + 3; if (det == 4) { double abeta = std::abs(eta); if (abeta < 1.479) lay = layer + 1; else lay--; if (lay < 3) lay = 3; if (lay == layer) lay++; if (lay > 20) lay = 20; if (lay != layer) nlayHE++; } else if (lay != layer) { if (lay >= 20) nlayHO++; else nlayHB++; } } LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos: Det " << det << " Layer " << lay << " Eta " << eta << " Phi " << phi/deg; } else if (layer == 1) { det = -1; lay = 2; } } if (det != 0) { if (lay != layer) { id = lay; layer = lay; } } if (id > idOld) { // edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Step at " << name; fillHisto(id-1); } } stepLen += step; radLen += step/radl; intLen += step/intl; if (fillHistos) { if (layer == 21 && det == 5) { if (!isItHF(aStep->GetPostStepPoint()->GetTouchable())) { LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos: After HF in " << aStep->GetPostStepPoint()->GetTouchable()->GetVolume(0)->GetName(); fillHisto(id); id++; layer = 0; } } } }
void MaterialBudgetHcalHistos::fillStartTrack | ( | const G4Track * | aTrack | ) |
Definition at line 100 of file MaterialBudgetHcalHistos.cc.
References dir, eta, intLen, intLength, layer, matList, nlayHB, nlayHE, nlayHF, nlayHO, phi, printSum, radLen, radLength, stepLen, stepLength, and steps.
Referenced by MaterialBudgetHcal::update().
{ id = layer = steps = 0; radLen = intLen = stepLen = 0; nlayHB = nlayHE = nlayHF = nlayHO = 0; const G4ThreeVector& dir = aTrack->GetMomentum() ; if (dir.theta() != 0 ) { eta = dir.eta(); } else { eta = -99; } phi = dir.phi(); double theEnergy = aTrack->GetTotalEnergy(); int theID = (int)(aTrack->GetDefinition()->GetPDGEncoding()); if (printSum) { matList.clear(); stepLength.clear(); radLength.clear(); intLength.clear(); } edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Track " << aTrack->GetTrackID() << " Code " << theID << " Energy " << theEnergy/GeV << " GeV; Eta " << eta << " Phi " << phi/deg << " PT " << dir.perp()/GeV << " GeV *****"; }
std::vector< double > MaterialBudgetHcalHistos::getDDDArray | ( | const std::string & | str, |
const DDsvalues_type & | sv | ||
) | [private] |
Definition at line 430 of file MaterialBudgetHcalHistos.cc.
References DDfetch(), DDValue::doubles(), Exception, LogDebug, and relativeConstraints::value.
Referenced by fillBeginJob().
{ LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos:getDDDArray called " << "for " << str; DDValue value(str); if (DDfetch(&sv,value)) { LogDebug("MaterialBudget") << value; const std::vector<double> & fvec = value.doubles(); int nval = fvec.size(); if (nval < 1) { edm::LogError("MaterialBudget") << "MaterialBudgetHcalHistos : # of " << str << " bins " << nval << " < 1 ==> illegal"; throw cms::Exception("Unknown", "MaterialBudgetHcalHistos") << "nval < 1 for array " << str << "\n"; } return fvec; } else { edm::LogError("MaterialBudget") << "MaterialBudgetHcalHistos : cannot get " << "array " << str; throw cms::Exception("Unknown", "MaterialBudgetHcalHistos") << "cannot get array " << str << "\n"; } }
std::vector< std::string > MaterialBudgetHcalHistos::getNames | ( | DDFilteredView & | fv | ) | [private] |
Definition at line 414 of file MaterialBudgetHcalHistos.cc.
References DDFilteredView::firstChild(), i, funct::log(), DDFilteredView::logicalPart(), DDName::name(), DDBase< N, C >::name(), DDFilteredView::next(), convertSQLiteXML::ok, and tmp.
Referenced by fillBeginJob().
{ std::vector<std::string> tmp; bool dodet = fv.firstChild(); while (dodet) { const DDLogicalPart & log = fv.logicalPart(); std::string namx = log.name().name(); bool ok = true; for (unsigned int i=0; i<tmp.size(); i++) if (namx == tmp[i]) ok = false; if (ok) tmp.push_back(namx); dodet = fv.next(); } return tmp; }
void MaterialBudgetHcalHistos::hend | ( | ) | [private] |
Definition at line 409 of file MaterialBudgetHcalHistos.cc.
Referenced by ~MaterialBudgetHcalHistos().
{ edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Save user " << "histos ==="; }
bool MaterialBudgetHcalHistos::isItEC | ( | std::string | name | ) | [private] |
Definition at line 478 of file MaterialBudgetHcalHistos.cc.
References sensitiveEC.
Referenced by fillPerStep().
{ std::vector<std::string>::const_iterator it = sensitiveEC.begin(); for (; it != sensitiveEC.end(); it++) if (name == *it) return true; return false; }
bool MaterialBudgetHcalHistos::isItHF | ( | const G4VTouchable * | touch | ) | [private] |
Definition at line 465 of file MaterialBudgetHcalHistos.cc.
References hfLevels, hfNames, and mergeVDriftHistosByStation::name.
Referenced by fillPerStep().
{ std::vector<std::string>::const_iterator it = hfNames.begin(); int levels = ((touch->GetHistoryDepth())+1); for (unsigned int it=0; it < hfNames.size(); it++) { if (levels >= hfLevels[it]) { std::string name = touch->GetVolume(levels-hfLevels[it])->GetName(); if (name == hfNames[it]) return true; } } return false; }
bool MaterialBudgetHcalHistos::isSensitive | ( | std::string | name | ) | [private] |
Definition at line 457 of file MaterialBudgetHcalHistos.cc.
References sensitives.
Referenced by fillPerStep().
{ std::vector<std::string>::const_iterator it = sensitives.begin(); for (; it != sensitives.end(); it++) if (name == *it) return true; return false; }
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] |
Definition at line 62 of file MaterialBudgetHcalHistos.h.
Referenced by fillHisto(), fillLayer(), fillPerStep(), and fillStartTrack().
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] |
Definition at line 49 of file MaterialBudgetHcalHistos.h.
Referenced by fillBeginJob(), fillEndTrack(), fillPerStep(), and MaterialBudgetHcalHistos().
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] |
Definition at line 60 of file MaterialBudgetHcalHistos.h.
Referenced by fillPerStep(), and fillStartTrack().
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 [static, private] |
Definition at line 46 of file MaterialBudgetHcalHistos.h.
Referenced by book(), fillEndTrack(), and fillHisto().
const int MaterialBudgetHcalHistos::maxSet2 = 9 [private] |
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] |
Definition at line 62 of file MaterialBudgetHcalHistos.h.
Referenced by fillHisto(), fillLayer(), fillPerStep(), and fillStartTrack().
bool MaterialBudgetHcalHistos::printSum [private] |
Definition at line 49 of file MaterialBudgetHcalHistos.h.
Referenced by fillEndTrack(), fillPerStep(), fillStartTrack(), and MaterialBudgetHcalHistos().
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().