CMS 3D CMS Logo

HistoAnalyzer.h
Go to the documentation of this file.
1 #ifndef UtilAlgos_HistoAnalyzer_h
2 #define UtilAlgos_HistoAnalyzer_h
3 
20 
21 template <typename C>
23 public:
27  ~HistoAnalyzer() override;
28 
29 protected:
31  void analyze(const edm::Event&, const edm::EventSetup&) override;
32 
33 private:
41  std::vector<ExpressionHisto<typename C::value_type>*> vhistograms;
42 };
43 
44 template <typename C>
46  : srcToken_(consumes<C>(par.template getParameter<edm::InputTag>("src"))),
47  usingWeights_(par.exists("weights")),
48  weightsToken_(
49  mayConsume<double>(par.template getUntrackedParameter<edm::InputTag>("weights", edm::InputTag("fake")))) {
51  std::vector<edm::ParameterSet> histograms = par.template getParameter<std::vector<edm::ParameterSet> >("histograms");
52  std::vector<edm::ParameterSet>::const_iterator it = histograms.begin();
53  std::vector<edm::ParameterSet>::const_iterator end = histograms.end();
54 
55  // create the histograms from the given parameter sets
56  for (; it != end; ++it) {
58  hist->initialize(fs->tFileDirectory());
59  vhistograms.push_back(hist);
60  }
61 }
62 
63 template <typename C>
65  // delete all histograms and clear the vector of pointers
66  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator it = vhistograms.begin();
67  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator end = vhistograms.end();
68  for (; it != end; ++it) {
69  (*it)->~ExpressionHisto<typename C::value_type>();
70  }
71  vhistograms.clear();
72 }
73 
74 template <typename C>
76  edm::Handle<C> coll;
77  iEvent.getByToken(srcToken_, coll);
78  double weight = 1.0;
79  if (usingWeights_) {
80  edm::Handle<double> weightColl;
81  iEvent.getByToken(weightsToken_, weightColl);
82  weight = *weightColl;
83  }
84 
85  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator it = vhistograms.begin();
86  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator end = vhistograms.end();
87 
88  for (; it != end; ++it) {
89  uint32_t i = 0;
90  for (typename C::const_iterator elem = coll->begin(); elem != coll->end(); ++elem, ++i) {
91  if (!(*it)->fill(*elem, weight, i)) {
92  break;
93  }
94  }
95  }
96 }
97 
98 #endif
edm::EDGetTokenT< C > srcToken_
label of the collection to be read in
Definition: HistoAnalyzer.h:35
Definition: weight.py:1
std::vector< ExpressionHisto< typename C::value_type > * > vhistograms
vector of the histograms
Definition: HistoAnalyzer.h:41
int iEvent
Definition: GenABIO.cc:224
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
__shared__ Hist hist
~HistoAnalyzer() override
destructor
Definition: HistoAnalyzer.h:64
HLT enums.
edm::EDGetTokenT< double > weightsToken_
label of the weight collection (can be null for weights = 1)
Definition: HistoAnalyzer.h:39
void analyze(const edm::Event &, const edm::EventSetup &) override
process an event
Definition: HistoAnalyzer.h:75
HistoAnalyzer(const edm::ParameterSet &)
constructor from parameter set
Definition: HistoAnalyzer.h:45
bool usingWeights_
Do we weight events?
Definition: HistoAnalyzer.h:37