CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/JetMET/src/HLTMhtHtFilter.cc

Go to the documentation of this file.
00001 
00008 #include "HLTrigger/JetMET/interface/HLTMhtHtFilter.h"
00009 
00010 #include "DataFormats/Common/interface/Handle.h"
00011 
00012 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00013 
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 
00016 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00017 
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019 #include "FWCore/Framework/interface/EventSetup.h"
00020 
00021 #include "DataFormats/Math/interface/deltaPhi.h"
00022 
00023 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00024 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00025 #include "FWCore/Utilities/interface/InputTag.h"
00026 #include <vector>
00027 
00028 
00029 //
00030 // constructors and destructor
00031 //
00032 HLTMhtHtFilter::HLTMhtHtFilter(const edm::ParameterSet& iConfig)
00033 {
00034   inputJetTag_ = iConfig.getParameter< edm::InputTag > ("inputJetTag");
00035   saveTag_     = iConfig.getUntrackedParameter<bool>("saveTag",false);
00036   minMht_= iConfig.getParameter<double> ("minMht");
00037   minPtJet_= iConfig.getParameter<std::vector<double> > ("minPtJet");
00038   minNJet_= iConfig.getParameter<int> ("minNJet");
00039   mode_= iConfig.getParameter<int>("mode");
00040   //----mode=1 for MHT only
00041   //----mode=2 for Meff
00042   //----mode=3 for PT12
00043   //----mode=4 for HT only
00044   //----mode=5 for HT and AlphaT cross trigger (ALWAYS uses jet ET, not pT)
00045   etaJet_= iConfig.getParameter<std::vector<double> > ("etaJet");
00046   usePt_= iConfig.getParameter<bool>("usePt");
00047   minPT12_= iConfig.getParameter<double> ("minPT12");
00048   minMeff_= iConfig.getParameter<double> ("minMeff");
00049   minHt_= iConfig.getParameter<double> ("minHt");
00050   minAlphaT_= iConfig.getParameter<double> ("minAlphaT");
00051 
00052   // sanity checks
00053   if (       (minPtJet_.size()    !=  etaJet_.size())
00054        || (  (minPtJet_.size()<1) || (etaJet_.size()<1) )
00055        || ( ((minPtJet_.size()<2) || (etaJet_.size()<2))
00056             && ( (mode_==1) || (mode_==2) || (mode_ == 5))) ) {
00057     edm::LogError("HLTMhtHtFilter") << "inconsistent module configuration!";
00058   }
00059 
00060   //register your products
00061   produces<trigger::TriggerFilterObjectWithRefs>();
00062 }
00063 
00064 HLTMhtHtFilter::~HLTMhtHtFilter(){}
00065 
00066 void HLTMhtHtFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00067   edm::ParameterSetDescription desc;
00068   desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
00069   desc.addUntracked<bool>("saveTag",false);
00070   desc.add<double>("minMht",0.0);
00071   {
00072     std::vector<double> temp1;
00073     temp1.reserve(2);
00074     temp1.push_back(20.0);
00075     temp1.push_back(20.0);
00076     desc.add<std::vector<double> >("minPtJet",temp1);
00077   }
00078   desc.add<int>("minNJet",0);
00079   desc.add<int>("mode",2);
00080   {
00081     std::vector<double> temp1;
00082     temp1.reserve(2);
00083     temp1.push_back(9999.0);
00084     temp1.push_back(9999.0);
00085     desc.add<std::vector<double> >("etaJet",temp1);
00086   }
00087   desc.add<bool>("usePt",true);
00088   desc.add<double>("minPT12",0.0);
00089   desc.add<double>("minMeff",180.0);
00090   desc.add<double>("minHt",0.0);
00091   desc.add<double>("minAlphaT",0.0);
00092   descriptions.add("hltMhtHtFilter",desc);
00093 }
00094 
00095 
00096 
00097 // ------------ method called to produce the data  ------------
00098 bool
00099   HLTMhtHtFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00100 {
00101   using namespace std;
00102   using namespace edm;
00103   using namespace reco;
00104   using namespace trigger;
00105   // The filter object
00106   auto_ptr<trigger::TriggerFilterObjectWithRefs> filterobject (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00107   if (saveTag_) filterobject->addCollectionTag(inputJetTag_);
00108 
00109   CaloJetRef ref;
00110   // Get the Candidates
00111 
00112   Handle<CaloJetCollection> recocalojets;
00113   iEvent.getByLabel(inputJetTag_,recocalojets);
00114 
00115   // look at all candidates,  check cuts and add to filter object
00116   int n(0), nj(0), flag(0);
00117   double ht=0.;
00118   double mhtx=0., mhty=0.;
00119   double jetVar;
00120   double dht = 0.;
00121   double aT = 0.;
00122   if(recocalojets->size() > 0){
00123     // events with at least one jet
00124     for (CaloJetCollection::const_iterator recocalojet = recocalojets->begin();
00125     recocalojet != recocalojets->end(); recocalojet++) {
00126       if (flag == 1){break;}
00127       jetVar = recocalojet->pt();
00128       if (!usePt_ || mode_==3 ) jetVar = recocalojet->et();
00129 
00130       if (mode_==1 || mode_==2 || mode_ == 5) {//---get MHT
00131         if (jetVar > minPtJet_.at(1) && fabs(recocalojet->eta()) < etaJet_.at(1)) {
00132           mhtx -= jetVar*cos(recocalojet->phi());
00133           mhty -= jetVar*sin(recocalojet->phi());
00134         }
00135       }
00136       if (mode_==2 || mode_==4 || mode_==5) {//---get HT
00137         if (jetVar > minPtJet_.at(0) && fabs(recocalojet->eta()) < etaJet_.at(0)) {
00138           ht += jetVar;
00139           nj++;
00140         }
00141       }
00142       if (mode_==3) {//---get PT12
00143         if (jetVar > minPtJet_.at(0) && fabs(recocalojet->eta()) < etaJet_.at(0)) {
00144           nj++;
00145           mhtx -= jetVar*cos(recocalojet->phi());
00146           mhty -= jetVar*sin(recocalojet->phi());
00147           if (nj==2) break;
00148         }
00149       }
00150       if(mode_ == 5){
00151         double mHT = sqrt( (mhtx*mhtx) + (mhty*mhty) );
00152         dht += ( nj < 2 ? jetVar : -1.* jetVar ); //@@ only use for njets < 4
00153         if ( nj == 2 || nj == 3 ) {
00154           aT = ( ht - fabs(dht) ) / ( 2. * sqrt( ( ht*ht ) - ( mHT*mHT  ) ) );
00155         } else if ( nj > 3 ) {
00156           aT = ht / ( 2.*sqrt( ( ht*ht ) - ( mHT*mHT  ) ) );
00157         }
00158         if(ht > minHt_ && aT > minAlphaT_){
00159   // put filter object into the Event
00160           flag = 1;
00161         }
00162       }
00163     }
00164 
00165   if( mode_==1 && sqrt(mhtx*mhtx + mhty*mhty) > minMht_) flag=1;
00166   if( mode_==2 && sqrt(mhtx*mhtx + mhty*mhty)+ht > minMeff_) flag=1;
00167   if( mode_==3 && sqrt(mhtx*mhtx + mhty*mhty) > minPT12_ && nj>1) flag=1;
00168   if( mode_==4 && ht > minHt_) flag=1;
00169 
00170   if (flag==1) {
00171     for (reco::CaloJetCollection::const_iterator recocalojet = recocalojets->begin(); recocalojet!=recocalojets->end(); recocalojet++) {
00172       jetVar = recocalojet->pt();
00173       if (!usePt_ || mode_==3) jetVar = recocalojet->et();
00174 
00175       if (jetVar > minPtJet_.at(0)) {
00176         ref = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
00177         filterobject->addObject(TriggerJet,ref);
00178         n++;
00179       }
00180     }
00181   }
00182 } // events with at least one jet
00183 
00184 
00185 
00186   // filter decision
00187 bool accept(n>0);
00188 
00189   // put filter object into the Event
00190 iEvent.put(filterobject);
00191 
00192 return accept;
00193 }