CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
ParticleDecayDrawer Class Reference
Inheritance diagram for ParticleDecayDrawer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 ParticleDecayDrawer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

bool accept (const reco::Candidate &, const std::list< const reco::Candidate * > &) const
 accept candidate More...
 
void analyze (const edm::Event &, const edm::EventSetup &) override
 
std::string decay (const reco::Candidate &, std::list< const reco::Candidate * > &) const
 
bool hasValidDaughters (const reco::Candidate &) const
 has valid daughters in the chain More...
 
std::string printP4 (const reco::Candidate &) const
 print 4 momenta More...
 
bool select (const reco::Candidate &) const
 select candidate More...
 

Private Attributes

edm::ESHandle< ParticleDataTablepdt_
 
bool printP4_
 print parameters More...
 
bool printPtEtaPhi_
 
bool printVertex_
 
edm::EDGetTokenT< edm::View
< reco::Candidate > > 
srcToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 12 of file ParticleDecayDrawer.cc.

Constructor & Destructor Documentation

ParticleDecayDrawer::ParticleDecayDrawer ( const edm::ParameterSet cfg)

Definition at line 44 of file ParticleDecayDrawer.cc.

44  :
46  printP4_( cfg.getUntrackedParameter<bool>( "printP4", false ) ),
47  printPtEtaPhi_( cfg.getUntrackedParameter<bool>( "printPtEtaPhi", false ) ),
48  printVertex_( cfg.getUntrackedParameter<bool>( "printVertex", false ) ) {
49 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool printP4_
print parameters
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_

Member Function Documentation

bool ParticleDecayDrawer::accept ( const reco::Candidate ,
const std::list< const reco::Candidate * > &   
) const
private

accept candidate

Definition at line 51 of file ParticleDecayDrawer.cc.

References EnergyCorrector::c, spr::find(), and select().

Referenced by Vispa.Gui.BoxContentDialog.BoxContentDialog::apply(), Vispa.Plugins.ConfigEditor.ToolDialog.ToolDialog::apply(), decay(), and esMonitoring.FDJsonServer::handle_accept().

51  {
52  if( find( skip.begin(), skip.end(), & c ) != skip.end() ) return false;
53  return select( c );
54 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
bool select(const reco::Candidate &) const
select candidate
void ParticleDecayDrawer::analyze ( const edm::Event event,
const edm::EventSetup es 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 68 of file ParticleDecayDrawer.cc.

References gather_cfg::cout, TauDecayModes::dec, decay(), edm::EventSetup::getData(), j, visualization-live-secondInstance_cfg::m, reco::Candidate::mother(), gen::n, class-composition::nodes, AlCaHLTBitMon_ParallelJobs::p, pdt_, select(), runGlobalFakeInputProducer::skip, and srcToken_.

68  {
69  es.getData( pdt_ );
70  Handle<View<Candidate> > particles;
71  event.getByToken( srcToken_, particles );
72  list<const Candidate *> skip;
73  vector<const Candidate *> nodes, moms;
74  for( View<Candidate>::const_iterator p = particles->begin();
75  p != particles->end(); ++ p ) {
76  if( p->numberOfMothers() > 1 ) {
77  if ( select( * p ) ) {
78  skip.push_back( & * p );
79  nodes.push_back( & * p );
80  for( size_t j = 0; j < p->numberOfMothers(); ++ j ) {
81  const Candidate * mom = p->mother( j );
82  const Candidate * grandMom;
83  while ( ( grandMom = mom->mother() ) != 0 )
84  mom = grandMom;
85  if ( select( * mom ) ) {
86  moms.push_back( mom );
87  }
88  }
89  }
90  }
91  }
92  cout << "-- decay: --" << endl;
93  if( moms.size() > 0 ) {
94  if ( moms.size() > 1 )
95  for( size_t m = 0; m < moms.size(); ++ m ) {
96  string dec = decay( * moms[ m ], skip );
97  if ( ! dec.empty() )
98  cout << "{ " << dec << " } ";
99  }
100  else
101  cout << decay( * moms[ 0 ], skip );
102  }
103  if ( nodes.size() > 0 ) {
104  cout << "-> ";
105  if ( nodes.size() > 1 ) {
106  for( size_t n = 0; n < nodes.size(); ++ n ) {
107  skip.remove( nodes[ n ] );
108  string dec = decay( * nodes[ n ], skip );
109  if ( ! dec.empty() ) {
110  if ( dec.find( "->", 0 ) != string::npos )
111  cout << " ( " << dec << " )";
112  else
113  cout << " " << dec;
114  }
115  }
116  } else {
117  skip.remove( nodes[ 0 ] );
118  cout << decay( * nodes[ 0 ], skip );
119  }
120  }
121  cout << endl;
122 }
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
edm::ESHandle< ParticleDataTable > pdt_
void getData(T &iHolder) const
Definition: EventSetup.h:78
std::string decay(const reco::Candidate &, std::list< const reco::Candidate * > &) const
bool select(const reco::Candidate &) const
select candidate
int j
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
tuple cout
Definition: gather_cfg.py:121
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_
string ParticleDecayDrawer::decay ( const reco::Candidate ,
std::list< const reco::Candidate * > &   
) const
private

Definition at line 132 of file ParticleDecayDrawer.cc.

References accept(), assert(), EnergyCorrector::c, ztail::d, reco::Candidate::daughter(), TauDecayModes::dec, spr::find(), i, reco::Candidate::numberOfDaughters(), dbtoconf::out, reco::Candidate::pdgId(), pdt_, printP4(), and runGlobalFakeInputProducer::skip.

Referenced by analyze().

133  {
134  string out;
135  if ( find( skip.begin(), skip.end(), & c ) != skip.end() )
136  return out;
137  skip.push_back( & c );
138 
139 
140  int id = c.pdgId();
141  const ParticleData * pd = pdt_->particle( id );
142  assert( pd != 0 );
143  out += ( pd->name() + printP4( c ) );
144 
145  size_t validDau = 0, ndau = c.numberOfDaughters();
146  for( size_t i = 0; i < ndau; ++ i )
147  if ( accept( * c.daughter( i ), skip ) )
148  ++ validDau;
149  if ( validDau == 0 ) return out;
150 
151  out += " ->";
152 
153  for( size_t i = 0; i < ndau; ++ i ) {
154  const Candidate * d = c.daughter( i );
155  if ( accept( * d, skip ) ) {
156  string dec = decay( * d, skip );
157  if ( dec.find( "->", 0 ) != string::npos )
158  out += ( " ( " + dec + " )" );
159  else
160  out += ( " " + dec );
161  }
162  }
163  return out;
164 }
int i
Definition: DBlmapReader.cc:9
bool accept(const reco::Candidate &, const std::list< const reco::Candidate * > &) const
accept candidate
std::string printP4(const reco::Candidate &) const
print 4 momenta
assert(m_qm.get())
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
edm::ESHandle< ParticleDataTable > pdt_
tuple d
Definition: ztail.py:151
std::string decay(const reco::Candidate &, std::list< const reco::Candidate * > &) const
HepPDT::ParticleData ParticleData
tuple out
Definition: dbtoconf.py:99
bool ParticleDecayDrawer::hasValidDaughters ( const reco::Candidate c) const
private

has valid daughters in the chain

Definition at line 60 of file ParticleDecayDrawer.cc.

References reco::Candidate::daughter(), i, reco::Candidate::numberOfDaughters(), and select().

60  {
61  size_t ndau = c.numberOfDaughters();
62  for( size_t i = 0; i < ndau; ++ i )
63  if ( select( * c.daughter( i ) ) )
64  return true;
65  return false;
66 }
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual size_type numberOfDaughters() const =0
number of daughters
bool select(const reco::Candidate &) const
select candidate
string ParticleDecayDrawer::printP4 ( const reco::Candidate c) const
private

print 4 momenta

Definition at line 124 of file ParticleDecayDrawer.cc.

References gather_cfg::cout, reco::Candidate::energy(), reco::Candidate::eta(), reco::Candidate::phi(), printP4_, printPtEtaPhi_, printVertex_, reco::Candidate::pt(), reco::Candidate::px(), reco::Candidate::py(), reco::Candidate::pz(), reco::Candidate::vx(), reco::Candidate::vy(), and reco::Candidate::vz().

Referenced by decay().

124  {
125  ostringstream cout;
126  if ( printP4_ ) cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")";
127  if ( printPtEtaPhi_ ) cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
128  if ( printVertex_ ) cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
129  return cout.str();
130 }
virtual double energy() const =0
energy
virtual double pt() const =0
transverse momentum
virtual double pz() const =0
z coordinate of momentum vector
virtual double vx() const =0
x coordinate of vertex position
virtual double vy() const =0
y coordinate of vertex position
bool printP4_
print parameters
virtual double py() const =0
y coordinate of momentum vector
virtual double px() const =0
x coordinate of momentum vector
virtual double vz() const =0
z coordinate of vertex position
tuple cout
Definition: gather_cfg.py:121
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
bool ParticleDecayDrawer::select ( const reco::Candidate c) const
private

select candidate

Definition at line 56 of file ParticleDecayDrawer.cc.

References reco::Candidate::status().

Referenced by accept(), analyze(), and hasValidDaughters().

56  {
57  return c.status() == 3;
58 }
virtual int status() const =0
status word

Member Data Documentation

edm::ESHandle<ParticleDataTable> ParticleDecayDrawer::pdt_
private

Definition at line 19 of file ParticleDecayDrawer.cc.

Referenced by analyze(), and decay().

bool ParticleDecayDrawer::printP4_
private

print parameters

Definition at line 21 of file ParticleDecayDrawer.cc.

Referenced by printP4().

bool ParticleDecayDrawer::printPtEtaPhi_
private

Definition at line 21 of file ParticleDecayDrawer.cc.

Referenced by printP4().

bool ParticleDecayDrawer::printVertex_
private

Definition at line 21 of file ParticleDecayDrawer.cc.

Referenced by printP4().

edm::EDGetTokenT<edm::View<reco::Candidate> > ParticleDecayDrawer::srcToken_
private

Definition at line 17 of file ParticleDecayDrawer.cc.

Referenced by analyze().