CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/HLTrigger/HLTfilters/src/HLTDoubletDZ.cc

Go to the documentation of this file.
00001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00002 #include "HLTrigger/HLTfilters/interface/HLTDoubletDZ.h"
00003 
00004 #include "DataFormats/Common/interface/Handle.h"
00005 
00006 #include "DataFormats/Common/interface/Ref.h"
00007 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00008 
00009 #include "DataFormats/Candidate/interface/Particle.h"
00010 
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 
00013 #include "DataFormats/Math/interface/deltaR.h"
00014 
00015 #include<cmath>
00016 
00017 //
00018 // constructors and destructor
00019 //
00020 template<typename T1, typename T2>
00021 HLTDoubletDZ<T1,T2>::HLTDoubletDZ(const edm::ParameterSet& iConfig) : HLTFilter(iConfig),
00022   originTag1_(iConfig.template getParameter<edm::InputTag>("originTag1")),
00023   originTag2_(iConfig.template getParameter<edm::InputTag>("originTag2")),
00024   inputTag1_(iConfig.template getParameter<edm::InputTag>("inputTag1")),
00025   inputTag2_(iConfig.template getParameter<edm::InputTag>("inputTag2")),
00026   triggerType1_(iConfig.template getParameter<int>("triggerType1")),
00027   triggerType2_(iConfig.template getParameter<int>("triggerType2")),
00028   minDR_ (iConfig.template getParameter<double>("MinDR")),
00029   maxDZ_ (iConfig.template getParameter<double>("MaxDZ")),
00030   min_N_    (iConfig.template getParameter<int>("MinN")),
00031   checkSC_  (iConfig.template getParameter<bool>("checkSC")),
00032   label_    (iConfig.getParameter<std::string>("@module_label")),
00033   coll1_(),
00034   coll2_()
00035 {
00036    // same collections to be compared?
00037    same_ = (inputTag1_.encode()==inputTag2_.encode());
00038 }
00039 
00040 template<typename T1, typename T2>
00041 HLTDoubletDZ<T1,T2>::~HLTDoubletDZ()
00042 {
00043 }
00044 
00045 template<typename T1, typename T2>
00046 void
00047 HLTDoubletDZ<T1,T2>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00048   edm::ParameterSetDescription desc;
00049   makeHLTFilterDescription(desc);
00050   desc.add<edm::InputTag>("originTag1",edm::InputTag("hltOriginal1"));
00051   desc.add<edm::InputTag>("originTag2",edm::InputTag("hltOriginal2"));
00052   desc.add<edm::InputTag>("inputTag1",edm::InputTag("hltFiltered1"));
00053   desc.add<edm::InputTag>("inputTag2",edm::InputTag("hltFiltered2"));
00054   desc.add<int>("triggerType1",0);
00055   desc.add<int>("triggerType2",0);
00056   desc.add<double>("MinDR",-1.0);
00057   desc.add<double>("MaxDZ",0.2);
00058   desc.add<bool>("checkSC",false);
00059   desc.add<int>("MinN",1);
00060   descriptions.add(std::string("hlt")+std::string(typeid(HLTDoubletDZ<T1,T2>).name()),desc);
00061 }
00062 
00063 // ------------ method called to produce the data  ------------
00064 template<typename T1, typename T2>
00065 bool
00066 HLTDoubletDZ<T1,T2>::hltFilter(edm::Event& iEvent, const edm::EventSetup& iSetup, trigger::TriggerFilterObjectWithRefs & filterproduct)
00067 {
00068    using namespace std;
00069    using namespace edm;
00070    using namespace reco;
00071    using namespace trigger;
00072 
00073    // All HLT filters must create and fill an HLT filter object,
00074    // recording any reconstructed physics objects satisfying (or not)
00075    // this HLT filter, and place it in the Event.
00076 
00077    bool accept(false);
00078 
00079    LogVerbatim("HLTDoubletDZ") << " XXX " << label_ << " 0 " << std::endl;
00080 
00081    // get hold of pre-filtered object collections
00082    Handle<TriggerFilterObjectWithRefs> coll1,coll2;
00083    if (iEvent.getByLabel (inputTag1_,coll1) && iEvent.getByLabel (inputTag2_,coll2)) {
00084      coll1_.clear();
00085      coll1->getObjects(triggerType1_,coll1_);
00086      const size_type n1(coll1_.size());
00087      coll2_.clear();
00088      coll2->getObjects(triggerType2_,coll2_);
00089      const size_type n2(coll2_.size());
00090 
00091      if (saveTags()) {
00092        InputTag tagOld;
00093        filterproduct.addCollectionTag(originTag1_);
00094        LogVerbatim("HLTDoubletDZ") << " XXX " << label_ << " 1a " << originTag1_.encode() << std::endl;
00095        tagOld=InputTag();
00096        for (size_type i1=0; i1!=n1; ++i1) {
00097          const ProductID pid(coll1_[i1].id());
00098          const string&    label(iEvent.getProvenance(pid).moduleLabel());
00099          const string& instance(iEvent.getProvenance(pid).productInstanceName());
00100          const string&  process(iEvent.getProvenance(pid).processName());
00101          InputTag tagNew(InputTag(label,instance,process));
00102          if (tagOld.encode()!=tagNew.encode()) {
00103            filterproduct.addCollectionTag(tagNew);
00104            tagOld=tagNew;
00105            LogVerbatim("HLTDoubletDZ") << " XXX " << label_ << " 1b " << tagNew.encode() << std::endl;
00106          }
00107        }
00108        filterproduct.addCollectionTag(originTag2_);
00109        LogVerbatim("HLTDoubletDZ") << " XXX " << label_ << " 2a " << originTag2_.encode() << std::endl;
00110        tagOld=InputTag();
00111        for (size_type i2=0; i2!=n2; ++i2) {
00112          const ProductID pid(coll2_[i2].id());
00113          const string&    label(iEvent.getProvenance(pid).moduleLabel());
00114          const string& instance(iEvent.getProvenance(pid).productInstanceName());
00115          const string&  process(iEvent.getProvenance(pid).processName());
00116          InputTag tagNew(InputTag(label,instance,process));
00117          if (tagOld.encode()!=tagNew.encode()) {
00118            filterproduct.addCollectionTag(tagNew);
00119            tagOld=tagNew;
00120            LogVerbatim("HLTDoubletDZ") << " XXX " << label_ << " 2b " << tagNew.encode() << std::endl;
00121          }
00122        }
00123      }
00124 
00125      int n(0);
00126      T1Ref r1;
00127      T2Ref r2;
00128      Particle::LorentzVector p1,p2,p;
00129      for (unsigned int i1=0; i1!=n1; i1++) {
00130        r1=coll1_[i1];
00131        const reco::Candidate& candidate1(*r1);
00132        unsigned int I(0);
00133        if (same_) {I=i1+1;}
00134        for (unsigned int i2=I; i2!=n2; i2++) {
00135          r2=coll2_[i2];
00136          if (checkSC_) {
00137            if (r1->superCluster().isNonnull() && r2->superCluster().isNonnull()) {
00138              if (r1->superCluster() == r2->superCluster()) continue;
00139            }
00140          }
00141          const reco::Candidate& candidate2(*r2);
00142          if ( reco::deltaR(candidate1, candidate2) < minDR_ ) continue;
00143          if ( std::abs(candidate1.vz()-candidate2.vz()) > maxDZ_ ) continue;
00144          n++;
00145          filterproduct.addObject(triggerType1_,r1);
00146          filterproduct.addObject(triggerType2_,r2);
00147        }
00148      }
00149      // filter decision
00150      accept = accept || (n>=min_N_);
00151    }
00152 
00153    return accept;
00154 }