CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
ParticleDecayDrawer Class Reference
Inheritance diagram for ParticleDecayDrawer:
edm::one::EDAnalyzer<> edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 ParticleDecayDrawer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

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_
 
edm::ESGetToken< ParticleDataTable, edm::DefaultRecordpdtToken_
 
bool printP4_
 print parameters More...
 
bool printPtEtaPhi_
 
bool printVertex_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 13 of file ParticleDecayDrawer.cc.

Constructor & Destructor Documentation

◆ ParticleDecayDrawer()

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

Definition at line 47 of file ParticleDecayDrawer.cc.

48  : srcToken_(consumes<edm::View<reco::Candidate> >(cfg.getParameter<InputTag>("src"))),
49  pdtToken_(esConsumes<ParticleDataTable, edm::DefaultRecord>()),
50  printP4_(cfg.getUntrackedParameter<bool>("printP4", false)),
51  printPtEtaPhi_(cfg.getUntrackedParameter<bool>("printPtEtaPhi", false)),
52  printVertex_(cfg.getUntrackedParameter<bool>("printVertex", false)) {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool printP4_
print parameters
edm::ESGetToken< ParticleDataTable, edm::DefaultRecord > pdtToken_
edm::EDGetTokenT< edm::View< reco::Candidate > > srcToken_

Member Function Documentation

◆ accept()

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

accept candidate

Definition at line 54 of file ParticleDecayDrawer.cc.

References HltBtagPostValidation_cff::c, spr::find(), select(), and optionsL1T::skip.

Referenced by decay(), and esMonitoring.FDJsonServer::handle_accept().

54  {
55  if (find(skip.begin(), skip.end(), &c) != skip.end())
56  return false;
57  return select(c);
58 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
bool select(const reco::Candidate &) const
select candidate

◆ analyze()

void ParticleDecayDrawer::analyze ( const edm::Event event,
const edm::EventSetup es 
)
overrideprivatevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 70 of file ParticleDecayDrawer.cc.

References gather_cfg::cout, TauDecayModes::dec, decay(), edm::EventSetup::getHandle(), dqmiolumiharvest::j, visualization-live-secondInstance_cfg::m, reco::Candidate::mother(), dqmiodumpmetadata::n, class-composition::nodes, AlCaHLTBitMon_ParallelJobs::p, ecalTrigSettings_cff::particles, pdt_, pdtToken_, select(), optionsL1T::skip, and srcToken_.

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

◆ decay()

string ParticleDecayDrawer::decay ( const reco::Candidate ,
std::list< const reco::Candidate *> &   
) const
private

Definition at line 136 of file ParticleDecayDrawer.cc.

References accept(), cms::cuda::assert(), HltBtagPostValidation_cff::c, ztail::d, TauDecayModes::dec, spr::find(), mps_fire::i, MillePedeFileConverter_cfg::out, pdt_, printP4(), and optionsL1T::skip.

Referenced by analyze().

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

◆ hasValidDaughters()

bool ParticleDecayDrawer::hasValidDaughters ( const reco::Candidate c) const
private

has valid daughters in the chain

Definition at line 62 of file ParticleDecayDrawer.cc.

References HltBtagPostValidation_cff::c, mps_fire::i, and select().

62  {
63  size_t ndau = c.numberOfDaughters();
64  for (size_t i = 0; i < ndau; ++i)
65  if (select(*c.daughter(i)))
66  return true;
67  return false;
68 }
bool select(const reco::Candidate &) const
select candidate

◆ printP4()

string ParticleDecayDrawer::printP4 ( const reco::Candidate c) const
private

print 4 momenta

Definition at line 125 of file ParticleDecayDrawer.cc.

References HltBtagPostValidation_cff::c, gather_cfg::cout, printP4_, printPtEtaPhi_, and printVertex_.

Referenced by decay().

125  {
126  ostringstream cout;
127  if (printP4_)
128  cout << " (" << c.px() << ", " << c.py() << ", " << c.pz() << "; " << c.energy() << ")";
129  if (printPtEtaPhi_)
130  cout << " [" << c.pt() << ": " << c.eta() << ", " << c.phi() << "]";
131  if (printVertex_)
132  cout << " {" << c.vx() << ", " << c.vy() << ", " << c.vz() << "}";
133  return cout.str();
134 }
bool printP4_
print parameters

◆ select()

bool ParticleDecayDrawer::select ( const reco::Candidate c) const
private

select candidate

Definition at line 60 of file ParticleDecayDrawer.cc.

References HltBtagPostValidation_cff::c.

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

60 { return c.status() == 3; }

Member Data Documentation

◆ pdt_

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

Definition at line 21 of file ParticleDecayDrawer.cc.

Referenced by analyze(), and decay().

◆ pdtToken_

edm::ESGetToken<ParticleDataTable, edm::DefaultRecord> ParticleDecayDrawer::pdtToken_
private

Definition at line 22 of file ParticleDecayDrawer.cc.

Referenced by analyze().

◆ printP4_

bool ParticleDecayDrawer::printP4_
private

print parameters

Definition at line 24 of file ParticleDecayDrawer.cc.

Referenced by printP4().

◆ printPtEtaPhi_

bool ParticleDecayDrawer::printPtEtaPhi_
private

Definition at line 24 of file ParticleDecayDrawer.cc.

Referenced by printP4().

◆ printVertex_

bool ParticleDecayDrawer::printVertex_
private

Definition at line 24 of file ParticleDecayDrawer.cc.

Referenced by printP4().

◆ srcToken_

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

Definition at line 19 of file ParticleDecayDrawer.cc.

Referenced by analyze().