CMS 3D CMS Logo

SimAnalyzerMinbias.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <string>
4 #include <iostream>
5 #include <fstream>
6 #include <sstream>
7 #include <vector>
8 #include <map>
9 
10 // user include files
20 
25 
31 
33 
34 #include "TFile.h"
35 #include "TTree.h"
36 
37 // class declaration
39 public:
40  explicit SimAnalyzerMinbias(const edm::ParameterSet&);
41  ~SimAnalyzerMinbias() override;
42  void analyze(const edm::Event&, const edm::EventSetup&) override;
43  void beginJob() override;
44  void endJob() override;
45  void beginRun(const edm::Run& r, const edm::EventSetup& iSetup) override;
46  void endRun(const edm::Run& r, const edm::EventSetup& iSetup) override;
47 
48 private:
49  // ----------member data ---------------------------
51  double timeCut_;
52  TTree* myTree_;
53 
54  // Root tree members
57  struct myInfo {
59  myInfo() { theMB0 = theMB1 = theMB2 = theMB3 = theMB4 = 0; }
60  };
61  std::map<HcalDetId, myInfo> myMap_;
64 };
65 
66 // constructors and destructor
67 
69  timeCut_ = iConfig.getUntrackedParameter<double>("TimeCut", 500);
70 
71  // get token names of modules, producing object collections
72  tok_evt_ = consumes<edm::HepMCProduct>(edm::InputTag("generator"));
73  tok_hcal_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", "HcalHits"));
74 
75  edm::LogInfo("AnalyzerMB") << "Use Time cut of " << timeCut_ << " ns";
76 }
77 
79  // do anything here that needs to be done at desctruction time
80  // (e.g. close files, deallocate resources etc.)
81 }
82 
84 
85 void SimAnalyzerMinbias::endRun(const edm::Run& r, const edm::EventSetup& iSetup) {}
86 
88  myTree_ = fs_->make<TTree>("SimJet", "SimJet Tree");
89  myTree_->Branch("mydet", &mydet, "mydet/I");
90  myTree_->Branch("mysubd", &mysubd, "mysubd/I");
91  myTree_->Branch("cells", &cells, "cells");
92  myTree_->Branch("depth", &depth, "depth/I");
93  myTree_->Branch("ieta", &ieta, "ieta/I");
94  myTree_->Branch("iphi", &iphi, "iphi/I");
95  myTree_->Branch("mom0_MB", &mom0_MB, "mom0_MB/F");
96  myTree_->Branch("mom1_MB", &mom1_MB, "mom1_MB/F");
97  myTree_->Branch("mom2_MB", &mom2_MB, "mom2_MB/F");
98  myTree_->Branch("mom3_MB", &mom3_MB, "mom3_MB/F");
99  myTree_->Branch("mom4_MB", &mom4_MB, "mom4_MB/F");
100 
101  myMap_.clear();
102 }
103 
104 // EndJob
105 //
107  cells = 0;
108  for (std::map<HcalDetId, myInfo>::const_iterator itr = myMap_.begin(); itr != myMap_.end(); ++itr) {
109  mysubd = itr->first.subdet();
110  depth = itr->first.depth();
111  iphi = itr->first.iphi();
112  ieta = itr->first.ieta();
113  myInfo info = itr->second;
114  if (info.theMB0 > 0) {
115  mom0_MB = info.theMB0;
116  mom1_MB = info.theMB1;
117  mom2_MB = info.theMB2;
118  mom3_MB = info.theMB3;
119  mom4_MB = info.theMB4;
120  cells++;
121 
122  edm::LogInfo("AnalyzerMB") << " Result= " << mysubd << " " << ieta << " " << iphi << " mom0 " << mom0_MB
123  << " mom1 " << mom1_MB << " mom2 " << mom2_MB << " mom3 " << mom3_MB << " mom4 "
124  << mom4_MB;
125  myTree_->Fill();
126  }
127  }
128  edm::LogInfo("AnalyzerMB") << "cells " << cells;
129 }
130 
131 //
132 // member functions
133 //
134 
135 // ------------ method called to produce the data ------------
136 
138  edm::LogInfo("AnalyzerMB") << " Start SimAnalyzerMinbias::analyze " << iEvent.id().run() << ":"
139  << iEvent.id().event();
140 
142  iEvent.getByToken(tok_evt_, evtMC);
143  if (!evtMC.isValid()) {
144  edm::LogWarning("AnalyzerMB") << "no HepMCProduct found";
145  } else {
146  const HepMC::GenEvent* myGenEvent = evtMC->GetEvent();
147  edm::LogInfo("AnalyzerMB") << "Event with " << myGenEvent->particles_size() << " particles + "
148  << myGenEvent->vertices_size() << " vertices";
149  }
150 
152  iEvent.getByToken(tok_hcal_, hcalHits);
153  if (!hcalHits.isValid()) {
154  edm::LogWarning("AnalyzerMB") << "Error! can't get HcalHits product!";
155  return;
156  }
157 
158  const edm::PCaloHitContainer* HitHcal = hcalHits.product();
159  std::map<HcalDetId, double> hitMap;
160  for (std::vector<PCaloHit>::const_iterator hcalItr = HitHcal->begin(); hcalItr != HitHcal->end(); ++hcalItr) {
161  double time = hcalItr->time();
162  if (time < timeCut_) {
163  double energyhit = hcalItr->energy();
164  HcalDetId hid = HcalDetId(hcalItr->id());
165  std::map<HcalDetId, double>::iterator itr1 = hitMap.find(hid);
166  if (itr1 == hitMap.end()) {
167  hitMap[hid] = 0;
168  itr1 = hitMap.find(hid);
169  }
170  itr1->second += energyhit;
171  }
172  }
173  edm::LogInfo("AnalyzerMB") << "extract information of " << hitMap.size() << " towers from " << HitHcal->size()
174  << " hits";
175 
176  for (std::map<HcalDetId, double>::const_iterator hcalItr = hitMap.begin(); hcalItr != hitMap.end(); ++hcalItr) {
177  HcalDetId hid = hcalItr->first;
178  double energyhit = hcalItr->second;
179  std::map<HcalDetId, myInfo>::iterator itr1 = myMap_.find(hid);
180  if (itr1 == myMap_.end()) {
181  myInfo info;
182  myMap_[hid] = info;
183  itr1 = myMap_.find(hid);
184  }
185  itr1->second.theMB0++;
186  itr1->second.theMB1 += energyhit;
187  itr1->second.theMB2 += (energyhit * energyhit);
188  itr1->second.theMB3 += (energyhit * energyhit * energyhit);
189  itr1->second.theMB4 += (energyhit * energyhit * energyhit * energyhit);
190  edm::LogInfo("AnalyzerMB") << "ID " << hid << " with energy " << energyhit;
191  }
192 }
193 
194 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
void beginJob() override
std::vector< PCaloHit > PCaloHitContainer
static const TGPicture * info(bool iBackgroundIsBlack)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
~SimAnalyzerMinbias() override
void endRun(const edm::Run &r, const edm::EventSetup &iSetup) override
void beginRun(const edm::Run &r, const edm::EventSetup &iSetup) override
edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void endJob() override
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hcal_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
T const * product() const
Definition: Handle.h:74
edm::Service< TFileService > fs_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EventID id() const
Definition: EventBase.h:59
std::map< HcalDetId, myInfo > myMap_
SimAnalyzerMinbias(const edm::ParameterSet &)
Definition: Run.h:45