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  * $Id $
11  *
12  * =====================================================================================
13  */
14 #include <boost/foreach.hpp>
15 #include <sstream>
16 
20 
23 
24 // Methods to write the different types
25 namespace {
26 
27 void write(std::ostringstream& output, const reco::PFTau& tau) {
28  output << " ------------------------------------" << std::endl;
29  output << tau << std::endl;
30  tau.dump(output);
31  if (tau.pfTauTagInfoRef().isNonnull()) {
32  output << " TTInfoJetRefID: "
33  << tau.pfTauTagInfoRef()->pfjetRef().id() << ":"
34  << tau.pfTauTagInfoRef()->pfjetRef().key() << std::endl;
35  output << " TTInfoJetRef: " << *(tau.pfTauTagInfoRef()->pfjetRef());
36  }
37  if (tau.jetRef().isNonnull()) {
38  output << " JetRefID: "
39  << tau.jetRef().id() << ":"
40  << tau.jetRef().key() << std::endl;
41  output << " JetRef: " << *(tau.jetRef());
42 
43  }
44  output << std::endl;
45 }
46 
47 void write(std::ostringstream& output, const reco::PFJet& jet) {
48  output << " ------------------------------------" << std::endl;
49  output << jet << std::endl;
50  BOOST_FOREACH(const reco::PFCandidatePtr& cand, jet.getPFConstituents()) {
51  output << " --> " << *cand << std::endl;
52  }
53  output << std::endl;
54 }
55 
56 void write(std::ostringstream& output, const reco::Candidate& cand) {
57  output << " ------------------------------------" << std::endl;
58  output <<
59  " candidate (pt/eta/phi): (" << cand.pt() << "/"
60  << cand.eta() << "/"
61  << cand.phi() << ")" << std::endl;
62  output << std::endl;
63 }
64 
65 }
66 
67 template<typename T>
69  public:
70  explicit CollectionDumper(const edm::ParameterSet& pset):
71  src_(pset.getParameter<edm::InputTag>("src")),
72  moduleName_(pset.getParameter<std::string>("@module_label")){}
73  virtual ~CollectionDumper() {}
74  virtual void analyze(const edm::Event& evt, const edm::EventSetup& es);
75  private:
77  std::string moduleName_;
78 };
79 
80 template<typename T> void
82  typedef edm::View<T> TView;
83  edm::Handle<TView> view;
84  evt.getByLabel(src_, view);
85 
86  std::ostringstream output;
87  output << " * * * <" << moduleName_
88  << "> Dump - source: [" << src_ << "]" << std::endl;
89 
90  BOOST_FOREACH(const T& obj, *view) {
91  write(output, obj);
92  }
93  std::cout << output.str();
94 }
95 
99 
101 DEFINE_FWK_MODULE(RecoTauDumper);
102 DEFINE_FWK_MODULE(PFJetDumper);
103 DEFINE_FWK_MODULE(CandidateDumper);
std::string moduleName_
const PFJetRef & jetRef() const
Definition: PFTau.cc:50
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
Jets made from PFObjects.
Definition: PFJet.h:22
edm::InputTag src_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
void dump(std::ostream &out=std::cout) const
prints information on this PFTau
Definition: PFTau.cc:209
CollectionDumper< reco::PFTau > RecoTauDumper
virtual void analyze(const edm::Event &evt, const edm::EventSetup &es)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
const PFTauTagInfoRef & pfTauTagInfoRef() const
Definition: PFTau.cc:53
CollectionDumper< reco::PFJet > PFJetDumper
key_type key() const
Accessor for product key.
Definition: Ref.h:266
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:53
tuple cout
Definition: gather_cfg.py:121
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
virtual ~CollectionDumper()
long double T
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity