ATLAS Offline Software
OriginalAodCounts.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
9 
10 #include "xAODRootAccess/TEvent.h"
11 
12 #include "TFile.h"
13 
15  // add default values here
16 {
17 }
18 
20  nEventsProcessed(0),
21  sumOfWeights(0),
22  sumOfWeightsSquared(0),
23  nIncomplete(0)
24 {
25 }
26 
28  nEventsProcessed += a.nEventsProcessed;
29  sumOfWeights += a.sumOfWeights;
30  sumOfWeightsSquared += a.sumOfWeightsSquared;
31  nIncomplete += a.nIncomplete;
32  return *this;
33 }
34 
35 
37  OriginalAodCounts ret = a;
38  ret += b;
39  return ret;
40 }
41 
42 // this is copied from
43 // https://twiki.cern.ch/twiki/bin/view/AtlasProtected/AnalysisMetadata#Analysis_Metadata_Root_or_athena
44 
45 
47  const AodCountsConfig&) {
48 
49  OriginalAodCounts counts;
50 
51  // check for corruption
52  const xAOD::CutBookkeeperContainer* incompleteCBC = nullptr;
53  if(!event.retrieveMetaInput(incompleteCBC, "IncompleteCutBookkeepers").isSuccess()){
54  throw std::runtime_error("Failed to retrieve IncompleteCutBookkeepers from MetaData!");
55  }
56  counts.nIncomplete = incompleteCBC->size();
57 
58  // Now, let's find the actual information
59  const xAOD::CutBookkeeperContainer* completeCBC = 0;
60  if(!event.retrieveMetaInput(completeCBC, "CutBookkeepers").isSuccess()){
61  throw std::runtime_error("Failed to retrieve CutBookkeepers from MetaData!");
62  }
63 
64  // Now, let's actually find the right one that contains all the needed info...
65  const xAOD::CutBookkeeper* allEventsCBK = 0;
66  int maxCycle = -1;
67  for (const auto *cbk: *completeCBC) {
68  if (cbk->cycle() > maxCycle &&
69  cbk->name() == "AllExecutedEvents" &&
70  cbk->inputStream() == "StreamAOD") {
71  allEventsCBK = cbk;
72  maxCycle = cbk->cycle();
73  }
74  }
75  counts.nEventsProcessed = allEventsCBK->nAcceptedEvents();
76  counts.sumOfWeights = allEventsCBK->sumOfEventWeights();
77  counts.sumOfWeightsSquared = allEventsCBK->sumOfEventWeightsSquared();
78  return counts;
79 }
80 
POOL::TEvent::retrieveMetaInput
StatusCode retrieveMetaInput(const T *&obj, const std::string &key)
Definition: PhysicsAnalysis/POOLRootAccess/POOLRootAccess/TEvent.h:90
OriginalAodCounts::nEventsProcessed
unsigned long long nEventsProcessed
Definition: OriginalAodCounts.h:21
CutBookkeeper.h
xAOD::CutBookkeeper_v1::cycle
int cycle() const
Get the skimming cycle that this CutBookkeeper was running in.
Definition: CutBookkeeper_v1.cxx:205
operator+
OriginalAodCounts operator+(const OriginalAodCounts &a, const OriginalAodCounts &b)
Definition: OriginalAodCounts.cxx:36
OriginalAodCounts::nIncomplete
int nIncomplete
Definition: OriginalAodCounts.h:24
xAOD::CutBookkeeper_v1
Description of the class that is used to keep track of event counts.
Definition: CutBookkeeper_v1.h:29
xAOD::CutBookkeeper_v1::sumOfEventWeightsSquared
double sumOfEventWeightsSquared() const
Get the sum-of-(event-weights-squared) that this CutBookkeeper has seen.
Definition: CutBookkeeper_v1.cxx:327
OriginalAodCounts.h
getOriginalAodCounts
OriginalAodCounts getOriginalAodCounts(xAOD::TEvent &event, const AodCountsConfig &)
Definition: OriginalAodCounts.cxx:46
AodCountsConfig
Definition: OriginalAodCounts.h:14
OriginalAodCounts::sumOfWeightsSquared
double sumOfWeightsSquared
Definition: OriginalAodCounts.h:23
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
TEvent.h
CutBookkeeperContainer.h
AodCountsConfig::AodCountsConfig
AodCountsConfig()
Definition: OriginalAodCounts.cxx:14
xAOD::CutBookkeeperContainer_v1
Container that holds the Container of all CutBookkeepers.
Definition: CutBookkeeperContainer_v1.h:27
OriginalAodCounts::OriginalAodCounts
OriginalAodCounts()
Definition: OriginalAodCounts.cxx:19
xAOD::CutBookkeeper_v1::nAcceptedEvents
uint64_t nAcceptedEvents() const
Get the number of accepted events that this CutBookkeeper has seen.
Definition: CutBookkeeper_v1.cxx:291
plotBeamSpotMon.b
b
Definition: plotBeamSpotMon.py:77
OriginalAodCounts::operator+=
OriginalAodCounts & operator+=(const OriginalAodCounts &a)
Definition: OriginalAodCounts.cxx:27
xAOD::CutBookkeeper_v1::sumOfEventWeights
double sumOfEventWeights() const
Get the sum-of-event-weights that this CutBookkeeper has seen.
Definition: CutBookkeeper_v1.cxx:309
a
TList * a
Definition: liststreamerinfos.cxx:10
OriginalAodCounts::sumOfWeights
double sumOfWeights
Definition: OriginalAodCounts.h:22
OriginalAodCounts
Definition: OriginalAodCounts.h:19
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
xAOD::TEvent
Tool for accessing xAOD files outside of Athena.
Definition: Control/xAODRootAccess/xAODRootAccess/TEvent.h:84