CMS 3D CMS Logo

SimpleSystematicsAnalyzer.cc
Go to the documentation of this file.
5 
7 public:
10  bool filter(edm::Event&, const edm::EventSetup&) override;
11  void beginJob() override;
12  void endJob() override;
13 
14 private:
16  std::vector<edm::InputTag> weightTags_;
17  std::vector<edm::EDGetTokenT<double> > weightTokens_;
19  unsigned int originalEvents_;
20  std::vector<double> weightedEvents_;
21  unsigned int selectedEvents_;
22  std::vector<double> weightedSelectedEvents_;
23  std::vector<double> weighted2SelectedEvents_;
24 };
25 
33 
35 
37 
40  : selectorPath_(pset.getUntrackedParameter<std::string>("SelectorPath", "")),
41  weightTags_(pset.getUntrackedParameter<std::vector<edm::InputTag> >("WeightTags")),
42  weightTokens_(
43  edm::vector_transform(weightTags_, [this](edm::InputTag const& tag) { return consumes<double>(tag); })),
44  triggerResultsToken_(consumes<edm::TriggerResults>(edm::InputTag("TriggerResults"))) {}
45 
48 
51  originalEvents_ = 0;
52  selectedEvents_ = 0;
53  edm::LogVerbatim("SimpleSystematicsAnalysis") << "Uncertainties will be determined for the following tags: ";
54  for (unsigned int i = 0; i < weightTags_.size(); ++i) {
55  edm::LogVerbatim("SimpleSystematicsAnalysis") << "\t" << weightTags_[i].encode();
56  weightedEvents_.push_back(0.);
57  weightedSelectedEvents_.push_back(0.);
58  weighted2SelectedEvents_.push_back(0.);
59  }
60 }
61 
64  if (originalEvents_ == 0) {
65  edm::LogVerbatim("SimpleSystematicsAnalysis") << "NO EVENTS => NO RESULTS";
66  return;
67  }
68  if (selectedEvents_ == 0) {
69  edm::LogVerbatim("SimpleSystematicsAnalysis") << "NO SELECTED EVENTS => NO RESULTS";
70  return;
71  }
72 
73  edm::LogVerbatim("SimpleSystematicsAnalysis") << "\n>>>> Begin of Weight systematics summary >>>>";
74  edm::LogVerbatim("SimpleSystematicsAnalysis") << "Total number of analyzed data: " << originalEvents_ << " [events]";
75  double originalAcceptance = double(selectedEvents_) / originalEvents_;
76  edm::LogVerbatim("SimpleSystematicsAnalysis")
77  << "Total number of selected data: " << selectedEvents_ << " [events], corresponding to acceptance: ["
78  << originalAcceptance * 100 << " +- "
79  << 100 * sqrt(originalAcceptance * (1. - originalAcceptance) / originalEvents_) << "] %";
80 
81  for (unsigned int i = 0; i < weightTags_.size(); ++i) {
82  edm::LogVerbatim("SimpleSystematicsAnalysis") << "Results for Weight Tag: " << weightTags_[i].encode() << " ---->";
83 
84  double acc_central = 0.;
85  double acc2_central = 0.;
86  if (weightedEvents_[i] > 0) {
87  acc_central = weightedSelectedEvents_[i] / weightedEvents_[i];
88  acc2_central = weighted2SelectedEvents_[i] / weightedEvents_[i];
89  }
90  double waverage = weightedEvents_[i] / originalEvents_;
91  edm::LogVerbatim("SimpleSystematicsAnalysis")
92  << "\tTotal Events after reweighting: " << weightedEvents_[i] << " [events]";
93  edm::LogVerbatim("SimpleSystematicsAnalysis")
94  << "\tEvents selected after reweighting: " << weightedSelectedEvents_[i] << " [events]";
95  edm::LogVerbatim("SimpleSystematicsAnalysis")
96  << "\tAcceptance after reweighting: [" << acc_central * 100 << " +- "
97  << 100 * sqrt((acc2_central / waverage - acc_central * acc_central) / originalEvents_) << "] %";
98  double xi = acc_central - originalAcceptance;
99  double deltaxi = (acc2_central - (originalAcceptance + 2 * xi + xi * xi)) / originalEvents_;
100  if (deltaxi > 0)
101  deltaxi = sqrt(deltaxi);
102  else
103  deltaxi = 0.;
104  edm::LogVerbatim("SimpleSystematicsAnalysis")
105  << "\ti.e. [" << std::setprecision(4) << 100 * xi / originalAcceptance << " +- " << std::setprecision(4)
106  << 100 * deltaxi / originalAcceptance << "] % relative variation with respect to the original acceptance";
107  }
108  edm::LogVerbatim("SimpleSystematicsAnalysis") << ">>>> End of Weight systematics summary >>>>";
109 }
110 
113  originalEvents_++;
114 
115  bool selectedEvent = false;
117  if (!ev.getByToken(triggerResultsToken_, triggerResults)) {
118  edm::LogError("SimpleSystematicsAnalysis") << ">>> TRIGGER collection does not exist !!!";
119  return false;
120  }
121 
122  const edm::TriggerNames& trigNames = ev.triggerNames(*triggerResults);
123  unsigned int pathIndex = trigNames.triggerIndex(selectorPath_);
124  bool pathFound = (pathIndex < trigNames.size()); // pathIndex >= 0, since pathIndex is unsigned
125  if (pathFound) {
126  if (triggerResults->accept(pathIndex))
127  selectedEvent = true;
128  }
129  //edm::LogVerbatim("SimpleSystematicsAnalysis") << ">>>> Path Name: " << selectorPath_ << ", selected? " << selectedEvent;
130 
131  if (selectedEvent)
132  selectedEvents_++;
133 
134  for (unsigned int i = 0; i < weightTags_.size(); ++i) {
135  edm::Handle<double> weightHandle;
136  ev.getByToken(weightTokens_[i], weightHandle);
137  weightedEvents_[i] += (*weightHandle);
138  if (selectedEvent) {
139  weightedSelectedEvents_[i] += (*weightHandle);
140  weighted2SelectedEvents_[i] += pow((*weightHandle), 2);
141  }
142  }
143 
144  return true;
145 }
146 
SimpleSystematicsAnalyzer::weightedEvents_
std::vector< double > weightedEvents_
Definition: SimpleSystematicsAnalyzer.cc:20
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
trigNames
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:57
ESHandle.h
TriggerResults.h
edm::EDGetTokenT< edm::TriggerResults >
edm
HLT enums.
Definition: AlignableModifier.h:19
hybridSuperClusters_cfi.xi
xi
Definition: hybridSuperClusters_cfi.py:10
SimpleSystematicsAnalyzer::filter
bool filter(edm::Event &, const edm::EventSetup &) override
Definition: SimpleSystematicsAnalyzer.cc:112
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
triggerResults
static const std::string triggerResults
Definition: EdmProvDump.cc:45
EDFilter.h
SimpleSystematicsAnalyzer::~SimpleSystematicsAnalyzer
~SimpleSystematicsAnalyzer() override
Definition: SimpleSystematicsAnalyzer.cc:47
watchdog.const
const
Definition: watchdog.py:83
edm::Handle< edm::TriggerResults >
SimpleSystematicsAnalyzer::selectedEvents_
unsigned int selectedEvents_
Definition: SimpleSystematicsAnalyzer.cc:21
SimpleSystematicsAnalyzer::weightTokens_
std::vector< edm::EDGetTokenT< double > > weightTokens_
Definition: SimpleSystematicsAnalyzer.cc:17
MakerMacros.h
SimpleSystematicsAnalyzer::endJob
void endJob() override
Definition: SimpleSystematicsAnalyzer.cc:63
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
GlobalPosition_Frontier_DevDB_cff.tag
tag
Definition: GlobalPosition_Frontier_DevDB_cff.py:11
SimpleSystematicsAnalyzer::weightedSelectedEvents_
std::vector< double > weightedSelectedEvents_
Definition: SimpleSystematicsAnalyzer.cc:22
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
SimpleSystematicsAnalyzer
Definition: SimpleSystematicsAnalyzer.cc:6
SimpleSystematicsAnalyzer::weighted2SelectedEvents_
std::vector< double > weighted2SelectedEvents_
Definition: SimpleSystematicsAnalyzer.cc:23
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::vector_transform
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
edm::ParameterSet
Definition: ParameterSet.h:47
SimpleSystematicsAnalyzer::SimpleSystematicsAnalyzer
SimpleSystematicsAnalyzer(const edm::ParameterSet &pset)
Definition: SimpleSystematicsAnalyzer.cc:39
Event.h
TriggerNames.h
SimpleSystematicsAnalyzer::weightTags_
std::vector< edm::InputTag > weightTags_
Definition: SimpleSystematicsAnalyzer.cc:16
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::EDFilter
Definition: EDFilter.h:38
edm::EventSetup
Definition: EventSetup.h:57
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
InputTag.h
SimpleSystematicsAnalyzer::triggerResultsToken_
edm::EDGetTokenT< edm::TriggerResults > triggerResultsToken_
Definition: SimpleSystematicsAnalyzer.cc:18
SimpleSystematicsAnalyzer::beginJob
void beginJob() override
Definition: SimpleSystematicsAnalyzer.cc:50
std
Definition: JetResolutionObject.h:76
Frameworkfwd.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
transform.h
SimpleSystematicsAnalyzer::selectorPath_
std::string selectorPath_
Definition: SimpleSystematicsAnalyzer.cc:15
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
edm::TriggerNames
Definition: TriggerNames.h:55
EventSetup.h
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
SimpleSystematicsAnalyzer::originalEvents_
unsigned int originalEvents_
Definition: SimpleSystematicsAnalyzer.cc:19
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27