CMS 3D CMS Logo

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