CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
reco::isodeposit::OtherJetConstituentsDeltaRVeto Class Reference

#include <EventDependentAbsVetos.h>

Inheritance diagram for reco::isodeposit::OtherJetConstituentsDeltaRVeto:
reco::isodeposit::EventDependentAbsVeto reco::isodeposit::AbsVeto

Public Member Functions

void centerOn (double eta, double phi) override
 Set axis for matching jets. More...
 
 OtherJetConstituentsDeltaRVeto (Direction dir, const edm::InputTag &jets, double dRjet, const edm::InputTag &pfCandAssocMap, double dRconstituent, edm::ConsumesCollector &iC)
 Create a veto specifying the input collection of the jets, the candidates, and the deltaR. More...
 
void setEvent (const edm::Event &evt, const edm::EventSetup &es) override
 Picks up the directions of the given candidates. More...
 
bool veto (double eta, double phi, float value) const override
 
 ~OtherJetConstituentsDeltaRVeto () override
 
- Public Member Functions inherited from reco::isodeposit::EventDependentAbsVeto
 ~EventDependentAbsVeto () override
 
- Public Member Functions inherited from reco::isodeposit::AbsVeto
virtual ~AbsVeto ()
 

Private Types

typedef edm::AssociationMap< edm::OneToMany< std::vector< reco::PFJet >, std::vector< reco::PFCandidate >, unsigned int > > JetToPFCandidateAssociation
 

Private Member Functions

void initialize ()
 

Private Attributes

double dR2constituent_
 
double dR2jet_
 
const edm::Eventevt_
 
std::vector< Directionitems_
 
edm::EDGetTokenT< reco::PFJetCollectionsrcJets_
 
edm::EDGetTokenT< JetToPFCandidateAssociationsrcPFCandAssocMap_
 
Direction vetoDir_
 

Detailed Description

Definition at line 64 of file EventDependentAbsVetos.h.

Member Typedef Documentation

◆ JetToPFCandidateAssociation

Definition at line 102 of file EventDependentAbsVetos.h.

Constructor & Destructor Documentation

◆ OtherJetConstituentsDeltaRVeto()

reco::isodeposit::OtherJetConstituentsDeltaRVeto::OtherJetConstituentsDeltaRVeto ( Direction  dir,
const edm::InputTag jets,
double  dRjet,
const edm::InputTag pfCandAssocMap,
double  dRconstituent,
edm::ConsumesCollector iC 
)
inline

Create a veto specifying the input collection of the jets, the candidates, and the deltaR.

Definition at line 67 of file EventDependentAbsVetos.h.

