3 #include "Math/VectorUtil.h"
10 : jetsToken_(consumes<std::
vector<pat::
Jet> >(cfg.getParameter<edm::
InputTag>(
"jets"))),
11 lepsToken_(consumes<edm::
View<
reco::RecoCandidate> >(cfg.getParameter<edm::
InputTag>(
"leps"))),
12 maxNJets_(cfg.getParameter<int>(
"maxNJets")),
13 useDeltaR_(cfg.getParameter<bool>(
"useDeltaR")),
14 useBTagging_(cfg.getParameter<bool>(
"useBTagging")),
15 bTagAlgorithm_(cfg.getParameter<std::
string>(
"bTagAlgorithm")),
16 minBDiscBJets_(cfg.getParameter<double>(
"minBDiscBJets")),
17 maxBDiscLightJets_(cfg.getParameter<double>(
"maxBDiscLightJets")) {
20 <<
"It has to be larger than 4 or can be set to -1 to take all jets.";
22 produces<std::vector<std::vector<int> > >();
23 produces<int>(
"NumberOfConsideredJets");
29 std::unique_ptr<std::vector<std::vector<int> > > pOut(
new std::vector<std::vector<int> >);
30 std::unique_ptr<int> pJetsConsidered(
new int);
32 std::vector<int>
match;
33 for (
unsigned int i = 0;
i < 4; ++
i)
45 if (leps->empty() || jets->size() < 4) {
46 pOut->push_back(match);
48 *pJetsConsidered = jets->size();
49 evt.
put(
std::move(pJetsConsidered),
"NumberOfConsideredJets");
55 maxNJets = jets->size();
56 *pJetsConsidered = maxNJets;
57 evt.
put(
std::move(pJetsConsidered),
"NumberOfConsideredJets");
59 std::vector<bool> isBJet;
60 std::vector<bool> isLJet;
63 for (
unsigned int idx = 0; idx < maxNJets; ++idx) {
78 for (
unsigned int idx = 0; idx < maxNJets; ++idx) {
79 if (
useBTagging_ && (!isLJet[idx] || (cntBJets <= 2 && isBJet[idx])))
81 for (
unsigned int jdx = (idx + 1); jdx < maxNJets; ++jdx) {
83 (!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) || (cntBJets == 3 && isBJet[idx] && isBJet[jdx])))
85 double dist =
distance((*jets)[idx].p4(), (*jets)[jdx].p4());
86 if (minDist < 0. || dist < minDist) {
96 wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4();
105 for (
unsigned int idx = 0; idx < maxNJets; ++idx) {
109 if ((
int)idx != lightQ && (
int)idx != lightQBar) {
110 double dist =
distance((*jets)[idx].p4(), wHad);
111 if (minDist < 0. || dist < minDist) {
125 for (
unsigned int idx = 0; idx < maxNJets; ++idx) {
129 if ((
int)idx != lightQ && (
int)idx != lightQBar && (
int)idx != hadB) {
130 double dist =
distance((*jets)[idx].p4(), (*leps)[0].p4());
131 if (minDist < 0. || dist < minDist) {
143 pOut->push_back(match);
152 return fabs(v1.theta() - v2.theta());
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
double maxBDiscLightJets_
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
TtSemiLepJetCombGeom(const edm::ParameterSet &)
double distance(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
void produce(edm::Event &evt, const edm::EventSetup &setup) override
std::string bTagAlgorithm_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
~TtSemiLepJetCombGeom() override
math::XYZTLorentzVector LorentzVector
Lorentz vector.