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 
10 
11 bool
13 {
14  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
15  if (::deltaR2(it->eta(), it->phi(), eta, phi) < deltaR2_) return true;
16  }
17  return false;
18 }
19 
20 void
22  items_.clear();
24  iEvent.getByLabel(src_, candidates);
25  for (edm::View<reco::Candidate>::const_iterator it = candidates->begin(), ed = candidates->end(); it != ed; ++it) {
26  items_.push_back(Direction(it->eta(), it->phi()));
27  }
28 }
29 
30 bool
31 reco::isodeposit::OtherCandVeto::veto(double eta, double phi, float value) const
32 {
33  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
34  veto_->centerOn(it->eta(), it->phi());
35  if ( veto_->veto(eta,phi,value) ) return true;
36  }
37  return false;
38 }
39 
40 void
42  items_.clear();
44  iEvent.getByLabel(src_, candidates);
45  for (edm::View<reco::Candidate>::const_iterator it = candidates->begin(), ed = candidates->end(); it != ed; ++it) {
46  items_.push_back(Direction(it->eta(), it->phi()));
47  }
48 }
49 
50 bool
52 {
53  for (std::vector<Direction>::const_iterator it = items_.begin(), ed = items_.end(); it != ed; ++it) {
54  if (::deltaR2(it->eta(), it->phi(), eta, phi) < dR2constituent_) return true;
55  }
56  return false;
57 }
58 
59 void
61 {
62  //std::cout << "<OtherJetConstituentsDeltaRVeto::setEvent>:" << std::endl;
63  evt_ = &evt;
64 }
65 void
67 {
68  //std::cout << "<OtherJetConstituentsDeltaRVeto::initialize>:" << std::endl;
69  //std::cout << " vetoDir: eta = " << vetoDir_.eta() << ", phi = " << vetoDir_.phi() << std::endl;
70  assert(evt_);
71  items_.clear();
73  evt_->getByLabel(srcJets_, jets);
74  typedef edm::AssociationMap<edm::OneToMany<std::vector<reco::PFJet>, std::vector<reco::PFCandidate>, unsigned int> > JetToPFCandidateAssociation;
76  evt_->getByLabel(srcPFCandAssocMap_, jetToPFCandMap);
77  double dR2min = dR2jet_;
78  reco::PFJetRef matchedJet;
79  size_t numJets = jets->size();
80  for ( size_t jetIndex = 0; jetIndex < numJets; ++jetIndex ) {
81  reco::PFJetRef jet(jets, jetIndex);
82  double dR2 = ::deltaR2(vetoDir_.eta(), vetoDir_.phi(), jet->eta(), jet->phi());
83  //std::cout << "jet #" << jetIndex << ": Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->phi() << " (dR = " << sqrt(dR2) << ")" << std::endl;
84  if ( dR2 < dR2min ) {
85  matchedJet = jet;
86  dR2min = dR2;
87  }
88  }
89  if ( matchedJet.isNonnull() ) {
90  edm::RefVector<reco::PFCandidateCollection> pfCandsMappedToJet = (*jetToPFCandMap)[matchedJet];
91  int idx = 0;
92  for ( edm::RefVector<reco::PFCandidateCollection>::const_iterator pfCand = pfCandsMappedToJet.begin();
93  pfCand != pfCandsMappedToJet.end(); ++pfCand ) {
94  //std::cout << "pfCand #" << idx << ": Pt = " << (*pfCand)->pt() << ", eta = " << (*pfCand)->eta() << ", phi = " << (*pfCand)->phi() << std::endl;
95  items_.push_back(Direction((*pfCand)->eta(), (*pfCand)->phi()));
96  ++idx;
97  }
98  }
99 }
100 
101 void
103 {
104  //std::cout << "<OtherJetConstituentsDeltaRVeto::centerOn>:" << std::endl;
105  //std::cout << " eta = " << eta << std::endl;
106  //std::cout << " phi = " << phi << std::endl;
107  vetoDir_ = Direction(eta,phi);
108  initialize();
109 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
virtual void setEvent(const edm::Event &evt, const edm::EventSetup &es)
Picks up the directions of the given candidates.
virtual bool veto(double eta, double phi, float value) const
virtual bool veto(double eta, double phi, float value) const
virtual bool veto(double eta, double phi, float value) const
T eta() const
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
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:244
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
int iEvent
Definition: GenABIO.cc:243
vector< PseudoJet > jets
virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Picks up the directions of the given candidates.
double deltaR2(double eta1, double phi1, double eta2, double phi2)
Definition: deltaR.h:13
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Picks up the directions of the given candidates.
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
Definition: DDAxes.h:10