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 
8 {
9  srcToken_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("src"));
10  r_ = cfg.exists("r") ? cfg.getParameter<double>("r") : 2.;
11 
12  produces<double>("thrust");
13  //produces<double>("oblateness");
14  produces<double>("isotropy");
15  produces<double>("circularity");
16  produces<double>("sphericity");
17  produces<double>("aplanarity");
18  produces<double>("C");
19  produces<double>("D");
20 
21 }
22 
23 void put(edm::Event& evt, double value, const char* instanceName)
24 {
25  evt.put(std::make_unique<double>(value), instanceName);
26 }
27 
29 {
30  //std::cout << "<EventShapeVarsProducer::produce>:" << std::endl;
31 
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  put(evt, eventShapeVarsAlgo.isotropy(), "isotropy");
41  put(evt, eventShapeVarsAlgo.circularity(), "circularity");
42  put(evt, eventShapeVarsAlgo.sphericity(r_), "sphericity");
43  put(evt, eventShapeVarsAlgo.aplanarity(r_), "aplanarity");
44  put(evt, eventShapeVarsAlgo.C(r_), "C");
45  put(evt, eventShapeVarsAlgo.D(r_), "D");
46 }
47 
49 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
EventShapeVarsProducer(const edm::ParameterSet &)
tuple cfg
Definition: looper.py:293
Class for the calculation of several event shape variables.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#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_
void produce(edm::Event &, const edm::EventSetup &)
Definition: Thrust.h:38