CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHemiDPhiFilter.cc
Go to the documentation of this file.
3 
5 
9 
12 
13 
16 
20 
21 //
22 // constructors and destructor
23 //
25  inputTag_ (iConfig.getParameter<edm::InputTag>("inputTag")),
26  min_dphi_ (iConfig.getParameter<double> ("minDPhi" )),
27  accept_NJ_ (iConfig.getParameter<bool> ("acceptNJ" ))
28 
29 {
30  m_theHemiToken = consumes<std::vector<math::XYZTLorentzVector>>(inputTag_);
31  LogDebug("") << "Inputs/minDphi/acceptNJ : "
32  << inputTag_.encode() << " "
33  << min_dphi_ << " "
34  << accept_NJ_ << ".";
35 }
36 
38 {
39 }
40 
41 void
45  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltRHemisphere"));
46  desc.add<double>("minDPhi",2.9415);
47  desc.add<bool>("acceptNJ",true);
48  descriptions.add("hltHemiDPhiFilter",desc);
49 }
50 
51 //
52 // member functions
53 //
54 
55 // ------------ method called to produce the data ------------
56 bool
58 {
59  using namespace std;
60  using namespace edm;
61  using namespace reco;
62  using namespace trigger;
63 
64  // get hold of collection of objects
67 
68  // check the the input collections are available
69  if (not hemispheres.isValid())
70  return false;
71 
72  if(hemispheres->size() ==0){ // the Hemisphere Maker will produce an empty collection of hemispheres if the number of jets in the
73  return accept_NJ_; // event is greater than the maximum number of jets
74  }
75 
76 
77 
78 
79  //***********************************
80  // Get the set of hemisphere axes
81 
82  TLorentzVector j1R(hemispheres->at(0).x(),hemispheres->at(0).y(),hemispheres->at(0).z(),hemispheres->at(0).t());
83  TLorentzVector j2R(hemispheres->at(1).x(),hemispheres->at(1).y(),hemispheres->at(1).z(),hemispheres->at(1).t());
84 
85  // compute the dPhi between them
86  double dphi = 50.;
87  dphi = deltaPhi(j1R.Phi(),j2R.Phi());
88 
89  // Dphi requirement
90 
91  if(dphi<=min_dphi_) return true;
92 
93  // filter decision
94  return false;
95 }
96 
97 
98 double HLTHemiDPhiFilter::deltaPhi(double v1, double v2)
99 {
100  // Computes the correctly normalized phi difference
101  // v1, v2 = phi of object 1 and 2
102  double diff = std::abs(v2 - v1);
103  return (diff < M_PI) ? diff : 2*M_PI - diff;
104 }
105 
#define LogDebug(id)
edm::InputTag inputTag_
HLTHemiDPhiFilter(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
std::string encode() const
Definition: InputTag.cc:164
int iEvent
Definition: GenABIO.cc:230
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< std::vector< math::XYZTLorentzVector > > m_theHemiToken
static double deltaPhi(double, double)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define M_PI
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void add(std::string const &label, ParameterSetDescription const &psetDescription)