#include <HLTDoublet.h>
Public Member Functions | |
HLTDoublet (const edm::ParameterSet &) | |
virtual bool | hltFilter (edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) |
~HLTDoublet () | |
Static Public Member Functions | |
static void | fillDescriptions (edm::ConfigurationDescriptions &descriptions) |
Private Types | |
typedef std::vector< T1 > | T1Collection |
typedef edm::Ref< T1Collection > | T1Ref |
typedef std::vector< T2 > | T2Collection |
typedef edm::Ref< T2Collection > | T2Ref |
Private Attributes | |
std::vector< T1Ref > | coll1_ |
std::vector< T2Ref > | coll2_ |
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_ |
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.
See header file for documentation
Definition at line 33 of file HLTDoublet.h.
typedef std::vector<T1> HLTDoublet< T1, T2 >::T1Collection [private] |
Definition at line 64 of file HLTDoublet.h.
typedef edm::Ref<T1Collection> HLTDoublet< T1, T2 >::T1Ref [private] |
Definition at line 65 of file HLTDoublet.h.
typedef std::vector<T2> HLTDoublet< T1, T2 >::T2Collection [private] |
Definition at line 67 of file HLTDoublet.h.
typedef edm::Ref<T2Collection> HLTDoublet< T1, T2 >::T2Ref [private] |
Definition at line 68 of file HLTDoublet.h.
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_.
: HLTFilter(iConfig), originTag1_(iConfig.template getParameter<edm::InputTag>("originTag1")), originTag2_(iConfig.template getParameter<edm::InputTag>("originTag2")), inputTag1_(iConfig.template getParameter<edm::InputTag>("inputTag1")), inputTag2_(iConfig.template getParameter<edm::InputTag>("inputTag2")), triggerType1_(iConfig.template getParameter<int>("triggerType1")), triggerType2_(iConfig.template getParameter<int>("triggerType2")), min_Dphi_ (iConfig.template getParameter<double>("MinDphi")), max_Dphi_ (iConfig.template getParameter<double>("MaxDphi")), min_Deta_ (iConfig.template getParameter<double>("MinDeta")), max_Deta_ (iConfig.template getParameter<double>("MaxDeta")), min_Minv_ (iConfig.template getParameter<double>("MinMinv")), max_Minv_ (iConfig.template getParameter<double>("MaxMinv")), min_DelR_ (iConfig.template getParameter<double>("MinDelR")), max_DelR_ (iConfig.template getParameter<double>("MaxDelR")), min_Pt_ (iConfig.template getParameter<double>("MinPt")), max_Pt_ (iConfig.template getParameter<double>("MaxPt")), min_N_ (iConfig.template getParameter<int>("MinN")), label_ (iConfig.getParameter<std::string>("@module_label")), coll1_(), coll2_() { // same collections to be compared? same_ = (inputTag1_.encode()==inputTag2_.encode()); cutdphi_ = (min_Dphi_ <= max_Dphi_); // cut active? cutdeta_ = (min_Deta_ <= max_Deta_); // cut active? cutminv_ = (min_Minv_ <= max_Minv_); // cut active? cutdelr_ = (min_DelR_ <= max_DelR_); // cut active? cutpt_ = (min_Pt_ <= max_Pt_ ); // cut active? LogDebug("") << "InputTags and cuts : " << inputTag1_.encode() << " " << inputTag2_.encode() << triggerType1_ << " " << triggerType2_ << " Dphi [" << min_Dphi_ << " " << max_Dphi_ << "]" << " Deta [" << min_Deta_ << " " << max_Deta_ << "]" << " Minv [" << min_Minv_ << " " << max_Minv_ << "]" << " DelR [" << min_DelR_ << " " << max_DelR_ << "]" << " Pt [" << min_Pt_ << " " << max_Pt_ << "]" << " MinN =" << min_N_ << " same/dphi/deta/minv/delr/pt " << same_ << cutdphi_ << cutdeta_ << cutminv_ << cutdelr_ << cutpt_; }
HLTDoublet< T1, T2 >::~HLTDoublet | ( | ) |
Definition at line 77 of file HLTDoublet.cc.
{ }
void HLTDoublet< T1, T2 >::fillDescriptions | ( | edm::ConfigurationDescriptions & | descriptions | ) | [static] |
Reimplemented from edm::EDFilter.
Definition at line 82 of file HLTDoublet.cc.
References edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::add(), and mergeVDriftHistosByStation::name.
{ edm::ParameterSetDescription desc; makeHLTFilterDescription(desc); desc.add<edm::InputTag>("originTag1",edm::InputTag("hltOriginal1")); desc.add<edm::InputTag>("originTag2",edm::InputTag("hltOriginal2")); desc.add<edm::InputTag>("inputTag1",edm::InputTag("hltFiltered1")); desc.add<edm::InputTag>("inputTag2",edm::InputTag("hltFiltered22")); desc.add<int>("triggerType1",0); desc.add<int>("triggerType2",0); desc.add<double>("MinDphi",+1.0); desc.add<double>("MaxDphi",-1.0); desc.add<double>("MinDeta",+1.0); desc.add<double>("MaxDeta",-1.0); desc.add<double>("MinMinv",+1.0); desc.add<double>("MaxMinv",-1.0); desc.add<double>("MinDelR",+1.0); desc.add<double>("MaxDelR",-1.0); desc.add<double>("MinPt" ,+1.0); desc.add<double>("MaxPt" ,-1.0); desc.add<int>("MinN",1); descriptions.add(std::string("hlt")+std::string(typeid(HLTDoublet<T1,T2>).name()),desc); }
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().
{ using namespace std; using namespace edm; using namespace reco; using namespace trigger; // All HLT filters must create and fill an HLT filter object, // recording any reconstructed physics objects satisfying (or not) // this HLT filter, and place it in the Event. bool accept(false); LogVerbatim("HLTDoublet") << " XXX " << label_ << " 0 " << std::endl; // get hold of pre-filtered object collections Handle<TriggerFilterObjectWithRefs> coll1,coll2; if (iEvent.getByLabel (inputTag1_,coll1) && iEvent.getByLabel (inputTag2_,coll2)) { coll1_.clear(); coll1->getObjects(triggerType1_,coll1_); const size_type n1(coll1_.size()); coll2_.clear(); coll2->getObjects(triggerType2_,coll2_); const size_type n2(coll2_.size()); if (saveTags()) { InputTag tagOld; filterproduct.addCollectionTag(originTag1_); LogVerbatim("HLTDoublet") << " XXX " << label_ << " 1a " << originTag1_.encode() << std::endl; tagOld=InputTag(); for (size_type i1=0; i1!=n1; ++i1) { const ProductID pid(coll1_[i1].id()); const string& label(iEvent.getProvenance(pid).moduleLabel()); const string& instance(iEvent.getProvenance(pid).productInstanceName()); const string& process(iEvent.getProvenance(pid).processName()); InputTag tagNew(InputTag(label,instance,process)); if (tagOld.encode()!=tagNew.encode()) { filterproduct.addCollectionTag(tagNew); tagOld=tagNew; LogVerbatim("HLTDoublet") << " XXX " << label_ << " 1b " << tagNew.encode() << std::endl; } } filterproduct.addCollectionTag(originTag2_); LogVerbatim("HLTDoublet") << " XXX " << label_ << " 2a " << originTag2_.encode() << std::endl; tagOld=InputTag(); for (size_type i2=0; i2!=n2; ++i2) { const ProductID pid(coll2_[i2].id()); const string& label(iEvent.getProvenance(pid).moduleLabel()); const string& instance(iEvent.getProvenance(pid).productInstanceName()); const string& process(iEvent.getProvenance(pid).processName()); InputTag tagNew(InputTag(label,instance,process)); if (tagOld.encode()!=tagNew.encode()) { filterproduct.addCollectionTag(tagNew); tagOld=tagNew; LogVerbatim("HLTDoublet") << " XXX " << label_ << " 2b " << tagNew.encode() << std::endl; } } } int n(0); T1Ref r1; T2Ref r2; Particle::LorentzVector p1,p2,p; for (unsigned int i1=0; i1!=n1; i1++) { r1=coll1_[i1]; p1=r1->p4(); unsigned int I(0); if (same_) {I=i1+1;} for (unsigned int i2=I; i2!=n2; i2++) { r2=coll2_[i2]; p2=r2->p4(); double Dphi(std::abs(p1.phi()-p2.phi())); if (Dphi>M_PI) Dphi=2.0*M_PI-Dphi; double Deta(std::abs(p1.eta()-p2.eta())); p=p1+p2; double Minv(std::abs(p.mass())); double Pt(p.pt()); double DelR(sqrt(Dphi*Dphi+Deta*Deta)); if ( ( (!cutdphi_) || ((min_Dphi_<=Dphi) && (Dphi<=max_Dphi_)) ) && ( (!cutdeta_) || ((min_Deta_<=Deta) && (Deta<=max_Deta_)) ) && ( (!cutminv_) || ((min_Minv_<=Minv) && (Minv<=max_Minv_)) ) && ( (!cutdelr_) || ((min_DelR_<=DelR) && (DelR<=max_DelR_)) ) && ( (!cutpt_ ) || ((min_Pt_ <=Pt ) && (Pt <=max_Pt_ )) ) ) { n++; filterproduct.addObject(triggerType1_,r1); filterproduct.addObject(triggerType2_,r2); } } } // filter decision accept = (n>=min_N_); } return accept; }
std::vector<T1Ref> HLTDoublet< T1, T2 >::coll1_ [private] |
Definition at line 66 of file HLTDoublet.h.
std::vector<T2Ref> HLTDoublet< T1, T2 >::coll2_ [private] |
Definition at line 69 of file HLTDoublet.h.
bool HLTDoublet< T1, T2 >::cutdelr_ [private] |
Definition at line 59 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
bool HLTDoublet< T1, T2 >::cutdeta_ [private] |
Definition at line 59 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
bool HLTDoublet< T1, T2 >::cutdphi_ [private] |
Definition at line 59 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
bool HLTDoublet< T1, T2 >::cutminv_ [private] |
Definition at line 59 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
bool HLTDoublet< T1, T2 >::cutpt_ [private] |
Definition at line 59 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
edm::InputTag HLTDoublet< T1, T2 >::inputTag1_ [private] |
Definition at line 46 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
edm::InputTag HLTDoublet< T1, T2 >::inputTag2_ [private] |
Definition at line 47 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
std::string HLTDoublet< T1, T2 >::label_ [private] |
Definition at line 61 of file HLTDoublet.h.
double HLTDoublet< T1, T2 >::max_DelR_ [private] |
Definition at line 53 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::max_Deta_ [private] |
Definition at line 51 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::max_Dphi_ [private] |
Definition at line 50 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::max_Minv_ [private] |
Definition at line 52 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::max_Pt_ [private] |
Definition at line 54 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::min_DelR_ [private] |
Definition at line 53 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::min_Deta_ [private] |
Definition at line 51 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::min_Dphi_ [private] |
Definition at line 50 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::min_Minv_ [private] |
Definition at line 52 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
int HLTDoublet< T1, T2 >::min_N_ [private] |
Definition at line 55 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
double HLTDoublet< T1, T2 >::min_Pt_ [private] |
Definition at line 54 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
edm::InputTag HLTDoublet< T1, T2 >::originTag1_ [private] |
Definition at line 44 of file HLTDoublet.h.
edm::InputTag HLTDoublet< T1, T2 >::originTag2_ [private] |
Definition at line 45 of file HLTDoublet.h.
bool HLTDoublet< T1, T2 >::same_ [private] |
Definition at line 58 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
int HLTDoublet< T1, T2 >::triggerType1_ [private] |
Definition at line 48 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().
int HLTDoublet< T1, T2 >::triggerType2_ [private] |
Definition at line 49 of file HLTDoublet.h.
Referenced by HLTDoublet< T1, T2 >::HLTDoublet().