CMS 3D CMS Logo

MixCollectionValidation.cc
Go to the documentation of this file.
1 // system include files
3 
4 // user include files
6 
9 
11 
13 #include "TFile.h"
14 
15 #include <memory>
16 #include <utility>
17 
18 using namespace edm;
19 
21  : minbunch_(iConfig.getParameter<int>("minBunch")),
22  maxbunch_(iConfig.getParameter<int>("maxBunch")),
23  verbose_(iConfig.getUntrackedParameter<bool>("verbose", false)),
24  nbin_(maxbunch_ - minbunch_ + 1) {
25  // Histograms will be defined according to the configuration
26  ParameterSet mixObjextsSet_ = iConfig.getParameter<ParameterSet>("mixObjects");
27 }
28 
30  // do anything here that needs to be done at desctruction time
31  // (e.g. close files, deallocate resources etc.)
32 }
33 
35  edm::Run const &iRun,
36  edm::EventSetup const & /* iSetup */) {
37  iBooker.setCurrentFolder("MixingV/Mixing");
38 
39  std::vector<std::string> names = mixObjextsSet_.getParameterNames();
40 
41  for (std::vector<std::string>::iterator it = names.begin(); it != names.end(); ++it) {
43  if (!pset.exists("type"))
44  continue; // to allow replacement by empty pset
45  std::string object = pset.getParameter<std::string>("type");
46  std::vector<InputTag> tags = pset.getParameter<std::vector<InputTag>>("input");
47 
48  if (object == "HepMCProduct") {
49  std::string title = "Log10 Number of GenParticle in " + object;
50  std::string name = "NumberOf" + object;
51  nrHepMCProductH_ = iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, 0., 40.);
52 
54  if (!HepMCProductTags_.empty()) {
56  consumes<CrossingFrame<HepMCProduct>>(edm::InputTag("mix", HepMCProductTags_[0].label()));
57  }
58  } else if (object == "SimTrack") {
59  std::string title = "Log10 Number of " + object;
60  std::string name = "NumberOf" + object;
61  nrSimTrackH_ = iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, 0., 40.);
62 
64  if (!SimTrackTags_.empty()) {
65  crossingFrame_SimTr_Token_ = consumes<CrossingFrame<SimTrack>>(edm::InputTag("mix", SimTrackTags_[0].label()));
66  }
67  } else if (object == "SimVertex") {
68  std::string title = "Log10 Number of " + object;
69  std::string name = "NumberOf" + object;
70  nrSimVertexH_ = iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, 0., 40.);
71 
73  if (!SimVertexTags_.empty()) {
75  consumes<CrossingFrame<SimVertex>>(edm::InputTag("mix", SimVertexTags_[0].label()));
76  }
77  } else if (object == "PSimHit") {
78  std::vector<std::string> subdets = pset.getParameter<std::vector<std::string>>("subdets");
79  for (unsigned int ii = 0; ii < subdets.size(); ii++) {
80  std::string title = "Log10 Number of " + subdets[ii];
81  std::string name = "NumberOf" + subdets[ii];
82  SimHitNrmap_[subdets[ii]] = iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, 0., 40.);
83 
84  title = "Time of " + subdets[ii];
85  name = "TimeOf" + subdets[ii];
87  iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, -125., 375.);
88  }
89 
91  for (auto const &it : PSimHitTags_)
93  consumes<CrossingFrame<PSimHit>>(edm::InputTag("mix", it.label() + it.instance())));
94  } else if (object == "PCaloHit") {
95  std::vector<std::string> subdets = pset.getParameter<std::vector<std::string>>("subdets");
96  for (unsigned int ii = 0; ii < subdets.size(); ii++) {
97  std::string title = "Log10 Number of " + subdets[ii];
98  std::string name = "NumberOf" + subdets[ii];
99  CaloHitNrmap_[subdets[ii]] = iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, 0., 40.);
100 
101  title = "Time of " + subdets[ii];
102  name = "TimeOf" + subdets[ii];
104  iBooker.bookProfile(name, title, nbin_, minbunch_, maxbunch_ + 1, 40, -125., 375.);
105  }
106 
108  for (auto const &it : PCaloHitTags_)
110  consumes<CrossingFrame<PCaloHit>>(edm::InputTag("mix", it.label() + it.instance())));
111  }
112  }
113 }
114 
116  using namespace edm;
117 
118  if (!HepMCProductTags_.empty()) {
119  bool gotHepMCProduct;
121  gotHepMCProduct = iEvent.getByToken(crossingFrame_Hep_Token_, crossingFrame);
122 
123  if (gotHepMCProduct) {
124  std::unique_ptr<MixCollection<HepMCProduct>> hepMCProduct(
125  new MixCollection<HepMCProduct>(crossingFrame.product()));
127 
128  fillGenParticleMulti(hitItr, hepMCProduct, nrHepMCProductH_);
129  }
130  }
131 
132  if (!SimTrackTags_.empty()) {
133  bool gotSimTrack;
135  gotSimTrack = iEvent.getByToken(crossingFrame_SimTr_Token_, crossingFrame);
136 
137  if (gotSimTrack) {
138  std::unique_ptr<MixCollection<SimTrack>> simTracks(new MixCollection<SimTrack>(crossingFrame.product()));
140 
142  }
143  }
144 
145  if (!SimVertexTags_.empty()) {
146  bool gotSimVertex;
148  std::string SimVertexLabel = SimVertexTags_[0].label();
149  gotSimVertex = iEvent.getByToken(crossingFrame_SimVtx_Token_, crossingFrame);
150 
151  if (gotSimVertex) {
152  std::unique_ptr<MixCollection<SimVertex>> simVerteces(new MixCollection<SimVertex>(crossingFrame.product()));
154 
155  fillMultiplicity(hitItr, simVerteces, nrSimVertexH_);
156  }
157  }
158 
159  if (!PSimHitTags_.empty()) {
160  edm::Handle<CrossingFrame<PSimHit>> crossingFrame;
161 
162  for (int i = 0; i < (int)PSimHitTags_.size(); i++) {
163  bool gotPSimHit;
164  gotPSimHit = iEvent.getByToken(crossingFrame_PSimHit_Tokens_[i], crossingFrame);
165 
166  if (gotPSimHit) {
167  std::unique_ptr<MixCollection<PSimHit>> simHits(new MixCollection<PSimHit>(crossingFrame.product()));
168 
170 
172 
174  }
175  }
176  }
177 
178  if (!PCaloHitTags_.empty()) {
180 
181  for (int i = 0; i < (int)PCaloHitTags_.size(); i++) {
182  bool gotPCaloHit;
183  std::string PCaloHitLabel = PCaloHitTags_[i].label() + PCaloHitTags_[i].instance();
184  gotPCaloHit = iEvent.getByToken(crossingFrame_PCaloHit_Tokens_[i], crossingFrame);
185 
186  if (gotPCaloHit) {
187  std::unique_ptr<MixCollection<PCaloHit>> caloHits(new MixCollection<PCaloHit>(crossingFrame.product()));
188 
190 
191  fillMultiplicity(hitItr, caloHits, CaloHitNrmap_[PCaloHitTags_[i].instance()]);
192 
193  fillCaloHitTime(hitItr, caloHits, CaloHitTimemap_[PCaloHitTags_[i].instance()]);
194  }
195  }
196  }
197 }
198 
199 template <class T1, class T2>
200 void MixCollectionValidation::fillMultiplicity(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_) {
201  std::vector<int> theMult(nbin_);
202 
203  for (theItr_ = theColl_->begin(); theItr_ != theColl_->end(); ++theItr_) {
204  int bunch = (*theItr_).eventId().bunchCrossing();
205  int index = bunch - minbunch_;
206  if (index >= 0 && index < nbin_) {
207  theMult[index] += 1;
208  } else {
209  edm::LogWarning("MixCollectionValidation") << "fillMultiplicity: bunch number " << bunch << " out of range";
210  }
211  }
212 
213  for (int i = 0; i < nbin_; i++) {
214  theProfile_->Fill(float(i + minbunch_ + 0.5), std::log10(std::max(float(0.1), float(theMult[i]))));
215  }
216 }
217 
218 template <class T1, class T2>
219 void MixCollectionValidation::fillGenParticleMulti(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_) {
220  std::vector<int> theMult(nbin_);
221 
222  for (theItr_ = theColl_->begin(); theItr_ != theColl_->end(); ++theItr_) {
223  int bunch = theItr_.bunch();
224  int index = bunch - minbunch_;
225  if (index >= 0 && index < nbin_) {
226  theMult[index] += (*theItr_).GetEvent()->particles_size();
227  } else {
228  edm::LogWarning("MixCollectionValidation") << "fillMultiplicity: bunch number " << bunch << " out of range";
229  }
230  }
231 
232  for (int i = 0; i < nbin_; i++) {
233  theProfile_->Fill(float(i + minbunch_ + 0.5), std::log10(std::max(float(0.1), float(theMult[i]))));
234  }
235 }
236 
237 template <class T1, class T2>
238 void MixCollectionValidation::fillSimHitTime(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_) {
239  for (theItr_ = theColl_->begin(); theItr_ != theColl_->end(); ++theItr_) {
240  int bunch = (*theItr_).eventId().bunchCrossing();
241  float time = (*theItr_).timeOfFlight();
242  int index = bunch - minbunch_;
243  if (index >= 0 && index < nbin_) {
244  theProfile_->Fill(float(bunch + 0.5), time);
245  } else {
246  edm::LogWarning("MixCollectionValidation") << "fillSimHitTime: bunch number " << bunch << " out of range";
247  }
248  }
249 }
250 
251 template <class T1, class T2>
252 void MixCollectionValidation::fillCaloHitTime(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_) {
253  for (theItr_ = theColl_->begin(); theItr_ != theColl_->end(); ++theItr_) {
254  int bunch = (*theItr_).eventId().bunchCrossing();
255  float time = (*theItr_).time();
256  int index = bunch - minbunch_;
257  if (index >= 0 && index < nbin_) {
258  theProfile_->Fill(float(bunch + 0.5), time);
259  } else {
260  edm::LogWarning("MixCollectionValidation") << "fillCaloHitTime: bunch number " << bunch << " out of range";
261  }
262  }
263 }
std::map< std::string, MonitorElement * > SimHitNrmap_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< edm::EDGetTokenT< CrossingFrame< PSimHit > > > crossingFrame_PSimHit_Tokens_
std::map< std::string, MonitorElement * > CaloHitTimemap_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::vector< edm::InputTag > HepMCProductTags_
MixCollectionValidation(const edm::ParameterSet &)
static PFTauRenderPlugin instance
std::vector< edm::EDGetTokenT< CrossingFrame< PCaloHit > > > crossingFrame_PCaloHit_Tokens_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::map< std::string, MonitorElement * > CaloHitNrmap_
T const * product() const
Definition: Handle.h:70
const std::string names[nVars_]
void Fill(long long x)
char const * label
int iEvent
Definition: GenABIO.cc:224
void fillCaloHitTime(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_)
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
edm::EDGetTokenT< CrossingFrame< SimVertex > > crossingFrame_SimVtx_Token_
void fillMultiplicity(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_)
std::vector< edm::InputTag > PSimHitTags_
void fillSimHitTime(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_)
std::map< std::string, MonitorElement * > SimHitTimemap_
edm::EDGetTokenT< CrossingFrame< SimTrack > > crossingFrame_SimTr_Token_
std::vector< edm::InputTag > PCaloHitTags_
ii
Definition: cuy.py:589
std::vector< edm::InputTag > SimTrackTags_
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > crossingFrame_Hep_Token_
std::vector< edm::InputTag > SimVertexTags_
HLT enums.
Log< level::Warning, false > LogWarning
static const std::string subdets[7]
Definition: TrackUtils.cc:60
std::vector< std::string > getParameterNames() const
Definition: Run.h:45
void analyze(const edm::Event &, const edm::EventSetup &) override
void fillGenParticleMulti(T1 &theItr_, T2 &theColl_, MonitorElement *theProfile_)