15 bool isValid(
const int& idx,
const std::vector<pat::Jet>&
jets)
const {
16 return (0 <= idx && idx < (
int)jets.size());
30 : jetsToken_(consumes<std::
vector<pat::
Jet>>(cfg.getParameter<edm::
InputTag>(
"jets"))),
31 lepsToken_(consumes<edm::
View<
reco::RecoCandidate>>(cfg.getParameter<edm::
InputTag>(
"leps"))),
32 maxNJets_(cfg.getParameter<int>(
"maxNJets")),
33 wMass_(cfg.getParameter<double>(
"wMass")),
34 useBTagging_(cfg.getParameter<bool>(
"useBTagging")),
35 bTagAlgorithm_(cfg.getParameter<std::
string>(
"bTagAlgorithm")),
36 minBDiscBJets_(cfg.getParameter<double>(
"minBDiscBJets")),
37 maxBDiscLightJets_(cfg.getParameter<double>(
"maxBDiscLightJets")) {
40 <<
"It has to be larger than 4 or can be set to -1 to take all jets.";
42 produces<std::vector<std::vector<int>>>();
43 produces<int>(
"NumberOfConsideredJets");
47 auto pOut = std::make_unique<std::vector<std::vector<int>>>();
48 auto pJetsConsidered = std::make_unique<int>(0);
50 std::vector<int>
match;
51 for (
unsigned int i = 0;
i < 4; ++
i)
61 if (leps.empty() ||
jets.size() < 4) {
62 pOut->push_back(match);
64 *pJetsConsidered =
jets.size();
65 evt.
put(
std::move(pJetsConsidered),
"NumberOfConsideredJets");
71 maxNJets =
jets.size();
72 *pJetsConsidered = maxNJets;
73 evt.
put(
std::move(pJetsConsidered),
"NumberOfConsideredJets");
75 std::vector<bool> isBJet;
76 std::vector<bool> isLJet;
79 for (
unsigned int idx = 0; idx < maxNJets; ++idx) {
92 std::vector<int> closestToWMassIndices;
93 closestToWMassIndices.push_back(-1);
94 closestToWMassIndices.push_back(-1);
95 for (
unsigned idx = 0; idx < maxNJets; ++idx) {
96 if (
useBTagging_ && (!isLJet[idx] || (cntBJets <= 2 && isBJet[idx])))
98 for (
unsigned jdx = (idx + 1); jdx < maxNJets; ++jdx) {
100 (!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) || (cntBJets == 3 && isBJet[idx] && isBJet[jdx])))
103 if (wDist < 0. || wDist > fabs(sum.mass() -
wMass_)) {
104 wDist = fabs(sum.mass() -
wMass_);
105 closestToWMassIndices.clear();
106 closestToWMassIndices.push_back(idx);
107 closestToWMassIndices.push_back(jdx);
119 for (
unsigned idx = 0; idx < maxNJets; ++idx) {
123 if ((
int)idx != closestToWMassIndices[0] && (
int)idx != closestToWMassIndices[1]) {
125 jets[closestToWMassIndices[0]].p4() +
jets[closestToWMassIndices[1]].p4() +
jets[idx].p4();
126 if (maxPt < 0. || maxPt < sum.pt()) {
141 for (
unsigned idx = 0; idx < maxNJets; ++idx) {
145 if ((
int)idx != closestToWMassIndices[0] && (
int)idx != closestToWMassIndices[1] && (
int)idx != hadB) {
147 if (maxPt < 0. || maxPt < sum.pt()) {
159 pOut->push_back(match);
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
TtSemiLepJetCombWMassMaxSumPt(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
std::string bTagAlgorithm_
bool get(ProductID const &oid, Handle< PROD > &result) const
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
bool isValid(const int &idx, const std::vector< pat::Jet > &jets) const
double maxBDiscLightJets_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &setup) const override