CMS 3D CMS Logo

DigiCollectionProfiler.h
Go to the documentation of this file.
1 #ifndef DPGAnalysis_SiStripTools_DigiCollectionProfile_H
2 #define DPGAnalysis_SiStripTools_DigiCollectionProfile_H
3 
4 #include <vector>
13 
14 #include "TH1.h"
15 #include "TH2.h"
16 #include "TProfile.h"
17 
18 template <class T>
20 public:
24 
25  void fill(edm::Handle<T> digis,
26  const std::vector<TH1F*>&,
27  const std::vector<TProfile*>&,
28  const std::vector<TH2F*>&) const;
29 
30 private:
31  bool m_folded;
35 
36  std::vector<DetIdSelector> m_selections;
37 };
38 
39 template <class T>
41  : m_folded(false), m_want1dHisto(false), m_wantProfile(false), m_want2dHisto(false), m_selections() {}
42 
43 template <class T>
45  : m_folded(iConfig.getUntrackedParameter<bool>("foldedStrips", false)),
46  m_want1dHisto(iConfig.getUntrackedParameter<bool>("want1dHisto", true)),
47  m_wantProfile(iConfig.getUntrackedParameter<bool>("wantProfile", true)),
48  m_want2dHisto(iConfig.getUntrackedParameter<bool>("want2dHisto", false))
49 
50 {
51  std::vector<edm::ParameterSet> selconfigs = iConfig.getParameter<std::vector<edm::ParameterSet> >("selections");
52 
53  for (std::vector<edm::ParameterSet>::const_iterator selconfig = selconfigs.begin(); selconfig != selconfigs.end();
54  ++selconfig) {
55  DetIdSelector selection(*selconfig);
56  m_selections.push_back(selection);
57  }
58 }
59 
60 template <class T>
62  const std::vector<TH1F*>& hist,
63  const std::vector<TProfile*>& hprof,
64  const std::vector<TH2F*>& hist2d) const {}
65 
66 template <>
69  const std::vector<TH1F*>& hist,
70  const std::vector<TProfile*>& hprof,
71  const std::vector<TH2F*>& hist2d) const {
72  for (edm::DetSetVector<SiStripDigi>::const_iterator mod = digis->begin(); mod != digis->end(); mod++) {
73  for (unsigned int isel = 0; isel < m_selections.size(); ++isel) {
74  if (m_selections[isel].isSelected(mod->detId())) {
75  TH1F* tobefilled1d = nullptr;
76  TProfile* tobefilledprof = nullptr;
77  TH2F* tobefilled2d = nullptr;
78 
79  if (m_want1dHisto)
80  tobefilled1d = hist[isel];
81  if (m_wantProfile)
82  tobefilledprof = hprof[isel];
83  if (m_want2dHisto)
84  tobefilled2d = hist2d[isel];
85 
86  for (edm::DetSet<SiStripDigi>::const_iterator digi = mod->begin(); digi != mod->end(); digi++) {
87  if (digi->adc() > 0) {
88  unsigned int strip = digi->strip();
89  if (m_folded)
90  strip = strip % 256;
91  if (tobefilled1d)
92  tobefilled1d->Fill(strip);
93  if (tobefilledprof)
94  tobefilledprof->Fill(strip, digi->adc());
95  if (tobefilled2d)
96  tobefilled2d->Fill(strip, digi->adc());
97  }
98  }
99  }
100  }
101  }
102 }
103 
104 template <>
107  const std::vector<TH1F*>& hist,
108  const std::vector<TProfile*>& hprof,
109  const std::vector<TH2F*>& hist2d) const {
110  for (edm::DetSetVector<SiStripRawDigi>::const_iterator mod = digis->begin(); mod != digis->end(); mod++) {
111  for (unsigned int isel = 0; isel < m_selections.size(); ++isel) {
112  if (m_selections[isel].isSelected(mod->detId())) {
113  TH1F* tobefilled1d = nullptr;
114  TProfile* tobefilledprof = nullptr;
115  TH2F* tobefilled2d = nullptr;
116 
117  if (m_want1dHisto)
118  tobefilled1d = hist[isel];
119  if (m_wantProfile)
120  tobefilledprof = hprof[isel];
121  if (m_want2dHisto)
122  tobefilled2d = hist2d[isel];
123 
124  unsigned int istrip = 0;
125  for (edm::DetSet<SiStripRawDigi>::const_iterator digi = mod->begin(); digi != mod->end(); digi++, ++istrip) {
126  if (digi->adc() > 0) {
127  unsigned int strip = istrip;
128  if (m_folded)
129  strip = strip % 256;
130  if (tobefilled1d)
131  tobefilled1d->Fill(strip);
132  if (tobefilledprof)
133  tobefilledprof->Fill(strip, digi->adc());
134  if (tobefilled2d)
135  tobefilled2d->Fill(strip, digi->adc());
136  }
137  }
138  }
139  }
140  }
141 }
142 
143 template <>
146  const std::vector<TH1F*>& hist,
147  const std::vector<TProfile*>& hprof,
148  const std::vector<TH2F*>& hist2d) const {
149  for (edmNew::DetSetVector<SiStripCluster>::const_iterator mod = digis->begin(); mod != digis->end(); mod++) {
150  for (unsigned int isel = 0; isel < m_selections.size(); ++isel) {
151  if (m_selections[isel].isSelected(mod->detId())) {
152  TH1F* tobefilled1d = nullptr;
153  TProfile* tobefilledprof = nullptr;
154  TH2F* tobefilled2d = nullptr;
155 
156  if (m_want1dHisto)
157  tobefilled1d = hist[isel];
158  if (m_wantProfile)
159  tobefilledprof = hprof[isel];
160  if (m_want2dHisto)
161  tobefilled2d = hist2d[isel];
162 
163  for (edmNew::DetSet<SiStripCluster>::const_iterator clus = mod->begin(); clus != mod->end(); clus++) {
164  for (unsigned int digi = 0; digi < clus->amplitudes().size(); ++digi) {
165  if (clus->amplitudes()[digi] > 0) {
166  unsigned int strip = clus->firstStrip() + digi;
167  if (m_folded)
168  strip = strip % 256;
169  if (tobefilled1d)
170  tobefilled1d->Fill(strip);
171  if (tobefilledprof)
172  tobefilledprof->Fill(strip, clus->amplitudes()[digi]);
173  if (tobefilled2d)
174  tobefilled2d->Fill(strip, clus->amplitudes()[digi]);
175  }
176  }
177  }
178  }
179  }
180  }
181 }
182 
183 #endif // DPGAnalysis_SiStripTools_DigiCollectionProfile_H
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< DetIdSelector > m_selections
selection
main part
Definition: corrVsCorr.py:100
data_type const * const_iterator
Definition: DetSetNew.h:31
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
bool isSelected(const std::vector< L1HPSPFTauQualityCut > &qualityCuts, const l1t::PFCandidate &pfCand, float_t primaryVertexZ)
void fill(edm::Handle< T > digis, const std::vector< TH1F *> &, const std::vector< TProfile *> &, const std::vector< TH2F *> &) const
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4