CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Types | Private Attributes
HLTDoublet< T1, T2 > Class Template Reference

#include <HLTDoublet.h>

Inheritance diagram for HLTDoublet< T1, T2 >:
HLTFilter edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 HLTDoublet (const edm::ParameterSet &)
 
virtual bool hltFilter (edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)
 
 ~HLTDoublet ()
 
- Public Member Functions inherited from HLTFilter
 HLTFilter (const edm::ParameterSet &config)
 
int module () const
 
const std::string * moduleLabel () const
 
int path () const
 
const std::string * pathName () const
 
std::pair< int, int > pmid () const
 
bool saveTags () const
 
virtual ~HLTFilter ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from HLTFilter
static void makeHLTFilterDescription (edm::ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 

Private Types

typedef std::vector< T1 > T1Collection
 
typedef edm::Ref< T1CollectionT1Ref
 
typedef std::vector< T2 > T2Collection
 
typedef edm::Ref< T2CollectionT2Ref
 

Private Attributes

std::vector< T1Refcoll1_
 
std::vector< T2Refcoll2_
 
bool cutdelr_
 
bool cutdeta_
 
bool cutdphi_
 
bool cutminv_
 
bool cutpt_
 
edm::InputTag inputTag1_
 
edm::InputTag inputTag2_
 
std::string label_
 
double max_DelR_
 
double max_Deta_
 
double max_Dphi_
 
double max_Minv_
 
double max_Pt_
 
double min_DelR_
 
double min_Deta_
 
double min_Dphi_
 
double min_Minv_
 
int min_N_
 
double min_Pt_
 
edm::InputTag originTag1_
 
edm::InputTag originTag2_
 
bool same_
 
int triggerType1_
 
int triggerType2_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
typedef WorkerT< EDFilterWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Protected Member Functions inherited from edm::EDFilter
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

template<typename T1, typename T2>
class HLTDoublet< T1, T2 >

This class is an HLTFilter (-> EDFilter) implementing a basic HLT trigger for pairs of object, evaluating all pairs with the first object from collection 1, and the second object from collection 2, cutting on variables relating to their 4-momentum representations. The object collections are assumed to be outputs of HLTSinglet single-object-type filters so that the access is thorugh RefToBases and polymorphic.

Date:
2012/02/24 13:34:20
Revision:
1.9
Author
Martin Grunewald

See header file for documentation

Date:
2012/03/27 18:04:53
Revision:
1.25
Author
Martin Grunewald

Definition at line 33 of file HLTDoublet.h.

Member Typedef Documentation

template<typename T1 , typename T2 >
typedef std::vector<T1> HLTDoublet< T1, T2 >::T1Collection
private

Definition at line 64 of file HLTDoublet.h.

template<typename T1 , typename T2 >
typedef edm::Ref<T1Collection> HLTDoublet< T1, T2 >::T1Ref
private

Definition at line 65 of file HLTDoublet.h.

template<typename T1 , typename T2 >
typedef std::vector<T2> HLTDoublet< T1, T2 >::T2Collection
private

Definition at line 67 of file HLTDoublet.h.

template<typename T1 , typename T2 >
typedef edm::Ref<T2Collection> HLTDoublet< T1, T2 >::T2Ref
private

Definition at line 68 of file HLTDoublet.h.

Constructor & Destructor Documentation

template<typename T1 , typename T2 >
HLTDoublet< T1, T2 >::HLTDoublet ( const edm::ParameterSet iConfig)
explicit

Definition at line 30 of file HLTDoublet.cc.

References HLTDoublet< T1, T2 >::cutdelr_, HLTDoublet< T1, T2 >::cutdeta_, HLTDoublet< T1, T2 >::cutdphi_, HLTDoublet< T1, T2 >::cutminv_, HLTDoublet< T1, T2 >::cutpt_, edm::InputTag::encode(), HLTDoublet< T1, T2 >::inputTag1_, HLTDoublet< T1, T2 >::inputTag2_, LogDebug, HLTDoublet< T1, T2 >::max_DelR_, HLTDoublet< T1, T2 >::max_Deta_, HLTDoublet< T1, T2 >::max_Dphi_, HLTDoublet< T1, T2 >::max_Minv_, HLTDoublet< T1, T2 >::max_Pt_, HLTDoublet< T1, T2 >::min_DelR_, HLTDoublet< T1, T2 >::min_Deta_, HLTDoublet< T1, T2 >::min_Dphi_, HLTDoublet< T1, T2 >::min_Minv_, HLTDoublet< T1, T2 >::min_N_, HLTDoublet< T1, T2 >::min_Pt_, HLTDoublet< T1, T2 >::same_, HLTDoublet< T1, T2 >::triggerType1_, and HLTDoublet< T1, T2 >::triggerType2_.

30  : HLTFilter(iConfig),
31  originTag1_(iConfig.template getParameter<edm::InputTag>("originTag1")),
32  originTag2_(iConfig.template getParameter<edm::InputTag>("originTag2")),
33  inputTag1_(iConfig.template getParameter<edm::InputTag>("inputTag1")),
34  inputTag2_(iConfig.template getParameter<edm::InputTag>("inputTag2")),
35  triggerType1_(iConfig.template getParameter<int>("triggerType1")),
36  triggerType2_(iConfig.template getParameter<int>("triggerType2")),
37  min_Dphi_ (iConfig.template getParameter<double>("MinDphi")),
38  max_Dphi_ (iConfig.template getParameter<double>("MaxDphi")),
39  min_Deta_ (iConfig.template getParameter<double>("MinDeta")),
40  max_Deta_ (iConfig.template getParameter<double>("MaxDeta")),
41  min_Minv_ (iConfig.template getParameter<double>("MinMinv")),
42  max_Minv_ (iConfig.template getParameter<double>("MaxMinv")),
43  min_DelR_ (iConfig.template getParameter<double>("MinDelR")),
44  max_DelR_ (iConfig.template getParameter<double>("MaxDelR")),
45  min_Pt_ (iConfig.template getParameter<double>("MinPt")),
46  max_Pt_ (iConfig.template getParameter<double>("MaxPt")),
47  min_N_ (iConfig.template getParameter<int>("MinN")),
48  label_ (iConfig.getParameter<std::string>("@module_label")),
49  coll1_(),
50  coll2_()
51 {
52 
53  // same collections to be compared?
55 
56  cutdphi_ = (min_Dphi_ <= max_Dphi_); // cut active?
57  cutdeta_ = (min_Deta_ <= max_Deta_); // cut active?
58  cutminv_ = (min_Minv_ <= max_Minv_); // cut active?
59  cutdelr_ = (min_DelR_ <= max_DelR_); // cut active?
60  cutpt_ = (min_Pt_ <= max_Pt_ ); // cut active?
61 
62  LogDebug("") << "InputTags and cuts : "
63  << inputTag1_.encode() << " " << inputTag2_.encode()
64  << triggerType1_ << " " << triggerType2_
65  << " Dphi [" << min_Dphi_ << " " << max_Dphi_ << "]"
66  << " Deta [" << min_Deta_ << " " << max_Deta_ << "]"
67  << " Minv [" << min_Minv_ << " " << max_Minv_ << "]"
68  << " DelR [" << min_DelR_ << " " << max_DelR_ << "]"
69  << " Pt [" << min_Pt_ << " " << max_Pt_ << "]"
70  << " MinN =" << min_N_
71  << " same/dphi/deta/minv/delr/pt "
72  << same_
73  << cutdphi_ << cutdeta_ << cutminv_ << cutdelr_ << cutpt_;
74 }
#define LogDebug(id)
double min_Pt_
Definition: HLTDoublet.h:54
T getParameter(std::string const &) const
edm::InputTag inputTag1_
Definition: HLTDoublet.h:46
double max_Dphi_
Definition: HLTDoublet.h:50
bool cutpt_
Definition: HLTDoublet.h:59
int triggerType2_
Definition: HLTDoublet.h:49
double max_Deta_
Definition: HLTDoublet.h:51
std::string label_
Definition: HLTDoublet.h:61
edm::InputTag originTag1_
Definition: HLTDoublet.h:44
double min_DelR_
Definition: HLTDoublet.h:53
edm::InputTag originTag2_
Definition: HLTDoublet.h:45
std::string encode() const
Definition: InputTag.cc:72
int triggerType1_
Definition: HLTDoublet.h:48
double min_Dphi_
Definition: HLTDoublet.h:50
bool cutminv_
Definition: HLTDoublet.h:59
double max_DelR_
Definition: HLTDoublet.h:53
double min_Deta_
Definition: HLTDoublet.h:51
double min_Minv_
Definition: HLTDoublet.h:52
bool same_
Definition: HLTDoublet.h:58
bool cutdelr_
Definition: HLTDoublet.h:59
HLTFilter(const edm::ParameterSet &config)
Definition: HLTFilter.cc:18
double max_Pt_
Definition: HLTDoublet.h:54
std::vector< T2Ref > coll2_
Definition: HLTDoublet.h:69
int min_N_
Definition: HLTDoublet.h:55
std::vector< T1Ref > coll1_
Definition: HLTDoublet.h:66
double max_Minv_
Definition: HLTDoublet.h:52
bool cutdphi_
Definition: HLTDoublet.h:59
edm::InputTag inputTag2_
Definition: HLTDoublet.h:47
bool cutdeta_
Definition: HLTDoublet.h:59
template<typename T1 , typename T2 >
HLTDoublet< T1, T2 >::~HLTDoublet ( )

Definition at line 77 of file HLTDoublet.cc.

78 {
79 }

Member Function Documentation

template<typename T1 , typename T2 >
void HLTDoublet< T1, T2 >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 82 of file HLTDoublet.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and mergeVDriftHistosByStation::name.

82  {
85  desc.add<edm::InputTag>("originTag1",edm::InputTag("hltOriginal1"));
86  desc.add<edm::InputTag>("originTag2",edm::InputTag("hltOriginal2"));
87  desc.add<edm::InputTag>("inputTag1",edm::InputTag("hltFiltered1"));
88  desc.add<edm::InputTag>("inputTag2",edm::InputTag("hltFiltered22"));
89  desc.add<int>("triggerType1",0);
90  desc.add<int>("triggerType2",0);
91  desc.add<double>("MinDphi",+1.0);
92  desc.add<double>("MaxDphi",-1.0);
93  desc.add<double>("MinDeta",+1.0);
94  desc.add<double>("MaxDeta",-1.0);
95  desc.add<double>("MinMinv",+1.0);
96  desc.add<double>("MaxMinv",-1.0);
97  desc.add<double>("MinDelR",+1.0);
98  desc.add<double>("MaxDelR",-1.0);
99  desc.add<double>("MinPt" ,+1.0);
100  desc.add<double>("MaxPt" ,-1.0);
101  desc.add<int>("MinN",1);
102  descriptions.add(std::string("hlt")+std::string(typeid(HLTDoublet<T1,T2>).name()),desc);
103 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:27
void add(std::string const &label, ParameterSetDescription const &psetDescription)
template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::hltFilter ( edm::Event iEvent,
const edm::EventSetup iSetup,
trigger::TriggerFilterObjectWithRefs filterproduct 
)
virtual

Implements HLTFilter.

Definition at line 112 of file HLTDoublet.cc.

References abs, accept(), trigger::TriggerFilterObjectWithRefs::addCollectionTag(), trigger::TriggerRefsCollections::addObject(), edm::HandleBase::clear(), edm::InputTag::encode(), edm::Event::getByLabel(), edm::Event::getProvenance(), Exhume::I, instance, diffTwoXMLs::label, M_PI, edm::Provenance::moduleLabel(), n, AlCaHLTBitMon_ParallelJobs::p, p1, p2, evf::utils::pid, LaserDQM_cfg::process, edm::Provenance::processName(), edm::Provenance::productInstanceName(), reco::tau::disc::Pt(), diffTwoXMLs::r1, diffTwoXMLs::r2, dt_dqm_sourceclient_common_cff::reco, CommPDSkim_cfg::saveTags, and mathSSE::sqrt().

113 {
114  using namespace std;
115  using namespace edm;
116  using namespace reco;
117  using namespace trigger;
118 
119  // All HLT filters must create and fill an HLT filter object,
120  // recording any reconstructed physics objects satisfying (or not)
121  // this HLT filter, and place it in the Event.
122 
123  bool accept(false);
124 
125  LogVerbatim("HLTDoublet") << " XXX " << label_ << " 0 " << std::endl;
126 
127  // get hold of pre-filtered object collections
129  if (iEvent.getByLabel (inputTag1_,coll1) && iEvent.getByLabel (inputTag2_,coll2)) {
130  coll1_.clear();
131  coll1->getObjects(triggerType1_,coll1_);
132  const size_type n1(coll1_.size());
133  coll2_.clear();
134  coll2->getObjects(triggerType2_,coll2_);
135  const size_type n2(coll2_.size());
136 
137  if (saveTags()) {
138  InputTag tagOld;
139  filterproduct.addCollectionTag(originTag1_);
140  LogVerbatim("HLTDoublet") << " XXX " << label_ << " 1a " << originTag1_.encode() << std::endl;
141  tagOld=InputTag();
142  for (size_type i1=0; i1!=n1; ++i1) {
143  const ProductID pid(coll1_[i1].id());
144  const string& label(iEvent.getProvenance(pid).moduleLabel());
145  const string& instance(iEvent.getProvenance(pid).productInstanceName());
146  const string& process(iEvent.getProvenance(pid).processName());
148  if (tagOld.encode()!=tagNew.encode()) {
149  filterproduct.addCollectionTag(tagNew);
150  tagOld=tagNew;
151  LogVerbatim("HLTDoublet") << " XXX " << label_ << " 1b " << tagNew.encode() << std::endl;
152  }
153  }
154  filterproduct.addCollectionTag(originTag2_);
155  LogVerbatim("HLTDoublet") << " XXX " << label_ << " 2a " << originTag2_.encode() << std::endl;
156  tagOld=InputTag();
157  for (size_type i2=0; i2!=n2; ++i2) {
158  const ProductID pid(coll2_[i2].id());
159  const string& label(iEvent.getProvenance(pid).moduleLabel());
160  const string& instance(iEvent.getProvenance(pid).productInstanceName());
161  const string& process(iEvent.getProvenance(pid).processName());
163  if (tagOld.encode()!=tagNew.encode()) {
164  filterproduct.addCollectionTag(tagNew);
165  tagOld=tagNew;
166  LogVerbatim("HLTDoublet") << " XXX " << label_ << " 2b " << tagNew.encode() << std::endl;
167  }
168  }
169  }
170 
171  int n(0);
172  T1Ref r1;
173  T2Ref r2;
175  for (unsigned int i1=0; i1!=n1; i1++) {
176  r1=coll1_[i1];
177  p1=r1->p4();
178  unsigned int I(0);
179  if (same_) {I=i1+1;}
180  for (unsigned int i2=I; i2!=n2; i2++) {
181  r2=coll2_[i2];
182  p2=r2->p4();
183 
184  double Dphi(std::abs(p1.phi()-p2.phi()));
185  if (Dphi>M_PI) Dphi=2.0*M_PI-Dphi;
186 
187  double Deta(std::abs(p1.eta()-p2.eta()));
188 
189  p=p1+p2;
190  double Minv(std::abs(p.mass()));
191  double Pt(p.pt());
192 
193  double DelR(sqrt(Dphi*Dphi+Deta*Deta));
194 
195  if ( ( (!cutdphi_) || ((min_Dphi_<=Dphi) && (Dphi<=max_Dphi_)) ) &&
196  ( (!cutdeta_) || ((min_Deta_<=Deta) && (Deta<=max_Deta_)) ) &&
197  ( (!cutminv_) || ((min_Minv_<=Minv) && (Minv<=max_Minv_)) ) &&
198  ( (!cutdelr_) || ((min_DelR_<=DelR) && (DelR<=max_DelR_)) ) &&
199  ( (!cutpt_ ) || ((min_Pt_ <=Pt ) && (Pt <=max_Pt_ )) ) ) {
200  n++;
201  filterproduct.addObject(triggerType1_,r1);
202  filterproduct.addObject(triggerType2_,r2);
203  }
204 
205  }
206  }
207  // filter decision
208  accept = (n>=min_N_);
209  }
210 
211  return accept;
212 }
double min_Pt_
Definition: HLTDoublet.h:54
edm::InputTag inputTag1_
Definition: HLTDoublet.h:46
double max_Dphi_
Definition: HLTDoublet.h:50
bool cutpt_
Definition: HLTDoublet.h:59
static PFTauRenderPlugin instance
int triggerType2_
Definition: HLTDoublet.h:49
edm::Ref< T2Collection > T2Ref
Definition: HLTDoublet.h:68
double max_Deta_
Definition: HLTDoublet.h:51
#define abs(x)
Definition: mlp_lapack.h:159
std::string label_
Definition: HLTDoublet.h:61
std::string const & processName() const
Definition: Provenance.h:63
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
edm::InputTag originTag1_
Definition: HLTDoublet.h:44
uint16_t size_type
double min_DelR_
Definition: HLTDoublet.h:53
edm::InputTag originTag2_
Definition: HLTDoublet.h:45
std::string encode() const
Definition: InputTag.cc:72
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
int triggerType1_
Definition: HLTDoublet.h:48
double min_Dphi_
Definition: HLTDoublet.h:50
bool cutminv_
Definition: HLTDoublet.h:59
edm::Ref< T1Collection > T1Ref
Definition: HLTDoublet.h:65
T sqrt(T t)
Definition: SSEVec.h:46
double max_DelR_
Definition: HLTDoublet.h:53
double min_Deta_
Definition: HLTDoublet.h:51
double min_Minv_
Definition: HLTDoublet.h:52
const std::complex< double > I
Definition: I.h:8
bool same_
Definition: HLTDoublet.h:58
double p2[4]
Definition: TauolaWrapper.h:90
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
bool cutdelr_
Definition: HLTDoublet.h:59
#define M_PI
Definition: BFit3D.cc:3
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
std::string const & moduleLabel() const
Definition: Provenance.h:62
double max_Pt_
Definition: HLTDoublet.h:54
bool saveTags() const
Definition: HLTFilter.h:45
double p1[4]
Definition: TauolaWrapper.h:89
std::vector< T2Ref > coll2_
Definition: HLTDoublet.h:69
int min_N_
Definition: HLTDoublet.h:55
std::vector< T1Ref > coll1_
Definition: HLTDoublet.h:66
std::string const & productInstanceName() const
Definition: Provenance.h:64
double max_Minv_
Definition: HLTDoublet.h:52
tuple process
Definition: LaserDQM_cfg.py:3
Provenance getProvenance(BranchID const &theID) const
Definition: Event.cc:60
bool cutdphi_
Definition: HLTDoublet.h:59
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
edm::InputTag inputTag2_
Definition: HLTDoublet.h:47
bool cutdeta_
Definition: HLTDoublet.h:59

Member Data Documentation

template<typename T1 , typename T2 >
std::vector<T1Ref> HLTDoublet< T1, T2 >::coll1_
private

Definition at line 66 of file HLTDoublet.h.

template<typename T1 , typename T2 >
std::vector<T2Ref> HLTDoublet< T1, T2 >::coll2_
private

Definition at line 69 of file HLTDoublet.h.

template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::cutdelr_
private

Definition at line 59 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::cutdeta_
private

Definition at line 59 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::cutdphi_
private

Definition at line 59 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::cutminv_
private

Definition at line 59 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::cutpt_
private

Definition at line 59 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
edm::InputTag HLTDoublet< T1, T2 >::inputTag1_
private

Definition at line 46 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
edm::InputTag HLTDoublet< T1, T2 >::inputTag2_
private

Definition at line 47 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
std::string HLTDoublet< T1, T2 >::label_
private
template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::max_DelR_
private

Definition at line 53 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::max_Deta_
private

Definition at line 51 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::max_Dphi_
private

Definition at line 50 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::max_Minv_
private

Definition at line 52 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::max_Pt_
private

Definition at line 54 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::min_DelR_
private

Definition at line 53 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::min_Deta_
private

Definition at line 51 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::min_Dphi_
private

Definition at line 50 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::min_Minv_
private

Definition at line 52 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
int HLTDoublet< T1, T2 >::min_N_
private

Definition at line 55 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
double HLTDoublet< T1, T2 >::min_Pt_
private

Definition at line 54 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
edm::InputTag HLTDoublet< T1, T2 >::originTag1_
private

Definition at line 44 of file HLTDoublet.h.

template<typename T1 , typename T2 >
edm::InputTag HLTDoublet< T1, T2 >::originTag2_
private

Definition at line 45 of file HLTDoublet.h.

template<typename T1 , typename T2 >
bool HLTDoublet< T1, T2 >::same_
private

Definition at line 58 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
int HLTDoublet< T1, T2 >::triggerType1_
private

Definition at line 48 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().

template<typename T1 , typename T2 >
int HLTDoublet< T1, T2 >::triggerType2_
private

Definition at line 49 of file HLTDoublet.h.

Referenced by HLTDoublet< T1, T2 >::HLTDoublet().