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>
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 <>
68  const std::vector<TH1F*>& hist,
69  const std::vector<TProfile*>& hprof,
70  const std::vector<TH2F*>& hist2d) const {
71  for (edm::DetSetVector<SiStripDigi>::const_iterator mod = digis->begin(); mod != digis->end(); mod++) {
72  for (unsigned int isel = 0; isel < m_selections.size(); ++isel) {
73  if (m_selections[isel].isSelected(mod->detId())) {
74  TH1F* tobefilled1d = nullptr;
75  TProfile* tobefilledprof = nullptr;
76  TH2F* tobefilled2d = nullptr;
77 
78  if (m_want1dHisto)
79  tobefilled1d = hist[isel];
80  if (m_wantProfile)
81  tobefilledprof = hprof[isel];
82  if (m_want2dHisto)
83  tobefilled2d = hist2d[isel];
84 
85  for (edm::DetSet<SiStripDigi>::const_iterator digi = mod->begin(); digi != mod->end(); digi++) {
86  if (digi->adc() > 0) {
87  unsigned int strip = digi->strip();
88  if (m_folded)
89  strip = strip % 256;
90  if (tobefilled1d)
91  tobefilled1d->Fill(strip);
92  if (tobefilledprof)
93  tobefilledprof->Fill(strip, digi->adc());
94  if (tobefilled2d)
95  tobefilled2d->Fill(strip, digi->adc());
96  }
97  }
98  }
99  }
100  }
101 }
102 
103 template <>
106  const std::vector<TH1F*>& hist,
107  const std::vector<TProfile*>& hprof,
108  const std::vector<TH2F*>& hist2d) const {
109  for (edm::DetSetVector<SiStripRawDigi>::const_iterator mod = digis->begin(); mod != digis->end(); mod++) {
110  for (unsigned int isel = 0; isel < m_selections.size(); ++isel) {
111  if (m_selections[isel].isSelected(mod->detId())) {
112  TH1F* tobefilled1d = nullptr;
113  TProfile* tobefilledprof = nullptr;
114  TH2F* tobefilled2d = nullptr;
115 
116  if (m_want1dHisto)
117  tobefilled1d = hist[isel];
118  if (m_wantProfile)
119  tobefilledprof = hprof[isel];
120  if (m_want2dHisto)
121  tobefilled2d = hist2d[isel];
122 
123  unsigned int istrip = 0;
124  for (edm::DetSet<SiStripRawDigi>::const_iterator digi = mod->begin(); digi != mod->end(); digi++, ++istrip) {
125  if (digi->adc() > 0) {
126  unsigned int strip = istrip;
127  if (m_folded)
128  strip = strip % 256;
129  if (tobefilled1d)
130  tobefilled1d->Fill(strip);
131  if (tobefilledprof)
132  tobefilledprof->Fill(strip, digi->adc());
133  if (tobefilled2d)
134  tobefilled2d->Fill(strip, digi->adc());
135  }
136  }
137  }
138  }
139  }
140 }
141 
142 template <>
145  const std::vector<TH1F*>& hist,
146  const std::vector<TProfile*>& hprof,
147  const std::vector<TH2F*>& hist2d) const {
148  for (edmNew::DetSetVector<SiStripCluster>::const_iterator mod = digis->begin(); mod != digis->end(); mod++) {
149  for (unsigned int isel = 0; isel < m_selections.size(); ++isel) {
150  if (m_selections[isel].isSelected(mod->detId())) {
151  TH1F* tobefilled1d = nullptr;
152  TProfile* tobefilledprof = nullptr;
153  TH2F* tobefilled2d = nullptr;
154 
155  if (m_want1dHisto)
156  tobefilled1d = hist[isel];
157  if (m_wantProfile)
158  tobefilledprof = hprof[isel];
159  if (m_want2dHisto)
160  tobefilled2d = hist2d[isel];
161 
162  for (edmNew::DetSet<SiStripCluster>::const_iterator clus = mod->begin(); clus != mod->end(); clus++) {
163  for (unsigned int digi = 0; digi < clus->amplitudes().size(); ++digi) {
164  if (clus->amplitudes()[digi] > 0) {
165  unsigned int strip = clus->firstStrip() + digi;
166  if (m_folded)
167  strip = strip % 256;
168  if (tobefilled1d)
169  tobefilled1d->Fill(strip);
170  if (tobefilledprof)
171  tobefilledprof->Fill(strip, clus->amplitudes()[digi]);
172  if (tobefilled2d)
173  tobefilled2d->Fill(strip, clus->amplitudes()[digi]);
174  }
175  }
176  }
177  }
178  }
179  }
180 }
181 
182 #endif // DPGAnalysis_SiStripTools_DigiCollectionProfile_H
T getParameter(std::string const &) const
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
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:32
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4