CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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:
28 
29 protected:
31  virtual 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_(mayConsume<double>(par.template getUntrackedParameter<edm::InputTag>("weights", edm::InputTag("fake"))))
49 {
51  std::vector<edm::ParameterSet> histograms =
52  par.template getParameter<std::vector<edm::ParameterSet> >("histograms");
53  std::vector<edm::ParameterSet>::const_iterator it = histograms.begin();
54  std::vector<edm::ParameterSet>::const_iterator end = histograms.end();
55 
56  // create the histograms from the given parameter sets
57  for (; it!=end; ++it)
58  {
60  hist->initialize(fs->tFileDirectory());
61  vhistograms.push_back(hist);
62  }
63 
64 }
65 
66 template<typename C>
68 {
69  // delete all histograms and clear the vector of pointers
70  typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator it = vhistograms.begin();
71  typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator end = vhistograms.end();
72  for (;it!=end; ++it){
73  (*it)->~ExpressionHisto<typename C::value_type>();
74  }
75  vhistograms.clear();
76 }
77 
78 template<typename C>
80 {
82  iEvent.getByToken( srcToken_, coll);
83  double weight = 1.0;
84  if (usingWeights_) {
85  edm::Handle<double> weightColl;
86  iEvent.getByToken( weightsToken_, weightColl );
87  weight = *weightColl;
88  }
89 
90  typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator it = vhistograms.begin();
91  typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator end = vhistograms.end();
92 
93  for (;it!=end; ++it){
94  uint32_t i = 0;
95  for( typename C::const_iterator elem=coll->begin(); elem!=coll->end(); ++elem, ++i ) {
96  if (!(*it)->fill( *elem, weight, i )) {
97  break;
98  }
99  }
100  }
101 }
102 
103 #endif
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< C > srcToken_
label of the collection to be read in
Definition: HistoAnalyzer.h:35
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
std::vector< ExpressionHisto< typename C::value_type > * > vhistograms
vector of the histograms
Definition: HistoAnalyzer.h:41
TFileDirectory & tFileDirectory()
Definition: TFileService.h:42
int iEvent
Definition: GenABIO.cc:230
#define end
Definition: vmac.h:37
JetCorrectorParametersCollection coll
Definition: classes.h:10
~HistoAnalyzer()
destructor
Definition: HistoAnalyzer.h:67
edm::EDGetTokenT< double > weightsToken_
label of the weight collection (can be null for weights = 1)
Definition: HistoAnalyzer.h:39
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
process an event
Definition: HistoAnalyzer.h:79
void initialize(TFileDirectory &fs)
int weight
Definition: histoStyle.py:50
HistoAnalyzer(const edm::ParameterSet &)
constructor from parameter set
Definition: HistoAnalyzer.h:45
bool usingWeights_
Do we weight events?
Definition: HistoAnalyzer.h:37
def template
Definition: svgfig.py:520
edm::Service< TFileService > fs