CMS 3D CMS Logo

MinBias.cc
Go to the documentation of this file.
1 // user include files
7 
8 namespace cms {
10  : caloGeometryESToken_(esConsumes<edm::Transition::BeginRun>()), geo_(nullptr) {
11  // get names of modules, producing object collections
12  hbheLabel_ = iConfig.getParameter<std::string>("hbheInput");
13  hoLabel_ = iConfig.getParameter<std::string>("hoInput");
14  hfLabel_ = iConfig.getParameter<std::string>("hfInput");
15  hbheToken_ = mayConsume<HBHERecHitCollection>(edm::InputTag(hbheLabel_));
16  hoToken_ = mayConsume<HORecHitCollection>(edm::InputTag(hoLabel_));
17  hfToken_ = mayConsume<HFRecHitCollection>(edm::InputTag(hfLabel_));
18  allowMissingInputs_ = iConfig.getUntrackedParameter<bool>("AllowMissingInputs", false);
19  }
20 
23  myTree_ = fs->make<TTree>("RecJet", "RecJet Tree");
24  myTree_->Branch("mydet", &mydet, "mydet/I");
25  myTree_->Branch("mysubd", &mysubd, "mysubd/I");
26  myTree_->Branch("depth", &depth, "depth/I");
27  myTree_->Branch("ieta", &ieta, "ieta/I");
28  myTree_->Branch("iphi", &iphi, "iphi/I");
29  myTree_->Branch("eta", &eta, "eta/F");
30  myTree_->Branch("phi", &phi, "phi/F");
31  myTree_->Branch("mom1", &mom1, "mom1/F");
32  myTree_->Branch("mom2", &mom2, "mom2/F");
33  myTree_->Branch("mom3", &mom3, "mom3/F");
34  myTree_->Branch("mom4", &mom4, "mom4/F");
35  }
36 
37  void MinBias::beginRun(edm::Run const&, edm::EventSetup const& iSetup) {
39  std::vector<DetId> did = geo_->getValidDetIds();
40 
41  for (auto const& id : did) {
42  if ((id).det() == DetId::Hcal) {
43  theFillDetMap0_[id] = 0.;
44  theFillDetMap1_[id] = 0.;
45  theFillDetMap2_[id] = 0.;
46  theFillDetMap3_[id] = 0.;
47  theFillDetMap4_[id] = 0.;
48  }
49  }
50  }
51 
52  void MinBias::endRun(edm::Run const&, edm::EventSetup const& iSetup) {}
53 
54  void MinBias::endJob() {
55  const HcalGeometry* hgeo = static_cast<const HcalGeometry*>(geo_->getSubdetectorGeometry(DetId::Hcal, 1));
56  const std::vector<DetId>& did = hgeo->getValidDetIds();
57  int i = 0;
58  for (const auto& id : did) {
59  // if( id.det() == DetId::Hcal ) {
60  GlobalPoint pos = hgeo->getPosition(id);
61  mydet = (int)(id.det());
62  mysubd = (id.subdetId());
63  depth = HcalDetId(id).depth();
64  ieta = HcalDetId(id).ieta();
65  iphi = HcalDetId(id).iphi();
66  phi = pos.phi();
67  eta = pos.eta();
68  if (theFillDetMap0_[id] > 0.) {
72  2. * pow(mom2, 3);
73  mom4 = (theFillDetMap4_[id] - 4. * mom1 * theFillDetMap3_[id] + 6. * pow(mom1, 2) * theFillDetMap2_[id]) /
74  theFillDetMap0_[id] -
75  3. * pow(mom1, 4);
76 
77  } else {
78  mom1 = 0.;
79  mom2 = 0.;
80  mom3 = 0.;
81  mom4 = 0.;
82  }
83  edm::LogWarning("MinBias") << " Detector " << id.rawId() << " mydet " << mydet << " " << mysubd << " " << depth
84  << " " << ieta << " " << iphi << " " << pos.eta() << " " << pos.phi();
85  edm::LogWarning("MinBias") << " Energy " << mom1 << " " << mom2 << std::endl;
86  myTree_->Fill();
87  i++;
88  // }
89  }
90  edm::LogWarning("MinBias") << " The number of CaloDet records " << did.size();
91  edm::LogWarning("MinBias") << " The number of Hcal records " << i << std::endl;
92 
93  /*
94  std::cout << "===== Start writing user histograms =====" << std::endl;
95  hOutputFile->SetCompressionLevel(2);
96  hOutputFile->cd();
97  myTree->Write();
98  hOutputFile->Close() ;
99  std::cout << "===== End writing user histograms =======" << std::endl;
100  */
101  }
102 
103  //
104  // member functions
105  //
106 
107  // ------------ method called to produce the data ------------
109  if (!hbheLabel_.empty()) {
111  iEvent.getByToken(hbheToken_, hbhe);
112  if (!hbhe.isValid()) {
113  // can't find it!
114  if (!allowMissingInputs_) {
115  *hbhe; // will throw the proper exception
116  }
117  } else {
118  for (auto const& hbheItr : (HBHERecHitCollection)(*hbhe)) {
119  DetId id = (hbheItr).detid();
120  if (hbheItr.energy() > 0.)
121  edm::LogWarning("MinBias") << " Energy = " << hbheItr.energy();
122  theFillDetMap0_[id] += 1.;
123  theFillDetMap1_[id] += hbheItr.energy();
124  theFillDetMap2_[id] += pow(hbheItr.energy(), 2);
125  theFillDetMap3_[id] += pow(hbheItr.energy(), 3);
126  theFillDetMap4_[id] += pow(hbheItr.energy(), 4);
127  }
128  }
129  }
130 
131  if (!hoLabel_.empty()) {
133  iEvent.getByToken(hoToken_, ho);
134  if (!ho.isValid()) {
135  // can't find it!
136  if (!allowMissingInputs_) {
137  *ho; // will throw the proper exception
138  }
139  } else {
140  for (auto const& hoItr : (HORecHitCollection)(*ho)) {
141  DetId id = hoItr.detid();
142  theFillDetMap0_[id] += 1.;
143  theFillDetMap1_[id] += hoItr.energy();
144  theFillDetMap2_[id] += pow(hoItr.energy(), 2);
145  theFillDetMap3_[id] += pow(hoItr.energy(), 3);
146  theFillDetMap4_[id] += pow(hoItr.energy(), 4);
147  }
148  }
149  }
150 
151  if (!hfLabel_.empty()) {
153  iEvent.getByToken(hfToken_, hf);
154  if (!hf.isValid()) {
155  // can't find it!
156  if (!allowMissingInputs_) {
157  *hf; // will throw the proper exception
158  }
159  } else {
160  for (auto const hfItr : (HFRecHitCollection)(*hf)) {
161  DetId id = hfItr.detid();
162  theFillDetMap0_[id] += 1.;
163  theFillDetMap1_[id] += hfItr.energy();
164  theFillDetMap2_[id] += pow(hfItr.energy(), 2);
165  theFillDetMap3_[id] += pow(hfItr.energy(), 3);
166  theFillDetMap4_[id] += pow(hfItr.energy(), 4);
167  }
168  }
169  }
170  }
171 } // namespace cms
float mom2
Definition: MinBias.h:51
float mom1
Definition: MinBias.h:51
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float mom4
Definition: MinBias.h:51
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
MinBias(const edm::ParameterSet &)
Definition: MinBias.cc:9
edm::EDGetTokenT< HFRecHitCollection > hfToken_
Definition: MinBias.h:41
void beginJob() override
Definition: MinBias.cc:21
std::map< DetId, double > theFillDetMap4_
Definition: MinBias.h:58
std::map< DetId, double > theFillDetMap3_
Definition: MinBias.h:57
std::string hoLabel_
Definition: MinBias.h:38
TTree * myTree_
Definition: MinBias.h:47
constexpr int pow(int x)
Definition: conifer.h:24
void endRun(edm::Run const &, edm::EventSetup const &) override
Definition: MinBias.cc:52
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryESToken_
Definition: MinBias.h:42
std::map< DetId, double > theFillDetMap0_
Definition: MinBias.h:54
int mydet
Definition: MinBias.h:49
T getUntrackedParameter(std::string const &, T const &) const
int ieta
Definition: MinBias.h:49
edm::EDGetTokenT< HORecHitCollection > hoToken_
Definition: MinBias.h:40
const CaloGeometry * geo_
Definition: MinBias.h:52
bool allowMissingInputs_
Definition: MinBias.h:44
int iEvent
Definition: GenABIO.cc:224
void endJob() override
Definition: MinBias.cc:54
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
std::map< DetId, double > theFillDetMap2_
Definition: MinBias.h:56
Transition
Definition: Transition.h:12
edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
Definition: MinBias.h:39
Namespace of DDCMS conversion namespace.
int depth
Definition: MinBias.h:49
std::map< DetId, double > theFillDetMap1_
Definition: MinBias.h:55
Definition: DetId.h:17
float eta
Definition: MinBias.h:50
float phi
Definition: MinBias.h:50
float mom3
Definition: MinBias.h:51
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: MinBias.cc:108
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:75
void beginRun(edm::Run const &, edm::EventSetup const &) override
Definition: MinBias.cc:37
HLT enums.
int iphi
Definition: MinBias.h:49
GlobalPoint getPosition(const DetId &id) const
std::string hbheLabel_
Definition: MinBias.h:38
int mysubd
Definition: MinBias.h:49
std::string hfLabel_
Definition: MinBias.h:38
Log< level::Warning, false > LogWarning
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
Definition: Run.h:45
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164