CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MaterialBudgetVolumeAnalysis.cc
Go to the documentation of this file.
3 
13 
15 
16 #include "TProfile.h"
17 #include "TProfile2D.h"
18 
19 #include <iostream>
20 #include <string>
21 #include <vector>
22 
23 //#define EDM_ML_DEBUG
24 
25 using namespace geant_units::operators;
26 
27 class MaterialBudgetVolumeAnalysis : public edm::one::EDAnalyzer<edm::one::SharedResources> {
28 public:
30 
31  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
32 
33 private:
34  void analyze(edm::Event const&, edm::EventSetup const&) override;
35  void bookHisto();
36 
37  const std::vector<std::string> names_;
39  const int binEta_, binPhi_;
40  const double etaLow_, etaHigh_, phiLow_, phiHigh_;
42  std::vector<TProfile*> meStepEta_, meStepPhi_;
43  std::vector<TProfile*> meRadLEta_, meRadLPhi_;
44  std::vector<TProfile*> meIntLEta_, meIntLPhi_;
45  std::vector<TProfile2D*> meStepEtaPhi_, meRadLEtaPhi_, meIntLEtaPhi_;
46 };
47 
49  : names_(p.getParameter<std::vector<std::string> >("names")),
50  tag_(p.getParameter<edm::InputTag>("inputTag")),
51  binEta_(p.getParameter<int>("nBinEta")),
52  binPhi_(p.getParameter<int>("nBinPhi")),
53  etaLow_(p.getParameter<double>("etaLow")),
54  etaHigh_(p.getParameter<double>("etaHigh")),
55  phiLow_(-1._pi),
56  phiHigh_(1._pi) {
57  usesResource(TFileService::kSharedResource);
58  tok_info_ = consumes<edm::MaterialInformationContainer>(tag_);
59 
60  edm::LogVerbatim("MaterialBudget") << "MaterialBudgetVolumeAnalysis: Eta plot: NX " << binEta_ << " Range "
61  << -etaLow_ << ":" << etaHigh_ << " Phi plot: NX " << binPhi_ << " Range "
62  << -1._pi << ":" << 1._pi << " for " << names_.size() << " detectors from "
63  << tag_;
64  std::ostringstream st1;
65  for (unsigned int k = 0; k < names_.size(); ++k)
66  st1 << " [" << k << "] " << names_[k];
67  edm::LogVerbatim("MaterialBudget") << "MaterialBudgetVolume: " << st1.str();
68  bookHisto();
69 }
70 
73  std::vector<std::string> names = {
74  "BEAM", "BEAM1", "BEAM2", "BEAM3", "BEAM4", "Tracker", "ECAL", "HCal", "MUON", "VCAL", "MGNT", "OQUA", "CALOEC"};
75  desc.add<std::vector<std::string> >("names", names);
76  desc.add<edm::InputTag>("inputTag", edm::InputTag("g4SimHits", "MaterialInformation"));
77  desc.add<int>("nBinEta", 300);
78  desc.add<int>("nBinPhi", 180);
79  desc.add<double>("etaLow", -6.0);
80  desc.add<double>("etaHigh", 6.0);
81  descriptions.add("materialBudgetVolumeAnalysis", desc);
82 }
83 
85  edm::Handle<edm::MaterialInformationContainer> materialInformationContainer;
86  iEvent.getByToken(tok_info_, materialInformationContainer);
87 #ifdef EDM_ML_DEBUG
88  unsigned int nsize(0), ntot(0), nused(0);
89 #endif
90  if (materialInformationContainer.isValid()) {
91 #ifdef EDM_ML_DEBUG
92  nsize = materialInformationContainer->size();
93 #endif
94  for (const auto& it : *(materialInformationContainer.product())) {
95 #ifdef EDM_ML_DEBUG
96  ntot++;
97 #endif
98  if (std::find(names_.begin(), names_.end(), it.vname()) != names_.end()) {
99 #ifdef EDM_ML_DEBUG
100  nused++;
101 #endif
102  unsigned int k =
103  static_cast<unsigned int>(std::find(names_.begin(), names_.end(), it.vname()) - names_.begin());
104  meStepEta_[k]->Fill(it.trackEta(), it.stepLength());
105  meRadLEta_[k]->Fill(it.trackEta(), it.radiationLength());
106  meIntLEta_[k]->Fill(it.trackEta(), it.interactionLength());
107  meStepPhi_[k]->Fill(it.trackPhi(), it.stepLength());
108  meRadLPhi_[k]->Fill(it.trackPhi(), it.radiationLength());
109  meIntLPhi_[k]->Fill(it.trackPhi(), it.interactionLength());
110  meStepEtaPhi_[k]->Fill(it.trackEta(), it.trackPhi(), it.stepLength());
111  meRadLEtaPhi_[k]->Fill(it.trackEta(), it.trackPhi(), it.radiationLength());
112  meIntLEtaPhi_[k]->Fill(it.trackEta(), it.trackPhi(), it.interactionLength());
113  }
114  }
115  }
116 #ifdef EDM_ML_DEBUG
117  edm::LogVerbatim("MaterialBudget") << "MaterialInformation with " << nsize << ":" << ntot << " elements of which "
118  << nused << " are used";
119 #endif
120 }
121 
124  char name[40], title[100];
125  for (unsigned int k = 0; k < names_.size(); ++k) {
126  sprintf(name, "stepEta%s", names_[k].c_str());
127  sprintf(title, "MB(Step) vs #eta for %s", names_[k].c_str());
128  meStepEta_.emplace_back(fs->make<TProfile>(name, title, binEta_, etaLow_, etaHigh_));
129  sprintf(name, "radlEta%s", names_[k].c_str());
130  sprintf(title, "MB(X0) vs #eta for %s", names_[k].c_str());
131  meRadLEta_.emplace_back(fs->make<TProfile>(name, title, binEta_, etaLow_, etaHigh_));
132  sprintf(name, "intlEta%s", names_[k].c_str());
133  sprintf(title, "MB(L0) vs #eta for %s", names_[k].c_str());
134  meIntLEta_.emplace_back(fs->make<TProfile>(name, title, binEta_, etaLow_, etaHigh_));
135  sprintf(name, "stepPhi%s", names_[k].c_str());
136  sprintf(title, "MB(Step) vs #phi for %s", names_[k].c_str());
137  meStepPhi_.emplace_back(fs->make<TProfile>(name, title, binPhi_, phiLow_, phiHigh_));
138  sprintf(name, "radlPhi%s", names_[k].c_str());
139  sprintf(title, "MB(X0) vs #phi for %s", names_[k].c_str());
140  meRadLPhi_.emplace_back(fs->make<TProfile>(name, title, binPhi_, phiLow_, phiHigh_));
141  sprintf(name, "intlPhi%s", names_[k].c_str());
142  sprintf(title, "MB(L0) vs #phi for %s", names_[k].c_str());
143  meIntLPhi_.emplace_back(fs->make<TProfile>(name, title, binPhi_, phiLow_, phiHigh_));
144  sprintf(name, "stepEtaPhi%s", names_[k].c_str());
145  sprintf(title, "MB(Step) vs #eta and #phi for %s", names_[k].c_str());
146  meStepEtaPhi_.emplace_back(
147  fs->make<TProfile2D>(name, title, binEta_ / 2, etaLow_, etaHigh_, binPhi_ / 2, phiLow_, phiHigh_));
148  sprintf(name, "radlEtaPhi%s", names_[k].c_str());
149  sprintf(title, "MB(X0) vs #eta and #phi for %s", names_[k].c_str());
150  meRadLEtaPhi_.emplace_back(
151  fs->make<TProfile2D>(name, title, binEta_ / 2, etaLow_, etaHigh_, binPhi_ / 2, phiLow_, phiHigh_));
152  sprintf(name, "intlEtaPhi%s", names_[k].c_str());
153  sprintf(title, "MB(L0) vs #eta and #phi for %s", names_[k].c_str());
154  meIntLEtaPhi_.emplace_back(
155  fs->make<TProfile2D>(name, title, binEta_ / 2, etaLow_, etaHigh_, binPhi_ / 2, phiLow_, phiHigh_));
156  }
157 }
158 
159 // define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
std::vector< TProfile2D * > meStepEtaPhi_
std::vector< TProfile2D * > meIntLEtaPhi_
MaterialBudgetVolumeAnalysis(edm::ParameterSet const &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const std::string names[nVars_]
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void analyze(edm::Event const &, edm::EventSetup const &) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T const * product() const
Definition: Handle.h:70
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const std::vector< std::string > names_
std::vector< TProfile2D * > meRadLEtaPhi_
__shared__ uint32_t ntot
edm::EDGetTokenT< edm::MaterialInformationContainer > tok_info_