CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauDumper.cc
Go to the documentation of this file.
1 /*
2  * =====================================================================================
3  *
4  * Filename: RecoTauDumper.cc
5  *
6  * Description: Dump information about reco::PFTaus
7  *
8  * Author: Evan K. Friis, UC Davis
9  *
10  *
11  * =====================================================================================
12  */
13 #include <boost/foreach.hpp>
14 #include <sstream>
15 
19 
22 
23 // Methods to write the different types
24 namespace {
25 
26 void write(std::ostringstream& output, const reco::PFTau& tau) {
27  output << " ------------------------------------" << std::endl;
28  output << tau << std::endl;
29  tau.dump(output);
30  if (tau.pfTauTagInfoRef().isNonnull()) {
31  output << " TTInfoJetRefID: "
32  << tau.pfTauTagInfoRef()->pfjetRef().id() << ":"
33  << tau.pfTauTagInfoRef()->pfjetRef().key() << std::endl;
34  output << " TTInfoJetRef: " << *(tau.pfTauTagInfoRef()->pfjetRef());
35  }
36  if (tau.jetRef().isNonnull()) {
37  output << " JetRefID: "
38  << tau.jetRef().id() << ":"
39  << tau.jetRef().key() << std::endl;
40  output << " JetRef: " << *(tau.jetRef());
41 
42  }
43  output << std::endl;
44 }
45 
46 void write(std::ostringstream& output, const reco::PFJet& jet) {
47  output << " ------------------------------------" << std::endl;
48  output << jet << std::endl;
49  BOOST_FOREACH(const reco::PFCandidatePtr& cand, jet.getPFConstituents()) {
50  output << " --> " << *cand << std::endl;
51  }
52  output << std::endl;
53 }
54 
55 void write(std::ostringstream& output, const reco::Candidate& cand) {
56  output << " ------------------------------------" << std::endl;
57  output <<
58  " candidate (pt/eta/phi): (" << cand.pt() << "/"
59  << cand.eta() << "/"
60  << cand.phi() << ")" << std::endl;
61  output << std::endl;
62 }
63 
64 }
65 
66 template<typename T>
68  public:
69  explicit CollectionDumper(const edm::ParameterSet& pset):
70  src_(pset.getParameter<edm::InputTag>("src")),
71  moduleName_(pset.getParameter<std::string>("@module_label")){}
72  virtual ~CollectionDumper() {}
73  virtual void analyze(const edm::Event& evt, const edm::EventSetup& es) override;
74  private:
77 };
78 
79 template<typename T> void
81  typedef edm::View<T> TView;
82  edm::Handle<TView> view;
83  evt.getByLabel(src_, view);
84 
85  std::ostringstream output;
86  output << " * * * <" << moduleName_
87  << "> Dump - source: [" << src_ << "]" << std::endl;
88 
89  BOOST_FOREACH(const T& obj, *view) {
90  write(output, obj);
91  }
92  std::cout << output.str();
93 }
94 
98 
100 DEFINE_FWK_MODULE(RecoTauDumper);
101 DEFINE_FWK_MODULE(PFJetDumper);
102 DEFINE_FWK_MODULE(CandidateDumper);
std::string moduleName_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
const PFJetRef & jetRef() const
Definition: PFTau.cc:58
virtual void analyze(const edm::Event &evt, const edm::EventSetup &es) override
CollectionDumper(const edm::ParameterSet &pset)
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
CollectionDumper< reco::Candidate > CandidateDumper
key_type key() const
Accessor for product key.
Definition: Ref.h:264
Jets made from PFObjects.
Definition: PFJet.h:21
ProductID id() const
Accessor for product ID.
Definition: Ref.h:258
edm::InputTag src_
void dump(std::ostream &out=std::cout) const
prints information on this PFTau
Definition: PFTau.cc:276
CollectionDumper< reco::PFTau > RecoTauDumper
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
const PFTauTagInfoRef & pfTauTagInfoRef() const
Definition: PFTau.cc:61
CollectionDumper< reco::PFJet > PFJetDumper
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
tuple cout
Definition: gather_cfg.py:121
virtual ~CollectionDumper()
long double T
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity