CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EventDependentAbsVetos.cc
Go to the documentation of this file.
1 
6 
7 bool
9 {
10  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
11  if (::deltaR2(it->eta(), it->phi(), eta, phi) < deltaR2_) return true;
12  }
13  return false;
14 }
15 
16 void
18  items_.clear();
20  iEvent.getByToken(src_, candidates);
21  for (edm::View<reco::Candidate>::const_iterator it = candidates->begin(), ed = candidates->end(); it != ed; ++it) {
22  items_.push_back(Direction(it->eta(), it->phi()));
23  }
24 }
25 
26 bool
27 reco::isodeposit::OtherCandVeto::veto(double eta, double phi, float value) const
28 {
29  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
30  veto_->centerOn(it->eta(), it->phi());
31  if ( veto_->veto(eta,phi,value) ) return true;
32  }
33  return false;
34 }
35 
36 void
38  items_.clear();
40  iEvent.getByToken(src_, candidates);
41  for (edm::View<reco::Candidate>::const_iterator it = candidates->begin(), ed = candidates->end(); it != ed; ++it) {
42  items_.push_back(Direction(it->eta(), it->phi()));
43  }
44 }
45 
46 bool
48 {
49  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
50  if (::deltaR2(it->eta(), it->phi(), eta, phi) < dR2constituent_) return true;
51  }
52  return false;
53 }
54 
55 void
57 {
58  //std::cout << "<OtherJetConstituentsDeltaRVeto::setEvent>:" << std::endl;
59  evt_ = &evt;
60 }
61 void
63 {
64  //std::cout << "<OtherJetConstituentsDeltaRVeto::initialize>:" << std::endl;
65  //std::cout << " vetoDir: eta = " << vetoDir_.eta() << ", phi = " << vetoDir_.phi() << std::endl;
66  assert(evt_);
67  items_.clear();
69  evt_->getByToken(srcJets_, jets);
71  evt_->getByToken(srcPFCandAssocMap_, jetToPFCandMap);
72  double dR2min = dR2jet_;
73  reco::PFJetRef matchedJet;
74  size_t numJets = jets->size();
75  for ( size_t jetIndex = 0; jetIndex < numJets; ++jetIndex ) {
76  reco::PFJetRef jet(jets, jetIndex);
77  double dR2 = ::deltaR2(vetoDir_.eta(), vetoDir_.phi(), jet->eta(), jet->phi());
78  //std::cout << "jet #" << jetIndex << ": Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->phi() << " (dR = " << sqrt(dR2) << ")" << std::endl;
79  if ( dR2 < dR2min ) {
80  matchedJet = jet;
81  dR2min = dR2;
82  }
83  }
84  if ( matchedJet.isNonnull() ) {
85  edm::RefVector<reco::PFCandidateCollection> pfCandsMappedToJet = (*jetToPFCandMap)[matchedJet];
86  int idx = 0;
87  for ( edm::RefVector<reco::PFCandidateCollection>::const_iterator pfCand = pfCandsMappedToJet.begin();
88  pfCand != pfCandsMappedToJet.end(); ++pfCand ) {
89  //std::cout << "pfCand #" << idx << ": Pt = " << (*pfCand)->pt() << ", eta = " << (*pfCand)->eta() << ", phi = " << (*pfCand)->phi() << std::endl;
90  items_.push_back(Direction((*pfCand)->eta(), (*pfCand)->phi()));
91  ++idx;
92  }
93  }
94 }
95 
96 void
98 {
99  //std::cout << "<OtherJetConstituentsDeltaRVeto::centerOn>:" << std::endl;
100  //std::cout << " eta = " << eta << std::endl;
101  //std::cout << " phi = " << phi << std::endl;
102  vetoDir_ = Direction(eta,phi);
103  initialize();
104 }
static AlgebraicMatrix initialize()
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
virtual void setEvent(const edm::Event &evt, const edm::EventSetup &es)
Picks up the directions of the given candidates.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual bool veto(double eta, double phi, float value) const
virtual bool veto(double eta, double phi, float value) const
assert(m_qm.get())
virtual bool veto(double eta, double phi, float value) const
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
virtual void centerOn(double eta, double phi)
Set axis for matching jets.
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
int iEvent
Definition: GenABIO.cc:230
vector< PseudoJet > jets
virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Picks up the directions of the given candidates.
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Picks up the directions of the given candidates.
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
Geom::Phi< T > phi() const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85