CMS 3D CMS Logo

EventShapeVarsProducer.cc
Go to the documentation of this file.
1 
36 
37 #include <memory>
38 #include <vector>
39 
41 public:
43 
44 private:
46  double r_;
47  unsigned fwmax_;
48 
49  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
50 };
51 
53  srcToken_ = consumes<edm::View<reco::Candidate>>(cfg.getParameter<edm::InputTag>("src"));
54  r_ = cfg.exists("r") ? cfg.getParameter<double>("r") : 2.;
55  fwmax_ = cfg.exists("fwmax") ? cfg.getParameter<unsigned>("fwmax") : 0;
56 
57  produces<double>("thrust");
58  //produces<double>("oblateness");
59  produces<double>("isotropy");
60  produces<double>("circularity");
61  produces<double>("sphericity");
62  produces<double>("aplanarity");
63  produces<double>("C");
64  produces<double>("D");
65  if (fwmax_ > 0)
66  produces<std::vector<double>>("FWmoments");
67 }
68 
69 void put(edm::Event& evt, double value, const char* instanceName) {
70  evt.put(std::make_unique<double>(value), instanceName);
71 }
72 
76 
77  Thrust thrustAlgo(objects->begin(), objects->end());
78  put(evt, thrustAlgo.thrust(), "thrust");
79  //put(evt, thrustAlgo.oblateness(), "oblateness");
80 
81  EventShapeVariables eventShapeVarsAlgo(*objects);
82  eventShapeVarsAlgo.set_r(r_);
83  put(evt, eventShapeVarsAlgo.isotropy(), "isotropy");
84  put(evt, eventShapeVarsAlgo.circularity(), "circularity");
85  put(evt, eventShapeVarsAlgo.sphericity(), "sphericity");
86  put(evt, eventShapeVarsAlgo.aplanarity(), "aplanarity");
87  put(evt, eventShapeVarsAlgo.C(), "C");
88  put(evt, eventShapeVarsAlgo.D(), "D");
89  if (fwmax_ > 0) {
90  eventShapeVarsAlgo.setFWmax(fwmax_);
91  auto vfw = std::make_unique<std::vector<double>>(eventShapeVarsAlgo.getFWmoments());
92  evt.put(std::move(vfw), "FWmoments");
93  }
94 }
95 
97 
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
EventShapeVarsProducer(const edm::ParameterSet &)
Class for the calculation of several event shape variables.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
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
def move(src, dest)
Definition: eostools.py:511