File indexing completed on 2024-04-06 12:18:27
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include <cmath>
0011
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "HLTGlobalSums.h"
0014
0015 #include "DataFormats/Common/interface/Handle.h"
0016
0017 #include "DataFormats/Common/interface/Ref.h"
0018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0019
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021
0022 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0023
0024
0025
0026
0027 template <typename T>
0028 HLTGlobalSums<T>::HLTGlobalSums(const edm::ParameterSet& iConfig)
0029 : HLTFilter(iConfig),
0030 inputTag_(iConfig.template getParameter<edm::InputTag>("inputTag")),
0031 inputToken_(consumes<std::vector<T>>(inputTag_)),
0032 triggerType_(iConfig.template getParameter<int>("triggerType")),
0033 observable_(iConfig.template getParameter<std::string>("observable")),
0034 min_(iConfig.template getParameter<double>("Min")),
0035 max_(iConfig.template getParameter<double>("Max")),
0036 min_N_(iConfig.template getParameter<int>("MinN")),
0037 tid_(triggerType_) {
0038 LogDebug("") << "InputTags and cuts : " << inputTag_.encode() << " " << triggerType_ << " " << observable_
0039 << " Range [" << min_ << " " << max_ << "]"
0040 << " MinN =" << min_N_;
0041
0042 if (observable_ == "sumEt") {
0043 tid_ = triggerType_;
0044 } else if (observable_ == "mEtSig") {
0045 if (triggerType_ == trigger::TriggerTET) {
0046 tid_ = trigger::TriggerMETSig;
0047 } else if (triggerType_ == trigger::TriggerTHT) {
0048 tid_ = trigger::TriggerMHTSig;
0049 } else {
0050 tid_ = triggerType_;
0051 }
0052 } else if (observable_ == "e_longitudinal") {
0053 if (triggerType_ == trigger::TriggerTET) {
0054 tid_ = trigger::TriggerELongit;
0055 } else if (triggerType_ == trigger::TriggerTHT) {
0056 tid_ = trigger::TriggerHLongit;
0057 } else {
0058 tid_ = triggerType_;
0059 }
0060 } else {
0061 tid_ = triggerType_;
0062 }
0063 }
0064
0065 template <typename T>
0066 HLTGlobalSums<T>::~HLTGlobalSums() = default;
0067
0068 template <typename T>
0069 void HLTGlobalSums<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0070 edm::ParameterSetDescription desc;
0071 makeHLTFilterDescription(desc);
0072 desc.add<edm::InputTag>("inputTag", edm::InputTag("hltCollection"));
0073 desc.add<int>("triggerType", 0);
0074 desc.add<std::string>("observable", "");
0075 desc.add<double>("Min", -1e125);
0076 desc.add<double>("Max", +1e125);
0077 desc.add<int>("MinN", 1);
0078 descriptions.add(defaultModuleLabel<HLTGlobalSums<T>>(), desc);
0079 }
0080
0081
0082
0083
0084
0085
0086 template <typename T>
0087 bool HLTGlobalSums<T>::hltFilter(edm::Event& iEvent,
0088 const edm::EventSetup& iSetup,
0089 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0090 using namespace std;
0091 using namespace edm;
0092 using namespace reco;
0093 using namespace trigger;
0094
0095 typedef vector<T> TCollection;
0096 typedef Ref<TCollection> TRef;
0097
0098
0099
0100
0101
0102
0103 if (saveTags())
0104 filterproduct.addCollectionTag(inputTag_);
0105
0106 TRef ref;
0107
0108
0109 Handle<TCollection> objects;
0110 iEvent.getByToken(inputToken_, objects);
0111 if (!objects.isValid()) {
0112 LogDebug("") << inputTag_ << " collection not found!";
0113 return false;
0114 }
0115
0116 LogDebug("") << "Size of MET collection: " << objects->size();
0117 if (objects->empty()) {
0118 LogDebug("") << "MET collection does not contain a MET object!";
0119 } else if (objects->size() > 1) {
0120 LogDebug("") << "MET collection contains more than one MET object!";
0121 }
0122
0123 int n(0);
0124 double value(0.0);
0125 typename TCollection::const_iterator ibegin(objects->begin());
0126 typename TCollection::const_iterator iend(objects->end());
0127 typename TCollection::const_iterator iter;
0128 for (iter = ibegin; iter != iend; iter++) {
0129
0130 if ((tid_ == TriggerTET) || (tid_ == TriggerTHT)) {
0131 value = iter->sumEt();
0132 } else if ((tid_ == TriggerMETSig) || (tid_ == TriggerMHTSig)) {
0133 value = iter->mEtSig();
0134 } else if ((tid_ == TriggerELongit) || (tid_ == TriggerHLongit)) {
0135 value = iter->e_longitudinal();
0136 } else {
0137 value = 0.0;
0138 }
0139
0140 value = std::abs(value);
0141
0142 if (((min_ < 0.0) || (min_ <= value)) && ((max_ < 0.0) || (value <= max_))) {
0143 n++;
0144 ref = TRef(objects, distance(ibegin, iter));
0145 filterproduct.addObject(tid_, ref);
0146 }
0147 }
0148
0149
0150 const bool accept(n >= min_N_);
0151
0152 return accept;
0153 }