CMS 3D CMS Logo

EventDependentAbsVetos.cc
Go to the documentation of this file.
1 
6 
8  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
9  if (::deltaR2(it->eta(), it->phi(), eta, phi) < deltaR2_)
10  return true;
11  }
12  return false;
13 }
14 
16  items_.clear();
18  iEvent.getByToken(src_, candidates);
19  for (edm::View<reco::Candidate>::const_iterator it = candidates->begin(), ed = candidates->end(); it != ed; ++it) {
20  items_.push_back(Direction(it->eta(), it->phi()));
21  }
22 }
23 
24 bool reco::isodeposit::OtherCandVeto::veto(double eta, double phi, float value) const {
25  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
26  veto_->centerOn(it->eta(), it->phi());
27  if (veto_->veto(eta, phi, value))
28  return true;
29  }
30  return false;
31 }
32 
34  items_.clear();
36  iEvent.getByToken(src_, candidates);
37  for (edm::View<reco::Candidate>::const_iterator it = candidates->begin(), ed = candidates->end(); it != ed; ++it) {
38  items_.push_back(Direction(it->eta(), it->phi()));
39  }
40 }
41 
42 bool reco::isodeposit::OtherJetConstituentsDeltaRVeto::veto(double eta, double phi, float value) const {
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 }
49 
51  //std::cout << "<OtherJetConstituentsDeltaRVeto::setEvent>:" << std::endl;
52  evt_ = &evt;
53 }
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();
60  evt_->getByToken(srcJets_, jets);
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  int idx = 0;
78  for (edm::RefVector<reco::PFCandidateCollection>::const_iterator pfCand = pfCandsMappedToJet.begin();
79  pfCand != pfCandsMappedToJet.end();
80  ++pfCand) {
81  //std::cout << "pfCand #" << idx << ": Pt = " << (*pfCand)->pt() << ", eta = " << (*pfCand)->eta() << ", phi = " << (*pfCand)->phi() << std::endl;
82  items_.push_back(Direction((*pfCand)->eta(), (*pfCand)->phi()));
83  ++idx;
84  }
85  }
86 }
87 
89  //std::cout << "<OtherJetConstituentsDeltaRVeto::centerOn>:" << std::endl;
90  //std::cout << " eta = " << eta << std::endl;
91  //std::cout << " phi = " << phi << std::endl;
92  vetoDir_ = Direction(eta, phi);
93  initialize();
94 }
reco::isodeposit::OtherCandVeto::veto
bool veto(double eta, double phi, float value) const override
Definition: EventDependentAbsVetos.cc:24
PFCandidate.h
EventDependentAbsVetos.h
edm::RefVector::begin
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
cms::cuda::assert
assert(be >=bs)
reco::isodeposit::OtherJetConstituentsDeltaRVeto::initialize
void initialize()
Definition: EventDependentAbsVetos.cc:54
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
edm::RefVector
Definition: EDProductfwd.h:27
edm::Handle
Definition: AssociativeIterator.h:50
reco::isodeposit::OtherCandidatesDeltaRVeto::items_
std::vector< Direction > items_
Definition: EventDependentAbsVetos.h:36
edm::Ref< PFJetCollection >
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
edm::RefVector::end
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
reco::isodeposit::OtherCandidatesDeltaRVeto::deltaR2_
float deltaR2_
Definition: EventDependentAbsVetos.h:35
PVValHelper::eta
Definition: PVValidationHelpers.h:70
reco::isodeposit::OtherJetConstituentsDeltaRVeto::centerOn
void centerOn(double eta, double phi) override
Set axis for matching jets.
Definition: EventDependentAbsVetos.cc:88
reco::isodeposit::OtherJetConstituentsDeltaRVeto::setEvent
void setEvent(const edm::Event &evt, const edm::EventSetup &es) override
Picks up the directions of the given candidates.
Definition: EventDependentAbsVetos.cc:50
PVValHelper::phi
Definition: PVValidationHelpers.h:69
deltaR.h
edm::Ref::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
reco::isodeposit::OtherCandidatesDeltaRVeto::setEvent
void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
Picks up the directions of the given candidates.
Definition: EventDependentAbsVetos.cc:15
iEvent
int iEvent
Definition: GenABIO.cc:224
value
Definition: value.py:1
edm::EventSetup
Definition: EventSetup.h:58
reco::isodeposit::OtherCandidatesDeltaRVeto::veto
bool veto(double eta, double phi, float value) const override
Definition: EventDependentAbsVetos.cc:7
metsig::jet
Definition: SignAlgoResolutions.h:47
HLT_FULL_cff.candidates
candidates
Definition: HLT_FULL_cff.py:55017
reco::isodeposit::Direction
Definition: IsoDepositDirection.h:19
reco::isodeposit::OtherJetConstituentsDeltaRVeto::veto
bool veto(double eta, double phi, float value) const override
Definition: EventDependentAbsVetos.cc:42
edm::RefVectorIterator
Definition: EDProductfwd.h:33
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
initialize
static AlgebraicMatrix initialize()
Definition: BeamSpotTransientTrackingRecHit.cc:24
edm::Event
Definition: Event.h:73
reco::isodeposit::OtherCandVeto::setEvent
void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
Picks up the directions of the given candidates.
Definition: EventDependentAbsVetos.cc:33
PFCandidateFwd.h