DataFormats
HepMCCandidate
src
FlavorHistoryEvent.cc
Go to the documentation of this file.
1
#include "
DataFormats/HepMCCandidate/interface/FlavorHistoryEvent.h
"
2
#include "TLorentzVector.h"
3
4
#include <iostream>
5
6
using namespace
reco
;
7
using namespace
std
;
8
9
// Loop over flavor histories, count number of genjets with
10
// flavor b, c, or l
11
void
FlavorHistoryEvent::cache
() {
12
bool
verbose
=
false
;
13
14
if
(
verbose
)
15
cout
<<
"----- Caching Flavor History Event -----"
<< endl;
16
// Set cached to false
17
cached_ =
false
;
18
nb_ = nc_ = 0;
19
highestFlavor_ = 0;
20
dR_ = 0.0;
21
flavorSource_ =
FlavorHistory::FLAVOR_NULL
;
22
23
// get list of flavor --> type --> dR.
24
// Will sort this later to determine event classification
25
vector<helpers::FlavorHistoryEventHelper> classification;
26
27
// get iterators to the history vector
28
const_iterator
i
= histories_.begin(), ibegin = histories_.begin(), iend = histories_.end();
29
// loop over the history vector and count the number of
30
// partons of flavors "b" and "c" that have a matched genjet.
31
for
(;
i
!= iend; ++
i
) {
32
FlavorHistory
const
& flavHist = *
i
;
33
if
(
verbose
)
34
cout
<<
" Processing flavor history: "
<<
i
- ibegin <<
" = "
<< endl << flavHist << endl;
35
CandidatePtr
const
& parton = flavHist.
parton
();
36
flavor_type
flavorSource = flavHist.
flavorSource
();
37
38
// Now examine the matched jets to see what the classification should be.
39
int
pdgId
= -1;
40
if
(parton.
isNonnull
())
41
pdgId
=
abs
(parton->
pdgId
());
42
ShallowClonePtrCandidate
const
& matchedJet = flavHist.
matchedJet
();
43
// Only count events with a matched genjet
44
if
(matchedJet.
masterClonePtr
().
isNonnull
()) {
45
TLorentzVector p41(matchedJet.
px
(), matchedJet.
py
(), matchedJet.
pz
(), matchedJet.
energy
());
46
if
(
pdgId
== 5)
47
nb_++;
48
if
(
pdgId
== 4)
49
nc_++;
50
51
// Get the sister genjet
52
ShallowClonePtrCandidate
const
& sisterJet =
i
->sisterJet();
53
TLorentzVector p42(sisterJet.
px
(), sisterJet.
py
(), sisterJet.
pz
(), sisterJet.
energy
());
54
55
// Now check the source.
56
double
dR
= -1;
57
if
(sisterJet.
masterClonePtr
().
isNonnull
()) {
58
dR
= p41.DeltaR(p42);
59
}
60
// Add to the vector to be sorted later
61
if
(
verbose
)
62
cout
<<
"Adding classification: pdgId = "
<<
pdgId
<<
", flavorSource = "
<< flavorSource <<
", dR = "
<<
dR
63
<< endl;
64
classification.push_back(
helpers::FlavorHistoryEventHelper
(
pdgId
, flavorSource,
dR
));
65
}
else
{
66
if
(
verbose
)
67
cout
<<
"No matched jet found, not adding to classification list"
<< endl;
68
}
69
}
70
71
// Sort by priority
72
73
// Priority is:
74
//
75
// 1. flavor (5 > 4)
76
// 2. type:
77
// 2a. Flavor decay
78
// 2b. Matrix element
79
// 2c. Flavor excitation
80
// 2d. Gluon splitting
81
// 3. delta R (if applicable)
82
if
(!classification.empty()) {
83
std::sort
(classification.begin(), classification.end());
84
85
if
(
verbose
) {
86
cout
<<
"Writing out list of classifications"
<< endl;
87
copy
(classification.begin(), classification.end(), ostream_iterator<helpers::FlavorHistoryEventHelper>(
cout
,
""
));
88
}
89
90
helpers::FlavorHistoryEventHelper
const
& best = *(classification.rbegin());
91
dR_ = best.
dR
;
92
highestFlavor_ = best.
flavor
;
93
flavorSource_ = best.
flavorSource
;
94
}
else
{
95
dR_ = -1.0;
96
highestFlavor_ = 0;
97
flavorSource_ =
FlavorHistory::FLAVOR_NULL
;
98
}
99
100
// now we're cached, can return values quickly
101
cached_ =
true
;
102
}
reco::FlavorHistoryEvent::cache
void cache()
Definition:
FlavorHistoryEvent.cc:11
mps_fire.i
i
Definition:
mps_fire.py:428
reco::FlavorHistory
Definition:
FlavorHistory.h:45
filterCSVwithJSON.copy
copy
Definition:
filterCSVwithJSON.py:36
reco::FlavorHistory::flavorSource
FLAVOR_T flavorSource() const
Definition:
FlavorHistory.h:85
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:46
reco::helpers::FlavorHistoryEventHelper
Definition:
FlavorHistoryEvent.h:50
reco::LeafCandidate::py
double py() const final
y coordinate of momentum vector
Definition:
LeafCandidate.h:142
reco::FlavorHistory::FLAVOR_T
FLAVOR_T
Definition:
FlavorHistory.h:47
verbose
static constexpr int verbose
Definition:
HLTExoticaSubAnalysis.cc:25
reco::helpers::FlavorHistoryEventHelper::dR
double dR
Definition:
FlavorHistoryEvent.h:58
reco::FlavorHistory::parton
const reco::CandidatePtr & parton() const
Definition:
FlavorHistory.h:91
jetUpdater_cfi.sort
sort
Definition:
jetUpdater_cfi.py:29
EgammaValidation_cff.pdgId
pdgId
Definition:
EgammaValidation_cff.py:117
reco::helpers::FlavorHistoryEventHelper::flavorSource
FlavorHistory::FLAVOR_T flavorSource
Definition:
FlavorHistoryEvent.h:57
reco::FlavorHistory::FLAVOR_NULL
Definition:
FlavorHistory.h:48
reco::Candidate::pdgId
virtual int pdgId() const =0
PDG identifier.
edm::Ptr< Candidate >
std
Definition:
JetResolutionObject.h:76
reco::FlavorHistoryEvent::const_iterator
collection_type::const_iterator const_iterator
Definition:
FlavorHistoryEvent.h:91
reco::helpers::FlavorHistoryEventHelper::flavor
int flavor
Definition:
FlavorHistoryEvent.h:56
edm::Ptr::isNonnull
bool isNonnull() const
Checks for non-null.
Definition:
Ptr.h:146
reco::LeafCandidate::energy
double energy() const final
energy
Definition:
LeafCandidate.h:125
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition:
HGC3DClusterGenMatchSelector_cfi.py:7
reco::ShallowClonePtrCandidate
Definition:
ShallowClonePtrCandidate.h:15
reco::ShallowClonePtrCandidate::masterClonePtr
const CandidatePtr & masterClonePtr() const override
returns reference to master clone pointer
Definition:
ShallowClonePtrCandidate.cc:20
reco::LeafCandidate::px
double px() const final
x coordinate of momentum vector
Definition:
LeafCandidate.h:140
reco::LeafCandidate::pz
double pz() const final
z coordinate of momentum vector
Definition:
LeafCandidate.h:144
reco::FlavorHistory::matchedJet
const reco::ShallowClonePtrCandidate & matchedJet() const
Definition:
FlavorHistory.h:94
FlavorHistoryEvent.h
Generated for CMSSW Reference Manual by
1.8.16