73  : evt_(nullptr),
74  vetoDir_(dir),
76  dR2jet_(dRjet * dRjet),
78  dR2constituent_(dRconstituent * dRconstituent) {
79  //std::cout << "<OtherJetConstituentsDeltaRVeto::OtherJetConstituentsDeltaRVeto>:" << std::endl;
80  //std::cout << " vetoDir: eta = " << vetoDir_.eta() << ", phi = " << vetoDir_.phi() << std::endl;
81  //std::cout << " srcJets = " << srcJets_.label() << ":" << srcJets_.instance() << std::endl;
82  //std::cout << " dRjet = " << sqrt(dR2jet_) << std::endl;
83  //std::cout << " srcPFCandAssocMap = " << srcPFCandAssocMap_.label() << ":" << srcPFCandAssocMap_.instance() << std::endl;
84  //std::cout << " dRconstituent = " << sqrt(dR2constituent_) << std::endl;
85  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::AssociationMap< edm::OneToMany< std::vector< reco::PFJet >, std::vector< reco::PFCandidate >, unsigned int > > JetToPFCandidateAssociation
edm::EDGetTokenT< JetToPFCandidateAssociation > srcPFCandAssocMap_
edm::EDGetTokenT< reco::PFJetCollection > srcJets_
std::vector< PFJet > PFJetCollection
collection of PFJet objects

◆ ~OtherJetConstituentsDeltaRVeto()

reco::isodeposit::OtherJetConstituentsDeltaRVeto::~OtherJetConstituentsDeltaRVeto ( )
inlineoverride

Definition at line 88 of file EventDependentAbsVetos.h.

88 {}

Member Function Documentation

◆ centerOn()

void reco::isodeposit::OtherJetConstituentsDeltaRVeto::centerOn ( double  eta,
double  phi 
)
overridevirtual

Set axis for matching jets.

Reimplemented from reco::isodeposit::AbsVeto.

Definition at line 85 of file EventDependentAbsVetos.cc.

References PVValHelper::eta, and initialize().

85  {
86  //std::cout << "<OtherJetConstituentsDeltaRVeto::centerOn>:" << std::endl;
87  //std::cout << " eta = " << eta << std::endl;
88  //std::cout << " phi = " << phi << std::endl;
89  vetoDir_ = Direction(eta, phi);
90  initialize();
91 }

◆ initialize()

void reco::isodeposit::OtherJetConstituentsDeltaRVeto::initialize ( void  )
private

Definition at line 54 of file EventDependentAbsVetos.cc.

References cms::cuda::assert(), edm::RefVector< C, T, F >::begin(), reco::deltaR2(), ALPAKA_ACCELERATOR_NAMESPACE::dR2(), edm::RefVector< C, T, F >::end(), edm::Ref< C, T, F >::isNonnull(), metsig::jet, PDWG_EXODelayedJetMET_cff::jets, and edm::RefVector< C, T, F >::push_back().

54  {
55  //std::cout << "<OtherJetConstituentsDeltaRVeto::initialize>:" << std::endl;
56  //std::cout << " vetoDir: eta = " << vetoDir_.eta() << ", phi = " << vetoDir_.phi() << std::endl;
57  assert(evt_);
58  items_.clear();
62  evt_->getByToken(srcPFCandAssocMap_, jetToPFCandMap);
63  double dR2min = dR2jet_;
64  reco::PFJetRef matchedJet;
65  size_t numJets = jets->size();
66  for (size_t jetIndex = 0; jetIndex < numJets; ++jetIndex) {
67  reco::PFJetRef jet(jets, jetIndex);
68  double dR2 = ::deltaR2(vetoDir_.eta(), vetoDir_.phi(), jet->eta(), jet->phi());
69  //std::cout << "jet #" << jetIndex << ": Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->phi() << " (dR = " << sqrt(dR2) << ")" << std::endl;
70  if (dR2 < dR2min) {
71  matchedJet = jet;
72  dR2min = dR2;
73  }
74  }
75  if (matchedJet.isNonnull()) {
76  edm::RefVector<reco::PFCandidateCollection> pfCandsMappedToJet = (*jetToPFCandMap)[matchedJet];
77  for (edm::RefVector<reco::PFCandidateCollection>::const_iterator pfCand = pfCandsMappedToJet.begin();
78  pfCand != pfCandsMappedToJet.end();
79  ++pfCand) {
80  items_.push_back(Direction((*pfCand)->eta(), (*pfCand)->phi()));
81  }
82  }
83 }
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:526
assert(be >=bs)
ALPAKA_FN_ACC static ALPAKA_FN_INLINE float dR2(Position4 pos1, Position4 pos2)
edm::EDGetTokenT< JetToPFCandidateAssociation > srcPFCandAssocMap_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
edm::EDGetTokenT< reco::PFJetCollection > srcJets_
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ setEvent()

void reco::isodeposit::OtherJetConstituentsDeltaRVeto::setEvent ( const edm::Event evt,
const edm::EventSetup es 
)
overridevirtual

Picks up the directions of the given candidates.

Implements reco::isodeposit::EventDependentAbsVeto.

Definition at line 50 of file EventDependentAbsVetos.cc.

50  {
51  //std::cout << "<OtherJetConstituentsDeltaRVeto::setEvent>:" << std::endl;
52  evt_ = &evt;
53 }

◆ veto()

bool reco::isodeposit::OtherJetConstituentsDeltaRVeto::veto ( double  eta,
double  phi,
float  value 
) const
overridevirtual

Return "true" if a deposit at specific (eta,phi) with that value must be vetoed in the sum This is true if the deposit is within the stored AbsVeto of any item of the source collection

Implements reco::isodeposit::AbsVeto.

Definition at line 42 of file EventDependentAbsVetos.cc.

References reco::deltaR2(), PVValHelper::eta, and ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it.

42  {
43  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
44  if (::deltaR2(it->eta(), it->phi(), eta, phi) < dR2constituent_)
45  return true;
46  }
47  return false;
48 }
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16

Member Data Documentation

◆ dR2constituent_

double reco::isodeposit::OtherJetConstituentsDeltaRVeto::dR2constituent_
private

Definition at line 112 of file EventDependentAbsVetos.h.

◆ dR2jet_

double reco::isodeposit::OtherJetConstituentsDeltaRVeto::dR2jet_
private

Definition at line 110 of file EventDependentAbsVetos.h.

◆ evt_

const edm::Event* reco::isodeposit::OtherJetConstituentsDeltaRVeto::evt_
private

Definition at line 106 of file EventDependentAbsVetos.h.

◆ items_

std::vector<Direction> reco::isodeposit::OtherJetConstituentsDeltaRVeto::items_
private

Definition at line 113 of file EventDependentAbsVetos.h.

◆ srcJets_

edm::EDGetTokenT<reco::PFJetCollection> reco::isodeposit::OtherJetConstituentsDeltaRVeto::srcJets_
private

Definition at line 109 of file EventDependentAbsVetos.h.

◆ srcPFCandAssocMap_

edm::EDGetTokenT<JetToPFCandidateAssociation> reco::isodeposit::OtherJetConstituentsDeltaRVeto::srcPFCandAssocMap_
private

Definition at line 111 of file EventDependentAbsVetos.h.

◆ vetoDir_

Direction reco::isodeposit::OtherJetConstituentsDeltaRVeto::vetoDir_
private

Definition at line 108 of file EventDependentAbsVetos.h.