CMS 3D CMS Logo

AssociationMapAnalyzer.cc
Go to the documentation of this file.
1 
14 
15 #include <vector>
16 
17 namespace edm {
18  class EventSetup;
19  class StreamID;
20 } // namespace edm
21 
22 namespace edmtest {
23 
25  public:
29  typedef edm::AssociationMap<edm::OneToManyWithQuality<std::vector<int>, std::vector<int>, double>>
32 
34  void analyze(edm::StreamID, edm::Event const& event, edm::EventSetup const&) const override;
35 
40 
49  };
50 
52  inputToken1_ = consumes<std::vector<int>>(pset.getParameter<edm::InputTag>("inputTag1"));
53  inputToken2_ = consumes<std::vector<int>>(pset.getParameter<edm::InputTag>("inputTag2"));
54  inputToken1V_ = consumes<edm::View<int>>(pset.getParameter<edm::InputTag>("inputTag1"));
55  inputToken2V_ = consumes<edm::View<int>>(pset.getParameter<edm::InputTag>("inputTag2"));
56  associationMapToken1_ = consumes<AssocOneToOne>(pset.getParameter<edm::InputTag>("associationMapTag1"));
57  associationMapToken2_ = consumes<AssocOneToOne>(pset.getParameter<edm::InputTag>("associationMapTag2"));
58  associationMapToken3_ = consumes<AssocOneToValue>(pset.getParameter<edm::InputTag>("associationMapTag3"));
59  associationMapToken4_ = consumes<AssocOneToValue>(pset.getParameter<edm::InputTag>("associationMapTag4"));
60  associationMapToken5_ = consumes<AssocOneToMany>(pset.getParameter<edm::InputTag>("associationMapTag5"));
61  associationMapToken6_ = consumes<AssocOneToManyWithQuality>(pset.getParameter<edm::InputTag>("associationMapTag6"));
62  associationMapToken7_ = consumes<AssocOneToOneView>(pset.getParameter<edm::InputTag>("associationMapTag7"));
63  associationMapToken8_ = consumes<AssocOneToOneView>(pset.getParameter<edm::InputTag>("associationMapTag8"));
64  }
65 
67  edm::Handle<std::vector<int>> inputCollection1 = event.getHandle(inputToken1_);
68 
69  edm::Handle<std::vector<int>> inputCollection2 = event.getHandle(inputToken2_);
70 
71  // Readout some entries from some AssociationMaps and check that
72  // we readout the same content as was was put in. We know the values
73  // by looking at the hard coded values in AssociationMapProducer.
74  // The particular values used are arbitrary and have no meaning.
75 
76  AssocOneToOne const& associationMap1 = event.get(associationMapToken1_);
77 
78  if (*associationMap1[edm::Ref<std::vector<int>>(inputCollection1, 0)] != 22 ||
79  *associationMap1[edm::Ref<std::vector<int>>(inputCollection1, 2)] != 24 ||
80  *associationMap1[edm::Ptr<int>(inputCollection1, 0)] != 22 ||
81  *associationMap1[edm::Ptr<int>(inputCollection1, 2)] != 24) {
82  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 1";
83  }
84  AssocOneToOne::const_iterator iter = associationMap1.begin();
85  if (*iter->val != 22 || iter->key != edm::Ref<std::vector<int>>(inputCollection1, 0)) {
86  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 2";
87  }
88  ++iter;
89  if (*iter->val != 24 || iter->key != edm::Ref<std::vector<int>>(inputCollection1, 2)) {
90  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 3";
91  }
92  ++iter;
93  if (iter != associationMap1.end()) {
94  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 4";
95  }
96 
97  // Case where handle arguments were used creating the AssociationMap
98  AssocOneToOne const& associationMap2 = event.get(associationMapToken2_);
99 
100  if (*associationMap2[edm::Ref<std::vector<int>>(inputCollection1, 0)] != 22 ||
101  *associationMap2[edm::Ref<std::vector<int>>(inputCollection1, 2)] != 25 ||
102  *associationMap2[edm::Ptr<int>(inputCollection1, 0)] != 22 ||
103  *associationMap2[edm::Ptr<int>(inputCollection1, 2)] != 25) {
104  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 5";
105  }
106 
107  AssocOneToOne::const_iterator iter2 = associationMap2.begin();
108  ++iter2;
109  if (*iter2->val != 25 || iter2->key != edm::Ref<std::vector<int>>(inputCollection1, 2)) {
110  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 6";
111  }
112 
113  // One to Value case
114 
115  AssocOneToValue const& associationMap3 = event.get(associationMapToken3_);
116 
117  if (associationMap3[edm::Ref<std::vector<int>>(inputCollection1, 0)] != 11.0 ||
118  associationMap3[edm::Ref<std::vector<int>>(inputCollection1, 2)] != 12.0 ||
119  associationMap3[edm::Ptr<int>(inputCollection1, 0)] != 11.0 ||
120  associationMap3[edm::Ptr<int>(inputCollection1, 2)] != 12.0) {
121  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 7";
122  }
123  AssocOneToValue::const_iterator iter3 = associationMap3.begin();
124  ++iter3;
125  if (iter3->val != 12.0 || iter3->key != edm::Ref<std::vector<int>>(inputCollection1, 2)) {
126  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 8";
127  }
128 
129  // One to Value case with handle argument
130 
131  AssocOneToValue const& associationMap4 = event.get(associationMapToken4_);
132 
133  if (associationMap4[edm::Ref<std::vector<int>>(inputCollection1, 0)] != 21.0 ||
134  associationMap4[edm::Ref<std::vector<int>>(inputCollection1, 2)] != 22.0 ||
135  associationMap4[edm::Ptr<int>(inputCollection1, 0)] != 21.0 ||
136  associationMap4[edm::Ptr<int>(inputCollection1, 2)] != 22.0) {
137  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 9";
138  }
139  AssocOneToValue::const_iterator iter4 = associationMap4.begin();
140  ++iter4;
141  if (iter4->val != 22.0 || iter4->key != edm::Ref<std::vector<int>>(inputCollection1, 2)) {
142  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 10";
143  }
144 
145  // One to Many
146 
147  AssocOneToMany const& associationMap5 = event.get(associationMapToken5_);
148 
149  if (*associationMap5[edm::Ref<std::vector<int>>(inputCollection1, 0)].at(0) != 22 ||
150  *associationMap5[edm::Ref<std::vector<int>>(inputCollection1, 2)].at(1) != 27 ||
151  *associationMap5[edm::Ptr<int>(inputCollection1, 0)].at(0) != 22 ||
152  *associationMap5[edm::Ptr<int>(inputCollection1, 2)].at(1) != 27) {
153  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 11";
154  }
155  AssocOneToMany::const_iterator iter5 = associationMap5.begin();
156  ++iter5;
157  if (*iter5->val.at(1) != 27 || iter5->key != edm::Ref<std::vector<int>>(inputCollection1, 2)) {
158  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 12";
159  }
160 
161  // One to Many With Quality
162 
163  AssocOneToManyWithQuality const& associationMap6 = event.get(associationMapToken6_);
164  if (*associationMap6[edm::Ref<std::vector<int>>(inputCollection1, 0)].at(0).first != 22 ||
165  *associationMap6[edm::Ref<std::vector<int>>(inputCollection1, 2)].at(1).first != 25 ||
166  *associationMap6[edm::Ptr<int>(inputCollection1, 0)].at(0).first != 22 ||
167  *associationMap6[edm::Ptr<int>(inputCollection1, 2)].at(1).first != 25 ||
168  associationMap6[edm::Ref<std::vector<int>>(inputCollection1, 0)].at(0).second != 31.0 ||
169  associationMap6[edm::Ref<std::vector<int>>(inputCollection1, 2)].at(1).second != 32.0 ||
170  associationMap6[edm::Ptr<int>(inputCollection1, 0)].at(0).second != 31.0 ||
171  associationMap6[edm::Ptr<int>(inputCollection1, 2)].at(1).second != 32.0) {
172  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 13";
173  }
174  AssocOneToManyWithQuality::const_iterator iter6 = associationMap6.begin();
175  ++iter6;
176  if (*iter6->val.at(1).first != 25 || iter6->val.at(1).second != 32.0 ||
177  iter6->key != edm::Ref<std::vector<int>>(inputCollection1, 2)) {
178  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 14";
179  }
180 
181  // One to One View
182 
183  edm::View<int> const& inputView1 = event.get(inputToken1V_);
184 
185  edm::Handle<edm::View<int>> inputView2 = event.getHandle(inputToken2V_);
186 
187  AssocOneToOneView const& associationMap7 = event.get(associationMapToken7_);
188  if (*associationMap7[inputView1.refAt(0)] != 24 || *associationMap7[inputView1.refAt(2)] != 25) {
189  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 15";
190  }
191  AssocOneToOneView::const_iterator iter7 = associationMap7.begin();
192  ++iter7;
193  if (*iter7->val != 25 || iter7->key != inputView1.refAt(2)) {
194  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 16";
195  }
196 
197  // One to One View built with 2 arguments constructor
198 
199  AssocOneToOneView const& associationMap8 = event.get(associationMapToken8_);
200  if (*associationMap8[inputView1.refAt(0)] != 26 || *associationMap8[inputView1.refAt(2)] != 27) {
201  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 17";
202  }
203  AssocOneToOneView::const_iterator iter8 = associationMap8.begin();
204  ++iter8;
205  if (*iter8->val != 27 || iter8->key != inputView1.refAt(2)) {
206  throw cms::Exception("TestFailure") << "unexpected result after using AssociationMap 18";
207  }
208  }
209 } // namespace edmtest
edm::EDGetTokenT< AssocOneToValue > associationMapToken3_
RefToBase< value_type > refAt(size_type i) const
edm::EDGetTokenT< AssocOneToOneView > associationMapToken8_
edm::AssociationMap< edm::OneToManyWithQuality< std::vector< int >, std::vector< int >, double > > AssocOneToManyWithQuality
edm::EDGetTokenT< AssocOneToOneView > associationMapToken7_
edm::EDGetTokenT< std::vector< int > > inputToken1_
const_iterator end() const
last iterator over the map (read only)
edm::AssociationMap< edm::OneToOne< std::vector< int >, std::vector< int > > > AssocOneToOne
U second(std::pair< T, U > const &p)
void analyze(edm::StreamID, edm::Event const &event, edm::EventSetup const &) const override
edm::EDGetTokenT< edm::View< int > > inputToken1V_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< AssocOneToValue > associationMapToken4_
edm::EDGetTokenT< AssocOneToOne > associationMapToken1_
edm::EDGetTokenT< std::vector< int > > inputToken2_
edm::EDGetTokenT< AssocOneToOne > associationMapToken2_
const value_type & get(size_type i) const
return value_typeelement with key i
AssociationMapAnalyzer(edm::ParameterSet const &)
const_iterator begin() const
first iterator over the map (read only)
HLT enums.
edm::AssociationMap< edm::OneToValue< std::vector< int >, double > > AssocOneToValue
edm::EDGetTokenT< edm::View< int > > inputToken2V_
edm::AssociationMap< edm::OneToMany< std::vector< int >, std::vector< int > > > AssocOneToMany
edm::AssociationMap< edm::OneToOne< edm::View< int >, edm::View< int > > > AssocOneToOneView
edm::EDGetTokenT< AssocOneToMany > associationMapToken5_
edm::EDGetTokenT< AssocOneToManyWithQuality > associationMapToken6_
Definition: event.py:1