CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EventDependentAbsVetos.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_IsolationAlgos_EventDependentAbsVetos_h
2 #define PhysicsTools_IsolationAlgos_EventDependentAbsVetos_h
3 
11 
12 namespace reco {
13  namespace isodeposit {
15  public:
18  src_(iC.consumes<edm::View<reco::Candidate> >(candidates)), deltaR2_(deltaR*deltaR) { }
19 
20  // Virtual destructor (should always be there)
22 
25  virtual bool veto(double eta, double phi, float value) const ;
26 
28  virtual void centerOn(double eta, double phi) { }
29 
31  virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) ;
32 
33  private:
35  float deltaR2_;
36  std::vector<Direction> items_;
37  };
38 
40  public:
43  src_(iC.consumes<edm::View<reco::Candidate> >(candidates)), veto_(veto) { }
44 
45  // Virtual destructor (should always be there)
46  virtual ~OtherCandVeto() {}
47 
50  virtual bool veto(double eta, double phi, float value) const ;
51 
53  virtual void centerOn(double eta, double phi) { }
54 
56  virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) ;
57 
58  private:
60  std::vector<Direction> items_;
61  std::unique_ptr<AbsVeto> veto_;
62  };
63 
65  public:
67  OtherJetConstituentsDeltaRVeto(Direction dir, const edm::InputTag& jets, double dRjet, const edm::InputTag& pfCandAssocMap, double dRconstituent, edm::ConsumesCollector& iC)
68  : evt_(0),
69  vetoDir_(dir),
70  srcJets_(iC.consumes<reco::PFJetCollection>(jets)),
71  dR2jet_(dRjet*dRjet),
72  srcPFCandAssocMap_(iC.consumes<JetToPFCandidateAssociation>(pfCandAssocMap)),
73  dR2constituent_(dRconstituent*dRconstituent)
74  {
75  //std::cout << "<OtherJetConstituentsDeltaRVeto::OtherJetConstituentsDeltaRVeto>:" << std::endl;
76  //std::cout << " vetoDir: eta = " << vetoDir_.eta() << ", phi = " << vetoDir_.phi() << std::endl;
77  //std::cout << " srcJets = " << srcJets_.label() << ":" << srcJets_.instance() << std::endl;
78  //std::cout << " dRjet = " << sqrt(dR2jet_) << std::endl;
79  //std::cout << " srcPFCandAssocMap = " << srcPFCandAssocMap_.label() << ":" << srcPFCandAssocMap_.instance() << std::endl;
80  //std::cout << " dRconstituent = " << sqrt(dR2constituent_) << std::endl;
81  }
82 
83  // Virtual destructor (should always be there)
85 
88  virtual bool veto(double eta, double phi, float value) const;
89 
91  virtual void centerOn(double eta, double phi);
92 
94  virtual void setEvent(const edm::Event& evt, const edm::EventSetup& es);
95 
96  private:
97  typedef edm::AssociationMap<edm::OneToMany<std::vector<reco::PFJet>, std::vector<reco::PFCandidate>, unsigned int> > JetToPFCandidateAssociation;
98 
99  void initialize();
100 
101  const edm::Event* evt_;
102 
105  double dR2jet_;
108  std::vector<Direction> items_;
109  };
110  }
111 }
112 #endif
edm::AssociationMap< edm::OneToMany< std::vector< reco::PFJet >, std::vector< reco::PFCandidate >, unsigned int > > JetToPFCandidateAssociation
virtual void centerOn(double eta, double phi)
Nothing to do for this.
edm::EDGetTokenT< edm::View< reco::Candidate > > src_
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
edm::EDGetTokenT< edm::View< reco::Candidate > > src_
virtual void centerOn(double eta, double phi)
Set axis for matching jets.
int iEvent
Definition: GenABIO.cc:230
vector< PseudoJet > jets
edm::EDGetTokenT< JetToPFCandidateAssociation > srcPFCandAssocMap_
virtual void centerOn(double eta, double phi)
Nothing to do for this.
virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Picks up the directions of the given candidates.
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
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...
virtual void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup)
Picks up the directions of the given candidates.
std::unique_ptr< AbsVeto > veto_
edm::EDGetTokenT< reco::PFJetCollection > srcJets_
OtherCandidatesDeltaRVeto(const edm::InputTag &candidates, double deltaR, edm::ConsumesCollector &iC)
Create a veto specifying the input collection of the candidates, and the deltaR.
std::vector< PFJet > PFJetCollection
collection of PFJet objects
dbl *** dir
Definition: mlp_gen.cc:35
OtherCandVeto(const edm::InputTag &candidates, AbsVeto *veto, edm::ConsumesCollector &iC)
Create a veto specifying the input collection of the candidates, and the deltaR.