CMS 3D CMS Logo

EventShapeVarsProducer.cc
Go to the documentation of this file.
2 
4 
6 
7 #include <vector>
8 #include <memory>
9 
11 {
12  srcToken_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("src"));
13  r_ = cfg.exists("r") ? cfg.getParameter<double>("r") : 2.;
14  fwmax_ = cfg.exists("fwmax") ? cfg.getParameter<unsigned>("fwmax") : 0;
15 
16  produces<double>("thrust");
17  //produces<double>("oblateness");
18  produces<double>("isotropy");
19  produces<double>("circularity");
20  produces<double>("sphericity");
21  produces<double>("aplanarity");
22  produces<double>("C");
23  produces<double>("D");
24  if(fwmax_ > 0) produces<std::vector<double>>("FWmoments");
25 }
26 
27 void put(edm::Event& evt, double value, const char* instanceName)
28 {
29  evt.put(std::make_unique<double>(value), instanceName);
30 }
31 
33 {
34 
36  evt.getByToken(srcToken_, objects);
37 
38  Thrust thrustAlgo(objects->begin(), objects->end());
39  put(evt, thrustAlgo.thrust(), "thrust");
40  //put(evt, thrustAlgo.oblateness(), "oblateness");
41 
42  EventShapeVariables eventShapeVarsAlgo(*objects);
43  eventShapeVarsAlgo.set_r(r_);
44  put(evt, eventShapeVarsAlgo.isotropy(), "isotropy");
45  put(evt, eventShapeVarsAlgo.circularity(), "circularity");
46  put(evt, eventShapeVarsAlgo.sphericity(), "sphericity");
47  put(evt, eventShapeVarsAlgo.aplanarity(), "aplanarity");
48  put(evt, eventShapeVarsAlgo.C(), "C");
49  put(evt, eventShapeVarsAlgo.D(), "D");
50  if(fwmax_ > 0){
51  eventShapeVarsAlgo.setFWmax(fwmax_);
52  auto vfw = std::make_unique<std::vector<double>>(eventShapeVarsAlgo.getFWmoments());
53  evt.put(std::move(vfw),"FWmoments");
54  }
55 }
56 
58 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
EventShapeVarsProducer(const edm::ParameterSet &)
Class for the calculation of several event shape variables.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
void put(edm::Event &evt, double value, const char *instanceName)
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_
Definition: value.py:1
Definition: Thrust.h:38
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511