File indexing completed on 2024-04-06 12:01:15
0001 #ifndef UtilAlgos_HistoAnalyzer_h
0002 #define UtilAlgos_HistoAnalyzer_h
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0019 #include "CommonTools/UtilAlgos/interface/ExpressionHisto.h"
0020
0021 template <typename C>
0022 class HistoAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0023 public:
0024
0025 HistoAnalyzer(const edm::ParameterSet&);
0026
0027 ~HistoAnalyzer() override;
0028
0029 protected:
0030
0031 void analyze(const edm::Event&, const edm::EventSetup&) override;
0032
0033 private:
0034
0035 edm::EDGetTokenT<C> srcToken_;
0036
0037 bool usingWeights_;
0038
0039 edm::EDGetTokenT<double> weightsToken_;
0040
0041 std::vector<ExpressionHisto<typename C::value_type>*> vhistograms;
0042 };
0043
0044 template <typename C>
0045 HistoAnalyzer<C>::HistoAnalyzer(const edm::ParameterSet& par)
0046 : srcToken_(consumes<C>(par.template getParameter<edm::InputTag>("src"))),
0047 usingWeights_(par.exists("weights")),
0048 weightsToken_(
0049 mayConsume<double>(par.template getUntrackedParameter<edm::InputTag>("weights", edm::InputTag("fake")))) {
0050 usesResource(TFileService::kSharedResource);
0051 edm::Service<TFileService> fs;
0052 std::vector<edm::ParameterSet> histograms = par.template getParameter<std::vector<edm::ParameterSet> >("histograms");
0053 std::vector<edm::ParameterSet>::const_iterator it = histograms.begin();
0054 std::vector<edm::ParameterSet>::const_iterator end = histograms.end();
0055
0056
0057 for (; it != end; ++it) {
0058 ExpressionHisto<typename C::value_type>* hist = new ExpressionHisto<typename C::value_type>(*it);
0059 hist->initialize(fs->tFileDirectory());
0060 vhistograms.push_back(hist);
0061 }
0062 }
0063
0064 template <typename C>
0065 HistoAnalyzer<C>::~HistoAnalyzer() {
0066
0067 typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator it = vhistograms.begin();
0068 typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator end = vhistograms.end();
0069 for (; it != end; ++it) {
0070 (*it)->~ExpressionHisto<typename C::value_type>();
0071 }
0072 vhistograms.clear();
0073 }
0074
0075 template <typename C>
0076 void HistoAnalyzer<C>::analyze(const edm::Event& iEvent, const edm::EventSetup&) {
0077 edm::Handle<C> coll;
0078 iEvent.getByToken(srcToken_, coll);
0079 double weight = 1.0;
0080 if (usingWeights_) {
0081 edm::Handle<double> weightColl;
0082 iEvent.getByToken(weightsToken_, weightColl);
0083 weight = *weightColl;
0084 }
0085
0086 typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator it = vhistograms.begin();
0087 typename std::vector<ExpressionHisto<typename C::value_type>*>::iterator end = vhistograms.end();
0088
0089 for (; it != end; ++it) {
0090 uint32_t i = 0;
0091 for (typename C::const_iterator elem = coll->begin(); elem != coll->end(); ++elem, ++i) {
0092 if (!(*it)->fill(*elem, weight, i)) {
0093 break;
0094 }
0095 }
0096 }
0097 }
0098
0099 #endif