CMS 3D CMS Logo

OtherObjectVariableComputer.cc
Go to the documentation of this file.
1 //
2 
15 
19 
26 
27 template<typename T>
29  public:
30  explicit OtherObjectVariableComputer(const edm::ParameterSet & iConfig);
31  ~OtherObjectVariableComputer() override ;
32 
33  void produce(edm::Event & iEvent, const edm::EventSetup& iSetup) override;
34 
35  private:
39  double default_;
40  StringCutObjectSelector<T,true> objCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
41  bool doSort_;
43 };
44 
45 template<typename T>
47  probesToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("probes"))),
48  objectsToken_(consumes<edm::View<T> >(iConfig.getParameter<edm::InputTag>("objects"))),
49  objVar_(iConfig.getParameter<std::string>("expression")),
50  default_(iConfig.getParameter<double>("default")),
51  objCut_(iConfig.existsAs<std::string>("objectSelection") ? iConfig.getParameter<std::string>("objectSelection") : "", true),
52  doSort_(iConfig.existsAs<std::string>("objectSortDescendingBy")),
53  objSort_(doSort_ ? iConfig.getParameter<std::string>("objectSortDescendingBy") : "1", true)
54 {
55  produces<edm::ValueMap<float> >();
56 }
57 
58 
59 template<typename T>
61 {
62 }
63 
64 template<typename T>
65 void
67  using namespace edm;
68 
69  // read input
72  iEvent.getByToken(probesToken_, probes);
73  iEvent.getByToken(objectsToken_, objects);
74 
75  // fill
76  std::vector<std::pair<double, double> > selected;
77  typename View<T>::const_iterator object, endobjects = objects->end();
78  for (object = objects->begin(); object != endobjects; ++object) {
79  if (objCut_(*object)) {
80  selected.push_back(std::pair<double, double>(objSort_(*object), objVar_(*object)));
81  if (!doSort_) break; // if we take just the first one, there's no need of computing the others
82  }
83  }
84  if (doSort_ && selected.size() > 1) std::sort(selected.begin(), selected.end()); // sorts (ascending)
85 
86  // prepare vector for output
87  std::vector<float> values(probes->size(), (selected.empty() ? default_ : selected.back().second));
88 
89  // convert into ValueMap and store
90  auto valMap = std::make_unique<ValueMap<float>>();
92  filler.insert(probes, values.begin(), values.end());
93  filler.fill();
94  iEvent.put(std::move(valMap));
95 }
96 
97 
99 //typedef OtherObjectVariableComputer<reco::Track> OtherTrackVariableComputer;
100 //typedef OtherObjectVariableComputer<reco::Vertex> OtherVertexVariableComputer;
101 
104 //DEFINE_FWK_MODULE(OtherTrackVariableComputer);
105 //DEFINE_FWK_MODULE(OtherVertexVariableComputer);
StringCutObjectSelector< T, true > objCut_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
OtherObjectVariableComputer(const edm::ParameterSet &iConfig)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
Matcher of number of reconstructed objects in the event to probe.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
StringObjectFunction< T, true > objSort_
StringObjectFunction< T, true > objVar_
edm::EDGetTokenT< edm::View< reco::Candidate > > probesToken_
OtherObjectVariableComputer< reco::Candidate > OtherCandVariableComputer
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::View< T > > objectsToken_
int iEvent
Definition: GenABIO.cc:230
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
long double T
def move(src, dest)
Definition: eostools.py:511