CMS 3D CMS Logo

HLTGlobalSums.cc

Go to the documentation of this file.
00001 
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "HLTrigger/HLTfilters/interface/HLTGlobalSums.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 "FWCore/MessageLogger/interface/MessageLogger.h"
00021 
00022 #include<cmath>
00023 
00024 //
00025 // constructors and destructor
00026 //
00027 template<typename T, int Tid>
00028 HLTGlobalSums<T,Tid>::HLTGlobalSums(const edm::ParameterSet& iConfig) :
00029   inputTag_   (iConfig.template getParameter<edm::InputTag>("inputTag")),
00030   saveTag_    (iConfig.template getUntrackedParameter<bool> ("saveTag",false)),
00031   observable_ (iConfig.template getParameter<std::string>("observable")),
00032   min_        (iConfig.template getParameter<double>("Min")),
00033   max_        (iConfig.template getParameter<double>("Max")),
00034   min_N_      (iConfig.template getParameter<int>("MinN")),
00035   tid_()
00036 {
00037    LogDebug("") << "InputTags and cuts : " 
00038                 << inputTag_.encode() << " " << observable_
00039                 << " Range [" << min_ << " " << max_ << "]"
00040                 << " MinN =" << min_N_
00041      ;
00042 
00043    if (observable_=="sumEt") {
00044      tid_=trigger::TriggerHT;
00045    } else if (observable_=="mEtSig") {
00046      tid_=trigger::TriggerMETSig;
00047    } else if (observable_=="e_longitudinal") {
00048      tid_=trigger::TriggerELongit;
00049    } else {
00050      tid_=Tid;
00051    }
00052 
00053    //register your products
00054    produces<trigger::TriggerFilterObjectWithRefs>();
00055 }
00056 
00057 template<typename T, int Tid>
00058 HLTGlobalSums<T,Tid>::~HLTGlobalSums()
00059 {
00060 }
00061 
00062 //
00063 // member functions
00064 //
00065 
00066 // ------------ method called to produce the data  ------------
00067 template<typename T, int Tid> 
00068 bool
00069 HLTGlobalSums<T,Tid>::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00070 {
00071    using namespace std;
00072    using namespace edm;
00073    using namespace reco;
00074    using namespace trigger;
00075 
00076    typedef vector<T> TCollection;
00077    typedef Ref<TCollection> TRef;
00078 
00079    // All HLT filters must create and fill an HLT filter object,
00080    // recording any reconstructed physics objects satisfying (or not)
00081    // this HLT filter, and place it in the Event.
00082 
00083    // The filter object
00084    auto_ptr<TriggerFilterObjectWithRefs>
00085      filterobject (new TriggerFilterObjectWithRefs(path(),module()));
00086    if (saveTag_) filterobject->addCollectionTag(inputTag_);
00087    // Ref to Candidate object to be recorded in filter object
00088    TRef ref;
00089 
00090 
00091    // get hold of MET product from Event
00092    Handle<TCollection>   objects;
00093    iEvent.getByLabel(inputTag_,objects);
00094    if (!objects.isValid()) {
00095      LogDebug("") << inputTag_ << " collection not found!";
00096      iEvent.put(filterobject);
00097      return false;
00098    }
00099 
00100    LogDebug("") << "Size of MET collection: " << objects->size();
00101    if (objects->size()==0) {
00102      LogDebug("") << "MET collection does not contain a MET object!";
00103    } else if (objects->size()>1) {
00104      LogDebug("") << "MET collection contains more than one MET object!";
00105    }
00106 
00107    int n(0);
00108    double value(0.0);
00109    typename TCollection::const_iterator ibegin(objects->begin());
00110    typename TCollection::const_iterator iend(objects->end());
00111    typename TCollection::const_iterator iter;
00112    for (iter=ibegin; iter!=iend; iter++) {
00113 
00114      // get hold of value of observable to cut on
00115      if (tid_==TriggerHT) {
00116        value=iter->sumEt();
00117      } else if (tid_==TriggerMETSig) {
00118        value=iter->mEtSig();
00119      } else if (tid_==TriggerELongit) {
00120        value=iter->e_longitudinal();
00121      } else {
00122        value=0.0;
00123      }
00124 
00125      value=abs(value);
00126 
00127      if ( ( (min_<0.0) || (min_<=value) ) &&
00128           ( (max_<0.0) || (value<=max_) ) ) {
00129        n++;
00130        ref=TRef(objects,distance(ibegin,iter));
00131        filterobject->addObject(tid_,ref);
00132      }
00133 
00134    }
00135 
00136    // filter decision
00137    const bool accept(n>=min_N_);
00138 
00139    // put filter object into the Event
00140    iEvent.put(filterobject);
00141 
00142    return accept;
00143 }

Generated on Tue Jun 9 17:37:55 2009 for CMSSW by  doxygen 1.5.4