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 
21  public:
25 
26  void fill(edm::Handle<T> digis, const std::vector<TH1F*>&, const std::vector<TProfile*>&, const std::vector<TH2F*>&) const;
27 
28  private:
29 
30  bool m_folded;
34 
35  std::vector<DetIdSelector> m_selections;
36 
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 
52  std::vector<edm::ParameterSet> selconfigs = iConfig.getParameter<std::vector<edm::ParameterSet> >("selections");
53 
54  for(std::vector<edm::ParameterSet>::const_iterator selconfig=selconfigs.begin();selconfig!=selconfigs.end();++selconfig) {
55  DetIdSelector selection(*selconfig);
56  m_selections.push_back(selection);
57  }
58 
59 }
60 
61 template <class T>
62 void DigiCollectionProfiler<T>::fill(edm::Handle<T> digis, const std::vector<TH1F*>& hist, const std::vector<TProfile*>& hprof, const std::vector<TH2F*>& hist2d) const {
63 
64 }
65 
66 
67 template <>
68 void DigiCollectionProfiler<edm::DetSetVector<SiStripDigi> >::fill(edm::Handle<edm::DetSetVector<SiStripDigi> > digis, const std::vector<TH1F*>& hist, const std::vector<TProfile*>& hprof, const std::vector<TH2F*>& hist2d) const {
69 
70  for(edm::DetSetVector<SiStripDigi>::const_iterator mod = digis->begin();mod!=digis->end();mod++) {
71 
72  for(unsigned int isel=0;isel< m_selections.size(); ++isel) {
73 
74  if(m_selections[isel].isSelected(mod->detId())) {
75  TH1F* tobefilled1d=0;
76  TProfile* tobefilledprof=0;
77  TH2F* tobefilled2d=0;
78 
79  if(m_want1dHisto) tobefilled1d = hist[isel];
80  if(m_wantProfile) tobefilledprof = hprof[isel];
81  if(m_want2dHisto) tobefilled2d = hist2d[isel];
82 
83  for(edm::DetSet<SiStripDigi>::const_iterator digi=mod->begin();digi!=mod->end();digi++) {
84 
85  if(digi->adc()>0) {
86  unsigned int strip = digi->strip();
87  if(m_folded) strip = strip%256;
88  if(tobefilled1d) tobefilled1d->Fill(strip);
89  if(tobefilledprof) tobefilledprof->Fill(strip,digi->adc());
90  if(tobefilled2d) tobefilled2d->Fill(strip,digi->adc());
91  }
92  }
93  }
94  }
95  }
96 }
97 
98 template <>
99 void DigiCollectionProfiler<edm::DetSetVector<SiStripRawDigi> >::fill(edm::Handle<edm::DetSetVector<SiStripRawDigi> > digis, const std::vector<TH1F*>& hist, const std::vector<TProfile*>& hprof, const std::vector<TH2F*>& hist2d) const {
100 
101  for(edm::DetSetVector<SiStripRawDigi>::const_iterator mod = digis->begin();mod!=digis->end();mod++) {
102 
103  for(unsigned int isel=0;isel< m_selections.size(); ++isel) {
104 
105  if(m_selections[isel].isSelected(mod->detId())) {
106  TH1F* tobefilled1d=0;
107  TProfile* tobefilledprof=0;
108  TH2F* tobefilled2d=0;
109 
110  if(m_want1dHisto) tobefilled1d = hist[isel];
111  if(m_wantProfile) tobefilledprof = hprof[isel];
112  if(m_want2dHisto) tobefilled2d = hist2d[isel];
113 
114  unsigned int istrip=0;
115  for(edm::DetSet<SiStripRawDigi>::const_iterator digi=mod->begin();digi!=mod->end();digi++,++istrip) {
116 
117  if(digi->adc()>0) {
118  unsigned int strip = istrip;
119  if(m_folded) strip = strip%256;
120  if(tobefilled1d) tobefilled1d->Fill(strip);
121  if(tobefilledprof) tobefilledprof->Fill(strip,digi->adc());
122  if(tobefilled2d) tobefilled2d->Fill(strip,digi->adc());
123  }
124  }
125  }
126  }
127  }
128 }
129 
130 template <>
131 void DigiCollectionProfiler<edmNew::DetSetVector<SiStripCluster> >::fill(edm::Handle<edmNew::DetSetVector<SiStripCluster> > digis, const std::vector<TH1F*>& hist, const std::vector<TProfile*>& hprof, const std::vector<TH2F*>& hist2d) const {
132 
133  for(edmNew::DetSetVector<SiStripCluster>::const_iterator mod = digis->begin();mod!=digis->end();mod++) {
134 
135  for(unsigned int isel=0;isel< m_selections.size(); ++isel) {
136 
137  if(m_selections[isel].isSelected(mod->detId())) {
138  TH1F* tobefilled1d=0;
139  TProfile* tobefilledprof=0;
140  TH2F* tobefilled2d=0;
141 
142  if(m_want1dHisto) tobefilled1d = hist[isel];
143  if(m_wantProfile) tobefilledprof = hprof[isel];
144  if(m_want2dHisto) tobefilled2d = hist2d[isel];
145 
146  for(edmNew::DetSet<SiStripCluster>::const_iterator clus=mod->begin();clus!=mod->end();clus++) {
147 
148  for(unsigned int digi=0; digi < clus->amplitudes().size() ; ++digi) {
149 
150  if(clus->amplitudes()[digi]>0) {
151  unsigned int strip = clus->firstStrip()+digi;
152  if(m_folded) strip = strip%256;
153  if(tobefilled1d) tobefilled1d->Fill(strip);
154  if(tobefilledprof) tobefilledprof->Fill(strip,clus->amplitudes()[digi]);
155  if(tobefilled2d) tobefilled2d->Fill(strip,clus->amplitudes()[digi]);
156  }
157  }
158  }
159  }
160  }
161  }
162 }
163 
164 #endif // DPGAnalysis_SiStripTools_DigiCollectionProfile_H
T getParameter(std::string const &) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
std::vector< DetIdSelector > m_selections
selection
main part
Definition: corrVsCorr.py:98
data_type const * const_iterator
Definition: DetSetNew.h:30
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:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4