HistoAnalyzer

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99
#ifndef UtilAlgos_HistoAnalyzer_h
#define UtilAlgos_HistoAnalyzer_h
/** \class HistoAnalyzer
 *
 * Creates histograms defined in config file
 *
 * \author: Benedikt Hegner, DESY
 *
 * Template parameters:
 * - C : Concrete candidate collection type
 *
 */
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "CommonTools/UtilAlgos/interface/TFileService.h"
#include "CommonTools/UtilAlgos/interface/ExpressionHisto.h"

template <typename C>
class HistoAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
public:
  /// constructor from parameter set
  HistoAnalyzer(const edm::ParameterSet&);
  /// destructor
  ~HistoAnalyzer() override;

protected:
  /// process an event
  void analyze(const edm::Event&, const edm::EventSetup&) override;

private:
  /// label of the collection to be read in
  edm::EDGetTokenT<C> srcToken_;
  /// Do we weight events?
  bool usingWeights_;
  /// label of the weight collection (can be null for weights = 1)
  edm::EDGetTokenT<double> weightsToken_;
  /// vector of the histograms
  std::vector<ExpressionHisto<typename C::value_type>*> vhistograms;
};

template <typename C>
HistoAnalyzer<C>::HistoAnalyzer(const edm::ParameterSet& par)
    : srcToken_(consumes<C>(par.template getParameter<edm::InputTag>("src"))),
      usingWeights_(par.exists("weights")),
      weightsToken_(
          mayConsume<double>(par.template getUntrackedParameter<edm::InputTag>("weights", edm::InputTag("fake")))) {
  usesResource(TFileService::kSharedResource);
  edm::Service<TFileService> fs;
  std::vector<edm::ParameterSet> histograms = par.template getParameter<std::vector<edm::ParameterSet> >("histograms");
  std::vector<edm::ParameterSet>::const_iterator it = histograms.begin();
  std::vector<edm::ParameterSet>::const_iterator end = histograms.end();

  // create the histograms from the given parameter sets
  for (; it != end; ++it) {
    ExpressionHisto<typename C::value_type>* hist = new ExpressionHisto<typename C::value_type>(*it);
    hist->initialize(fs->tFileDirectory());
    vhistograms.push_back(hist);
  }
}

template <typename C>
HistoAnalyzer<C>::~HistoAnalyzer() {
  // delete all histograms and clear the vector of pointers
  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator it = vhistograms.begin();
  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator end = vhistograms.end();
  for (; it != end; ++it) {
    (*it)->~ExpressionHisto<typename C::value_type>();
  }
  vhistograms.clear();
}

template <typename C>
void HistoAnalyzer<C>::analyze(const edm::Event& iEvent, const edm::EventSetup&) {
  edm::Handle<C> coll;
  iEvent.getByToken(srcToken_, coll);
  double weight = 1.0;
  if (usingWeights_) {
    edm::Handle<double> weightColl;
    iEvent.getByToken(weightsToken_, weightColl);
    weight = *weightColl;
  }

  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator it = vhistograms.begin();
  typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator end = vhistograms.end();

  for (; it != end; ++it) {
    uint32_t i = 0;
    for (typename C::const_iterator elem = coll->begin(); elem != coll->end(); ++elem, ++i) {
      if (!(*it)->fill(*elem, weight, i)) {
        break;
      }
    }
  }
}

#endif