CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/HLTfilters/src/HLTDoublet.cc

Go to the documentation of this file.
00001 
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "HLTrigger/HLTfilters/interface/HLTDoublet.h"
00014 
00015 #include "DataFormats/Common/interface/Handle.h"
00016 
00017 #include "DataFormats/Common/interface/Ref.h"
00018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00019 
00020 #include "DataFormats/Candidate/interface/Particle.h"
00021 
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 
00024 #include<cmath>
00025 
00026 //
00027 // constructors and destructor
00028 //
00029 template<typename T1, int Tid1, typename T2, int Tid2>
00030 HLTDoublet<T1,Tid1,T2,Tid2>::HLTDoublet(const edm::ParameterSet& iConfig) :
00031   inputTag1_(iConfig.template getParameter<edm::InputTag>("inputTag1")),
00032   inputTag2_(iConfig.template getParameter<edm::InputTag>("inputTag2")),
00033   saveTags_ (iConfig.template getUntrackedParameter<bool>("saveTags",false)),
00034   min_Dphi_ (iConfig.template getParameter<double>("MinDphi")),
00035   max_Dphi_ (iConfig.template getParameter<double>("MaxDphi")),
00036   min_Deta_ (iConfig.template getParameter<double>("MinDeta")),
00037   max_Deta_ (iConfig.template getParameter<double>("MaxDeta")),
00038   min_Minv_ (iConfig.template getParameter<double>("MinMinv")),
00039   max_Minv_ (iConfig.template getParameter<double>("MaxMinv")),
00040   min_DelR_ (iConfig.template getParameter<double>("MinDelR")),
00041   max_DelR_ (iConfig.template getParameter<double>("MaxDelR")),
00042   min_N_    (iConfig.template getParameter<int>("MinN")),
00043   coll1_(),
00044   coll2_()
00045 {
00046    // same collections to be compared?
00047    same_ = (inputTag1_.encode()==inputTag2_.encode());
00048 
00049    cutdphi_ = (min_Dphi_ <= max_Dphi_); // cut active?
00050    cutdeta_ = (min_Deta_ <= max_Deta_); // cut active?
00051    cutminv_ = (min_Minv_ <= max_Minv_); // cut active?
00052    cutdelr_ = (min_DelR_ <= max_DelR_); // cut active?
00053 
00054    LogDebug("") << "InputTags and cuts : " 
00055                 << inputTag1_.encode() << " " << inputTag2_.encode()
00056                 << " Dphi [" << min_Dphi_ << " " << max_Dphi_ << "]"
00057                 << " Deta [" << min_Deta_ << " " << max_Deta_ << "]"
00058                 << " Minv [" << min_Minv_ << " " << max_Minv_ << "]"
00059                 << " DelR [" << min_DelR_ << " " << max_DelR_ << "]"
00060                 << " MinN =" << min_N_
00061                 << " same/dphi/deta/minv/delr "
00062                 << same_ << cutdphi_ << cutdeta_ << cutminv_ << cutdelr_;
00063 
00064    //register your products
00065    produces<trigger::TriggerFilterObjectWithRefs>();
00066 }
00067 
00068 template<typename T1, int Tid1, typename T2, int Tid2>
00069 HLTDoublet<T1,Tid1,T2,Tid2>::~HLTDoublet()
00070 {
00071 }
00072 
00073 //
00074 // member functions
00075 //
00076 
00077 // ------------ method called to produce the data  ------------
00078 template<typename T1, int Tid1, typename T2, int Tid2>
00079 bool
00080 HLTDoublet<T1,Tid1,T2,Tid2>::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00081 {
00082    using namespace std;
00083    using namespace edm;
00084    using namespace reco;
00085    using namespace trigger;
00086 
00087    // All HLT filters must create and fill an HLT filter object,
00088    // recording any reconstructed physics objects satisfying (or not)
00089    // this HLT filter, and place it in the Event.
00090 
00091    // The filter object
00092    auto_ptr<TriggerFilterObjectWithRefs>
00093      filterobject (new TriggerFilterObjectWithRefs(path(),module()));
00094    // Don't saveTag the TFOWRs, but rather the collections pointed to!
00095    //   if (saveTags_) {
00096    //     filterobject->addCollectionTag(inputTag1_);
00097    //     filterobject->addCollectionTag(inputTag2_);
00098    //   }
00099    bool accept(false);
00100 
00101    // get hold of pre-filtered object collections
00102    Handle<TriggerFilterObjectWithRefs> coll1,coll2;
00103    if (iEvent.getByLabel (inputTag1_,coll1) && iEvent.getByLabel (inputTag2_,coll2)) {
00104      coll1_.clear();
00105      coll1->getObjects(Tid1,coll1_);
00106      const size_type n1(coll1_.size());
00107      coll2_.clear();
00108      coll2->getObjects(Tid2,coll2_);
00109      const size_type n2(coll2_.size());
00110 
00111      if (saveTags_) {
00112        InputTag tagOld;
00113        tagOld=InputTag();
00114        for (size_type i1=0; i1!=n1; ++i1) {
00115          const ProductID pid(coll1_[i1].id());
00116          const string&    label(iEvent.getProvenance(pid).moduleLabel());
00117          const string& instance(iEvent.getProvenance(pid).productInstanceName());
00118          const string&  process(iEvent.getProvenance(pid).processName());
00119          InputTag tagNew(InputTag(label,instance,process));
00120          if (tagOld.encode()!=tagNew.encode()) {
00121            filterobject->addCollectionTag(tagNew);
00122            tagOld=tagNew;
00123          }
00124        }
00125        tagOld=InputTag();
00126        for (size_type i2=0; i2!=n2; ++i2) {
00127          const ProductID pid(coll2_[i2].id());
00128          const string&    label(iEvent.getProvenance(pid).moduleLabel());
00129          const string& instance(iEvent.getProvenance(pid).productInstanceName());
00130          const string&  process(iEvent.getProvenance(pid).processName());
00131          InputTag tagNew(InputTag(label,instance,process));
00132          if (tagOld.encode()!=tagNew.encode()) {
00133            filterobject->addCollectionTag(tagNew);
00134            tagOld=tagNew;
00135          }
00136        }
00137      }
00138 
00139      int n(0);
00140      T1Ref r1;
00141      T2Ref r2;
00142      Particle::LorentzVector p1,p2,p;
00143      for (unsigned int i1=0; i1!=n1; i1++) {
00144        r1=coll1_[i1];
00145        p1=r1->p4();
00146        unsigned int I(0);
00147        if (same_) {I=i1+1;}
00148        for (unsigned int i2=I; i2!=n2; i2++) {
00149          r2=coll2_[i2];
00150          p2=r2->p4();
00151 
00152          double Dphi(std::abs(p1.phi()-p2.phi()));
00153          if (Dphi>M_PI) Dphi=2.0*M_PI-Dphi;
00154          
00155          double Deta(std::abs(p1.eta()-p2.eta()));
00156          
00157          p=p1+p2;
00158          double Minv(std::abs(p.mass()));
00159 
00160          double DelR(sqrt(Dphi*Dphi+Deta*Deta));
00161          
00162          if ( ( (!cutdphi_) || ((min_Dphi_<=Dphi) && (Dphi<=max_Dphi_)) ) &&
00163               ( (!cutdeta_) || ((min_Deta_<=Deta) && (Deta<=max_Deta_)) ) &&
00164               ( (!cutminv_) || ((min_Minv_<=Minv) && (Minv<=max_Minv_)) ) &&
00165               ( (!cutdelr_) || ((min_DelR_<=DelR) && (DelR<=max_DelR_)) ) ) {
00166            n++;
00167            filterobject->addObject(Tid1,r1);
00168            filterobject->addObject(Tid2,r2);
00169          }
00170          
00171        }
00172      }
00173      // filter decision
00174      accept = accept || (n>=min_N_);
00175    }
00176 
00177    iEvent.put(filterobject);
00178    return accept;
00179 }