CMS 3D CMS Logo

HLTDoublet.cc
Go to the documentation of this file.
1 
12 
14 
17 
19 
21 
23 #include <cmath>
24 
25 //
26 // constructors and destructor
27 //
28 template<typename T1, typename T2>
30  originTag1_(iConfig.template getParameter<std::vector<edm::InputTag> >("originTag1")),
31  originTag2_(iConfig.template getParameter<std::vector<edm::InputTag> >("originTag2")),
32  inputTag1_(iConfig.template getParameter<edm::InputTag>("inputTag1")),
33  inputTag2_(iConfig.template getParameter<edm::InputTag>("inputTag2")),
34  inputToken1_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag1_)),
35  inputToken2_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag2_)),
36  triggerType1_(iConfig.template getParameter<int>("triggerType1")),
37  triggerType2_(iConfig.template getParameter<int>("triggerType2")),
38  min_Dphi_ (iConfig.template getParameter<double>("MinDphi")),
39  max_Dphi_ (iConfig.template getParameter<double>("MaxDphi")),
40  min_Deta_ (iConfig.template getParameter<double>("MinDeta")),
41  max_Deta_ (iConfig.template getParameter<double>("MaxDeta")),
42  min_Minv_ (iConfig.template getParameter<double>("MinMinv")),
43  max_Minv_ (iConfig.template getParameter<double>("MaxMinv")),
44  min_DelR_ (iConfig.template getParameter<double>("MinDelR")),
45  max_DelR_ (iConfig.template getParameter<double>("MaxDelR")),
46  min_Pt_ (iConfig.template getParameter<double>("MinPt")),
47  max_Pt_ (iConfig.template getParameter<double>("MaxPt")),
48  min_N_ (iConfig.template getParameter<int>("MinN")),
49  same_ (inputTag1_.encode() == inputTag2_.encode()), // same collections to be compared?
50  cutdphi_ (min_Dphi_ <= max_Dphi_), // cut active?
51  cutdeta_ (min_Deta_ <= max_Deta_), // cut active?
52  cutminv_ (min_Minv_ <= max_Minv_), // cut active?
53  cutdelr_ (min_DelR_ <= max_DelR_), // cut active?
54  cutpt_ (min_Pt_ <= max_Pt_ ) // cut active?
55 {
56  LogDebug("") << "InputTags and cuts : "
57  << inputTag1_.encode() << " " << inputTag2_.encode()
58  << triggerType1_ << " " << triggerType2_
59  << " Dphi [" << min_Dphi_ << " " << max_Dphi_ << "]"
60  << " Deta [" << min_Deta_ << " " << max_Deta_ << "]"
61  << " Minv [" << min_Minv_ << " " << max_Minv_ << "]"
62  << " DelR [" << min_DelR_ << " " << max_DelR_ << "]"
63  << " Pt [" << min_Pt_ << " " << max_Pt_ << "]"
64  << " MinN =" << min_N_
65  << " same/dphi/deta/minv/delr/pt "
66  << same_
67  << cutdphi_ << cutdeta_ << cutminv_ << cutdelr_ << cutpt_;
68 }
69 
70 template<typename T1, typename T2>
72 template<typename T1, typename T2>
73 void
77  std::vector<edm::InputTag> originTag1(1,edm::InputTag("hltOriginal1"));
78  std::vector<edm::InputTag> originTag2(1,edm::InputTag("hltOriginal2"));
79  desc.add<std::vector<edm::InputTag> >("originTag1",originTag1);
80  desc.add<std::vector<edm::InputTag> >("originTag2",originTag2);
81  desc.add<edm::InputTag>("inputTag1",edm::InputTag("hltFiltered1"));
82  desc.add<edm::InputTag>("inputTag2",edm::InputTag("hltFiltered22"));
83  desc.add<int>("triggerType1",0);
84  desc.add<int>("triggerType2",0);
85  desc.add<double>("MinDphi",+1.0);
86  desc.add<double>("MaxDphi",-1.0);
87  desc.add<double>("MinDeta",+1.0);
88  desc.add<double>("MaxDeta",-1.0);
89  desc.add<double>("MinMinv",+1.0);
90  desc.add<double>("MaxMinv",-1.0);
91  desc.add<double>("MinDelR",+1.0);
92  desc.add<double>("MaxDelR",-1.0);
93  desc.add<double>("MinPt" ,+1.0);
94  desc.add<double>("MaxPt" ,-1.0);
95  desc.add<int>("MinN",1);
96  descriptions.add(defaultModuleLabel<HLTDoublet<T1,T2>>(), desc);
97 }
98 
99 //
100 // member functions
101 //
102 
103 // ------------ method called to produce the data ------------
104 template<typename T1, typename T2>
105 bool
107 {
108  using namespace std;
109  using namespace edm;
110  using namespace reco;
111  using namespace trigger;
112 
113  // All HLT filters must create and fill an HLT filter object,
114  // recording any reconstructed physics objects satisfying (or not)
115  // this HLT filter, and place it in the Event.
116 
117  bool accept(false);
118 
119  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 0 " << std::endl;
120 
121  std::vector<T1Ref> coll1;
122  std::vector<T2Ref> coll2;
123 
124  // get hold of pre-filtered object collections
125  Handle<TriggerFilterObjectWithRefs> handle1, handle2;
126  if (iEvent.getByToken(inputToken1_,handle1) && iEvent.getByToken(inputToken2_,handle2)) {
127  handle1->getObjects(triggerType1_, coll1);
128  handle2->getObjects(triggerType2_, coll2);
129  const size_type n1(coll1.size());
130  const size_type n2(coll2.size());
131 
132  if (saveTags()) {
133  InputTag tagOld;
134  for (unsigned int i=0; i<originTag1_.size(); ++i) {
135  filterproduct.addCollectionTag(originTag1_[i]);
136  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 1a/" << i << " " << originTag1_[i].encode() << std::endl;
137  }
138  tagOld=InputTag();
139  for (size_type i1=0; i1!=n1; ++i1) {
140  const ProductID pid(coll1[i1].id());
141  const string& label(iEvent.getProvenance(pid).moduleLabel());
142  const string& instance(iEvent.getProvenance(pid).productInstanceName());
143  const string& process(iEvent.getProvenance(pid).processName());
145  if (tagOld.encode()!=tagNew.encode()) {
146  filterproduct.addCollectionTag(tagNew);
147  tagOld=tagNew;
148  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 1b " << tagNew.encode() << std::endl;
149  }
150  }
151  for (unsigned int i=0; i<originTag2_.size(); ++i) {
152  filterproduct.addCollectionTag(originTag2_[i]);
153  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 2a/" << i << " " << originTag2_[i].encode() << std::endl;
154  }
155  tagOld=InputTag();
156  for (size_type i2=0; i2!=n2; ++i2) {
157  const ProductID pid(coll2[i2].id());
158  const string& label(iEvent.getProvenance(pid).moduleLabel());
159  const string& instance(iEvent.getProvenance(pid).productInstanceName());
160  const string& process(iEvent.getProvenance(pid).processName());
162  if (tagOld.encode()!=tagNew.encode()) {
163  filterproduct.addCollectionTag(tagNew);
164  tagOld=tagNew;
165  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 2b " << tagNew.encode() << std::endl;
166  }
167  }
168  }
169 
170  int n(0);
171  T1Ref r1;
172  T2Ref r2;
174  for (unsigned int i1=0; i1!=n1; i1++) {
175  r1=coll1[i1];
176  p1=r1->p4();
177  unsigned int I(0);
178  if (same_) {I=i1+1;}
179  for (unsigned int i2=I; i2!=n2; i2++) {
180  r2=coll2[i2];
181  p2=r2->p4();
182 
183  double Dphi(std::abs(p1.phi()-p2.phi()));
184  if (Dphi>M_PI) Dphi=2.0*M_PI-Dphi;
185 
186  double Deta(std::abs(p1.eta()-p2.eta()));
187 
188  p=p1+p2;
189  double Minv(std::abs(p.mass()));
190  double Pt(p.pt());
191 
192  double DelR(sqrt(Dphi*Dphi+Deta*Deta));
193 
194  if ( ( (!cutdphi_) || ((min_Dphi_<=Dphi) && (Dphi<=max_Dphi_)) ) &&
195  ( (!cutdeta_) || ((min_Deta_<=Deta) && (Deta<=max_Deta_)) ) &&
196  ( (!cutminv_) || ((min_Minv_<=Minv) && (Minv<=max_Minv_)) ) &&
197  ( (!cutdelr_) || ((min_DelR_<=DelR) && (DelR<=max_DelR_)) ) &&
198  ( (!cutpt_ ) || ((min_Pt_ <=Pt ) && (Pt <=max_Pt_ )) ) ) {
199  n++;
200  filterproduct.addObject(triggerType1_,r1);
201  filterproduct.addObject(triggerType2_,r2);
202  }
203 
204  }
205  }
206  // filter decision
207  accept = (n>=min_N_);
208  }
209 
210  return accept;
211 }
#define LogDebug(id)
const edm::InputTag inputTag2_
Definition: HLTDoublet.h:48
std::string defaultModuleLabel()
const double min_Pt_
Definition: HLTDoublet.h:57
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
const double min_Deta_
Definition: HLTDoublet.h:54
HLTDoublet(const edm::ParameterSet &)
Definition: HLTDoublet.cc:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTDoublet.cc:74
const bool same_
Definition: HLTDoublet.h:61
const bool cutpt_
Definition: HLTDoublet.h:62
static PFTauRenderPlugin instance
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
const double min_Dphi_
Definition: HLTDoublet.h:53
const edm::InputTag inputTag1_
Definition: HLTDoublet.h:47
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken1_
Definition: HLTDoublet.h:49
const double min_Minv_
Definition: HLTDoublet.h:55
const double max_Deta_
Definition: HLTDoublet.h:54
std::string const & processName() const
Definition: Provenance.h:52
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
uint16_t size_type
std::string encode() const
Definition: InputTag.cc:165
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:230
const int triggerType2_
Definition: HLTDoublet.h:52
const bool cutdphi_
Definition: HLTDoublet.h:62
T sqrt(T t)
Definition: SSEVec.h:18
const double min_DelR_
Definition: HLTDoublet.h:56
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const bool cutdelr_
Definition: HLTDoublet.h:62
const std::complex< double > I
Definition: I.h:8
const int min_N_
Definition: HLTDoublet.h:58
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:520
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const double max_Minv_
Definition: HLTDoublet.h:55
const std::string * moduleLabel() const
Definition: HLTFilter.cc:66
double p2[4]
Definition: TauolaWrapper.h:90
#define M_PI
const double max_DelR_
Definition: HLTDoublet.h:56
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_
Definition: HLTDoublet.h:50
const double max_Pt_
Definition: HLTDoublet.h:57
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::string const & moduleLabel() const
Definition: Provenance.h:50
def encode(args, files)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
double p1[4]
Definition: TauolaWrapper.h:89
const std::vector< edm::InputTag > originTag2_
Definition: HLTDoublet.h:46
const std::vector< edm::InputTag > originTag1_
Definition: HLTDoublet.h:45
const int triggerType1_
Definition: HLTDoublet.h:51
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition: HLTDoublet.cc:106
const bool cutminv_
Definition: HLTDoublet.h:62
const double max_Dphi_
Definition: HLTDoublet.h:53
std::string const & productInstanceName() const
Definition: Provenance.h:53
Provenance getProvenance(BranchID const &theID) const
Definition: Event.cc:81
const bool cutdeta_
Definition: HLTDoublet.h:62
math::PtEtaPhiELorentzVectorF LorentzVector