ATLAS Offline Software
Loading...
Searching...
No Matches
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
11
12#include "TFile.h"
13
15 // add default values here
16{
17}
18
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
static Double_t a
OriginalAodCounts getOriginalAodCounts(xAOD::TEvent &event, const AodCountsConfig &)
OriginalAodCounts operator+(const OriginalAodCounts &a, const OriginalAodCounts &b)
size_type size() const noexcept
Returns the number of elements in the collection.
double sumOfEventWeightsSquared() const
Get the sum-of-(event-weights-squared) that this CutBookkeeper has seen.
double sumOfEventWeights() const
Get the sum-of-event-weights that this CutBookkeeper has seen.
uint64_t nAcceptedEvents() const
Get the number of accepted events that this CutBookkeeper has seen.
Tool for accessing xAOD files outside of Athena.
CutBookkeeper_v1 CutBookkeeper
Define the latest version of the CutBookkeeper class.
CutBookkeeperContainer_v1 CutBookkeeperContainer
Define the latest version of the CutBookkeeperContainer class.
OriginalAodCounts & operator+=(const OriginalAodCounts &a)
unsigned long long nEventsProcessed