CMS 3D CMS Logo

HLTJetCollectionsForElePlusJets.cc
Go to the documentation of this file.
1 #include <string>
2 #include <vector>
3 
4 #include "TVector3.h"
5 
14 
15 template <typename T>
17  : hltElectronTag(iConfig.getParameter<edm::InputTag>("HltElectronTag")),
18  sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
19  //minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
20  //maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
21  //minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
22  minDeltaR_(iConfig.getParameter<double>("minDeltaR"))
23 //Only for VBF
24 //minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
25 //minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
26 {
27  typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
28  TCollectionVector;
29  m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
30  m_theJetToken = consumes<std::vector<T>>(sourceJetTag);
31  produces<TCollectionVector>();
32 }
33 
34 template <typename T>
36  // do anything here that needs to be done at desctruction time
37  // (e.g. close files, deallocate resources etc.)
38 }
39 
40 template <typename T>
43  desc.add<edm::InputTag>("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
44  desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
45  // desc.add<double> ("MinJetPt", 30.);
46  // desc.add<double> ("MaxAbsJetEta", 2.6);
47  // desc.add<unsigned int> ("MinNJets", 1);
48  desc.add<double>("minDeltaR", 0.5);
49  //Only for VBF
50  // desc.add<double> ("MinSoftJetPt", 25.);
51  //desc.add<double> ("MinDeltaEta", -1.);
53 }
54 
55 //
56 // member functions
57 //
58 
59 // ------------ method called to produce the data ------------
60 template <typename T>
62  using namespace edm;
63  using namespace std;
64 
65  typedef vector<T> TCollection;
66  typedef Ref<TCollection> TRef;
67 
68  typedef edm::RefVector<TCollection> TRefVector;
69 
70  typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
71  TCollectionVector;
72 
74  iEvent.getByToken(m_theElectronToken, PrevFilterOutput);
75 
76  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
77  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
78  std::vector<edm::Ref<reco::RecoEcalCandidateCollection>> clusCands;
79  PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
80 
81  std::vector<edm::Ref<reco::ElectronCollection>> eleCands;
82  PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
83 
84  trigger::VRphoton photonCands;
85  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
86 
87  //prepare the collection of 3-D vector for electron momenta
88  std::vector<TVector3> ElePs;
89 
90  if (!clusCands.empty()) { //try trigger cluster
91  for (auto& clusCand : clusCands) {
92  TVector3 positionVector(clusCand->superCluster()->position().x(),
93  clusCand->superCluster()->position().y(),
94  clusCand->superCluster()->position().z());
95  ElePs.push_back(positionVector);
96  }
97  } else if (!eleCands.empty()) { // try trigger electrons
98  for (auto& eleCand : eleCands) {
99  TVector3 positionVector(eleCand->superCluster()->position().x(),
100  eleCand->superCluster()->position().y(),
101  eleCand->superCluster()->position().z());
102  ElePs.push_back(positionVector);
103  }
104  } else if (!photonCands.empty()) { // try trigger photons
105  for (auto& photonCand : photonCands) {
106  TVector3 positionVector(photonCand->superCluster()->position().x(),
107  photonCand->superCluster()->position().y(),
108  photonCand->superCluster()->position().z());
109  ElePs.push_back(positionVector);
110  }
111  }
112 
113  edm::Handle<TCollection> theJetCollectionHandle;
114  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
115 
116  const TCollection& theJetCollection = *theJetCollectionHandle;
117 
118  //std::unique_ptr< TCollection > theFilteredJetCollection(new TCollection);
119 
120  std::unique_ptr<TCollectionVector> allSelections(new TCollectionVector());
121 
122  //bool foundSolution(false);
123 
124  for (auto& EleP : ElePs) {
125  // bool VBFJetPair = false;
126  //std::vector<int> store_jet;
127  TRefVector refVector;
128 
129  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
130  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
131  double DR = EleP.DeltaR(JetP);
132 
133  if (DR > minDeltaR_)
134  refVector.push_back(TRef(theJetCollectionHandle, j));
135  }
136  allSelections->push_back(refVector);
137  }
138 
139  //iEvent.put(std::move(theFilteredJetCollection));
140  iEvent.put(std::move(allSelections));
141 
142  return;
143 }
defaultModuleLabel.h
HLTJetCollectionsForElePlusJets::HLTJetCollectionsForElePlusJets
HLTJetCollectionsForElePlusJets(const edm::ParameterSet &)
Definition: HLTJetCollectionsForElePlusJets.cc:16
Handle.h
trigger::TriggerElectron
Definition: TriggerTypeDefs.h:78
trigger::VRphoton
std::vector< reco::RecoEcalCandidateRef > VRphoton
Definition: TriggerRefsCollections.h:67
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
edm::RefVector
Definition: EDProductfwd.h:27
HLTJetCollectionsForElePlusJets::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: HLTJetCollectionsForElePlusJets.cc:61
edm::Handle< trigger::TriggerFilterObjectWithRefs >
RecoCandidate.h
trigger::TriggerRefsCollections::getObjects
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Definition: TriggerRefsCollections.h:590
edm::Ref
Definition: AssociativeIterator.h:58
HLTJetCollectionsForElePlusJets::~HLTJetCollectionsForElePlusJets
~HLTJetCollectionsForElePlusJets() override
Definition: HLTJetCollectionsForElePlusJets.cc:35
MakerMacros.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
HLTJetCollectionsForElePlusJets::sourceJetTag
edm::InputTag sourceJetTag
Definition: HLTJetCollectionsForElePlusJets.h:52
edm::ParameterSet
Definition: ParameterSet.h:47
defaultModuleLabel
std::string defaultModuleLabel()
Definition: defaultModuleLabel.h:16
HLTJetCollectionsForElePlusJets.h
iEvent
int iEvent
Definition: GenABIO.cc:224
HLTJetCollectionsForElePlusJets::hltElectronTag
edm::InputTag hltElectronTag
Definition: HLTJetCollectionsForElePlusJets.h:51
edm::EventSetup
Definition: EventSetup.h:58
trigger::TriggerCluster
Definition: TriggerTypeDefs.h:88
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
HLTJetCollectionsForElePlusJets
Definition: HLTJetCollectionsForElePlusJets.h:40
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
Electron.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
RecoEcalCandidate.h
T
long double T
Definition: Basic3DVectorLD.h:48
HLTJetCollectionsForElePlusJets::m_theElectronToken
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theElectronToken
Definition: HLTJetCollectionsForElePlusJets.h:49
SuperCluster.h
trigger::TriggerPhoton
HLT.
Definition: TriggerTypeDefs.h:77
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
HLTJetCollectionsForElePlusJets::m_theJetToken
edm::EDGetTokenT< std::vector< T > > m_theJetToken
Definition: HLTJetCollectionsForElePlusJets.h:50
edm::refhelper::FindUsingAdvance
Definition: RefTraits.h:14
HLTJetCollectionsForElePlusJets::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTJetCollectionsForElePlusJets.cc:41