File indexing completed on 2024-04-06 12:18:21
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "HLTPMDocaFilter.h"
0010
0011 #include "DataFormats/Common/interface/Handle.h"
0012
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016
0017 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0018 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
0019
0020
0021
0022
0023 HLTPMDocaFilter::HLTPMDocaFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0024 candTag_ = iConfig.getParameter<edm::InputTag>("candTag");
0025 docaDiffPerpCutHigh_ = iConfig.getParameter<double>("docaDiffPerpCutHigh");
0026 docaDiffPerpCutLow_ = iConfig.getParameter<double>("docaDiffPerpCutLow");
0027 nZcandcut_ = iConfig.getParameter<int>("nZcandcut");
0028 candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
0029 }
0030
0031 HLTPMDocaFilter::~HLTPMDocaFilter() = default;
0032
0033 void HLTPMDocaFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0034 edm::ParameterSetDescription desc;
0035 makeHLTFilterDescription(desc);
0036 desc.add<edm::InputTag>("candTag", edm::InputTag("HltZeePMMassFilter"));
0037 desc.add<double>("docaDiffPerpCutHigh", 0.055691);
0038 desc.add<double>("docaDiffPerpCutLow", 0.0);
0039 desc.add<int>("nZcandcut", 1);
0040 descriptions.add("hltPMDocaFilter", desc);
0041 }
0042
0043
0044 bool HLTPMDocaFilter::hltFilter(edm::Event& iEvent,
0045 const edm::EventSetup& iSetup,
0046 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0047 using namespace std;
0048 using namespace edm;
0049 using namespace reco;
0050 using namespace trigger;
0051
0052
0053 edm::Ref<reco::ElectronCollection> ref;
0054
0055 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0056 iEvent.getByToken(candToken_, PrevFilterOutput);
0057
0058 std::vector<edm::Ref<reco::ElectronCollection> > electrons;
0059 PrevFilterOutput->getObjects(TriggerElectron, electrons);
0060
0061 int n = 0;
0062
0063 unsigned int size = electrons.size();
0064 std::vector<double> vx(size);
0065 std::vector<double> vy(size);
0066
0067 for (unsigned int i = 0; i < size; i++) {
0068 ref = electrons[i];
0069 vx[i] = ref->vx();
0070 vy[i] = ref->vy();
0071 }
0072
0073 for (unsigned int jj = 0; jj < size; jj++) {
0074 for (unsigned int ii = jj + 1; ii < size; ii++) {
0075 double docaDiffPerp = sqrt((vx[jj] - vx[ii]) * (vx[jj] - vx[ii]) + (vy[jj] - vy[ii]) * (vy[jj] - vy[ii]));
0076
0077 if ((docaDiffPerp >= docaDiffPerpCutLow_) && (docaDiffPerp <= docaDiffPerpCutHigh_)) {
0078 n++;
0079 ref = electrons[ii];
0080 filterproduct.addObject(TriggerElectron, ref);
0081 ref = electrons[jj];
0082 filterproduct.addObject(TriggerElectron, ref);
0083 }
0084 }
0085 }
0086
0087
0088 bool accept(n >= nZcandcut_);
0089
0090 return accept;
0091 }
0092
0093
0094 #include "FWCore/Framework/interface/MakerMacros.h"
0095 DEFINE_FWK_MODULE(HLTPMDocaFilter);