CMS 3D CMS Logo

MultiToken.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_EgammaTools_MultiToken_H
2 #define RecoEgamma_EgammaTools_MultiToken_H
3 
10 
11 /*
12  * This class is a wrapper around a vector of EDM tokens, of which at least one
13  * is expected to yield a valid handle.
14  *
15  * The first time you call getValidHandle() or getHandle(), it will go over all
16  * tokens and try to get a valid handle. If no token yields a valid handle, it
17  * will either throw and exception or return the last non-valid handle.
18  *
19  * Once it found a valid handle, it will remember which token was used to get
20  * it and therefore doesn't need to loop over all tokens from there on.
21  *
22  * Example use case: auto-detection of AOD vs. MiniAOD.
23  *
24  * Created by Jonas Rembser on August 3, 2018.
25  */
26 
27 #include <memory>
28 
29 template <typename T>
30 class MultiTokenT {
31 
32  public:
33 
34  // Constructor which gets the input tags from a config to create the tokens
35  template <typename ... Tags>
36  MultiTokenT(edm::ConsumesCollector && cc, const edm::ParameterSet& pset, Tags && ... tags)
37  : isMaster_(true)
38  {
39  for (auto&& tag : { tags... }) {
40  tokens_.push_back(cc.mayConsume<T>(pset.getParameter<edm::InputTag>(tag)));
41  }
42  goodIndex_ = std::make_shared<int>(-1);
43  }
44 
45  // Constructor which gets the input tags from a config to create the tokens plus master token
46  template <typename S, typename ... Tags>
48  : isMaster_(false)
50  {
51  for (auto&& tag : { tags... }) {
52  tokens_.push_back(cc.mayConsume<T>(pset.getParameter<edm::InputTag>(tag)));
53  }
54  }
55 
56  // Get a handle on the event data, non-valid handle is allowed
58  {
60 
61  // If we already know which token works, take that one
62  if (*goodIndex_ >= 0) {
63  iEvent.getByToken(tokens_[*goodIndex_], handle);
64  return handle;
65  }
66 
67  if (!isMaster_) {
68  throw cms::Exception("MultiTokenTException") <<
69  "Trying to get a handle from a depending MultiToken before the master!";
70  }
71 
72  // If not, set the good token index parallel to getting the handle
73  handle = getInitialHandle(iEvent);
74 
75  if (*goodIndex_ == -1) {
76  *goodIndex_ = tokens_.size() - 1;
77  }
78  return handle;
79  }
80 
81  // Get a handle on the event data,
82  // throw exception if no token yields a valid handle
84  {
86 
87  // If we already know which token works, take that one
88  if (*goodIndex_ >= 0) {
89  iEvent.getByToken(tokens_[*goodIndex_], handle);
90  if (!handle.isValid())
91  throw cms::Exception("MultiTokenTException") <<
92  "Token gave valid handle in previously but not anymore!";
93  return handle;
94  }
95 
96  if (!isMaster_) {
97  throw cms::Exception("MultiTokenTException") <<
98  "Trying to get a handle from a depending MultiToken before the master!";
99  }
100 
101  // If not, set the good token index parallel to getting the handle
102  handle = getInitialHandle(iEvent);
103 
104  if (*goodIndex_ == -1) {
105  throw cms::Exception("MultiTokenTException") << "Neither handle is valid!";
106  }
107  return handle;
108  }
109 
110  // get the good token
112  {
113  // If we already know which token works, take that index
114  if (*goodIndex_ >= 0)
115  return tokens_[*goodIndex_];
116 
117  // If this is not a master MultiToken, just return what it got
118  if (!isMaster_) {
119  throw cms::Exception("MultiTokenTException") <<
120  "Trying to get a handle from a depending MultiToken before the master!";
121  }
122 
123  // Find which token is the good one by trying to get a handle
125  for (auto token:tokens_ ) {
126  iEvent.getByToken(token, handle);
127  if (handle.isValid()) {
128  return token;
129  }
130  }
131 
132  throw cms::Exception("MultiTokenTException") << "Neither token is valid!";
133  }
134 
135  int getGoodTokenIndex() const
136  {
137  return *goodIndex_;
138  }
139 
140  std::shared_ptr<int> getGoodTokenIndexPtr() const
141  {
142  return goodIndex_;
143  }
144 
145  private:
146 
148  {
149  // Try to retrieve the collection from the event. If we fail to
150  // retrieve the collection with one name, we next look for the one with
151  // the other name and so on.
153  for (size_t i = 0; i < tokens_.size(); ++i) {
154  iEvent.getByToken(tokens_[i], handle);
155  if (handle.isValid()) {
156  *goodIndex_ = i;
157  return handle;
158  }
159  }
160  return handle;
161  }
162 
163  const bool isMaster_;
164  std::vector<edm::EDGetTokenT<T>> tokens_;
165  std::shared_ptr<int> goodIndex_;
166 };
167 
168 #endif
T getParameter(std::string const &) const
Master< F > master(const F &f)
Definition: FunctClone.h:68
edm::Handle< T > getInitialHandle(const edm::Event &iEvent)
Definition: MultiToken.h:147
edm::Handle< T > getHandle(const edm::Event &iEvent)
Definition: MultiToken.h:57
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
MultiTokenT(MultiTokenT< S > &master, edm::ConsumesCollector &&cc, const edm::ParameterSet &pset, Tags &&...tags)
Definition: MultiToken.h:47
std::shared_ptr< int > getGoodTokenIndexPtr() const
Definition: MultiToken.h:140
std::vector< edm::EDGetTokenT< T > > tokens_
Definition: MultiToken.h:164
int iEvent
Definition: GenABIO.cc:230
const bool isMaster_
Definition: MultiToken.h:163
bool isValid() const
Definition: HandleBase.h:74
edm::Handle< T > getValidHandle(const edm::Event &iEvent)
Definition: MultiToken.h:83
MultiTokenT(edm::ConsumesCollector &&cc, const edm::ParameterSet &pset, Tags &&...tags)
Definition: MultiToken.h:36
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
std::shared_ptr< int > goodIndex_
Definition: MultiToken.h:165
int getGoodTokenIndex() const
Definition: MultiToken.h:135
long double T