CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
L1TrackFastJetProducer Class Reference
Inheritance diagram for L1TrackFastJetProducer:
edm::stream::EDProducer<>

Public Types

typedef std::vector< L1TTTrackTypeL1TTTrackCollectionType
 
typedef TTTrack< Ref_Phase2TrackerDigi_L1TTTrackType
 
- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Public Member Functions

 L1TrackFastJetProducer (const edm::ParameterSet &)
 
 ~L1TrackFastJetProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const double coneSize_
 
const double deltaZ0Cut_
 
const bool displaced_
 
const bool doTightChi2_
 
edm::EDGetTokenT< VertexCollectionpvToken_
 
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
 
const double trkBendChi2Max_
 
const float trkChi2dofMax_
 
const float trkChi2dofTightChi2_
 
const float trkEtaMax_
 
const int trkNPSStubMin_
 
const int trkNStubMin_
 
const float trkPtMin_
 
const float trkPtTightChi2_
 
const float trkZMax_
 
edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtTopoToken_
 

Detailed Description

Definition at line 51 of file L1TrackFastJetProducer.cc.

Member Typedef Documentation

◆ L1TTTrackCollectionType

Definition at line 54 of file L1TrackFastJetProducer.cc.

◆ L1TTTrackType

Definition at line 53 of file L1TrackFastJetProducer.cc.

Constructor & Destructor Documentation

◆ L1TrackFastJetProducer()

L1TrackFastJetProducer::L1TrackFastJetProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 87 of file L1TrackFastJetProducer.cc.

References displaced_.

88  : trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
89  trkChi2dofMax_((float)iConfig.getParameter<double>("trk_chi2dofMax")),
90  trkBendChi2Max_(iConfig.getParameter<double>("trk_bendChi2Max")),
91  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
92  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
93  trkNStubMin_((int)iConfig.getParameter<int>("trk_nStubMin")),
94  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
95  deltaZ0Cut_((float)iConfig.getParameter<double>("deltaZ0Cut")),
96  coneSize_((float)iConfig.getParameter<double>("coneSize")),
97  doTightChi2_(iConfig.getParameter<bool>("doTightChi2")),
98  trkPtTightChi2_((float)iConfig.getParameter<double>("trk_ptTightChi2")),
99  trkChi2dofTightChi2_((float)iConfig.getParameter<double>("trk_chi2dofTightChi2")),
100  displaced_(iConfig.getParameter<bool>("displaced")),
102  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
103  pvToken_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("L1PrimaryVertexTag"))),
104  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
105  if (displaced_)
106  produces<TkJetCollection>("L1TrackFastJetsExtended");
107  else
108  produces<TkJetCollection>("L1TrackFastJets");
109 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< VertexCollection > pvToken_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29

◆ ~L1TrackFastJetProducer()

L1TrackFastJetProducer::~L1TrackFastJetProducer ( )
override

Definition at line 112 of file L1TrackFastJetProducer.cc.

112 {}

Member Function Documentation

◆ fillDescriptions()

void L1TrackFastJetProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 220 of file L1TrackFastJetProducer.cc.

References edm::ConfigurationDescriptions::addDefault(), and submitPVResolutionJobs::desc.

220  {
221  //The following says we do not know what parameters are allowed so do no validation
222  // Please change this to state exactly what you do use, even if it is no parameters
224  desc.setUnknown();
225  descriptions.addDefault(desc);
226 }
void addDefault(ParameterSetDescription const &psetDescription)

◆ produce()

void L1TrackFastJetProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 115 of file L1TrackFastJetProducer.cc.

References funct::abs(), coneSize_, callgraph::cs, deltaZ0Cut_, displaced_, doTightChi2_, edm::EventSetup::getData(), mps_fire::i, iEvent, createfilelist::int, L1TrackFastJetProducer_cfi::L1TrackFastJets, TTTrack< T >::momentum(), eostools::move(), PV3DBase< T, PVType, FrameType >::perp(), pvToken_, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, StripSubdetector::TID, TrackerTopology::tidRing(), StripSubdetector::TOB, TrackerTopology::tobLayer(), align::Tracker, trackToken_, trkBendChi2Max_, trkChi2dofMax_, trkChi2dofTightChi2_, trkEtaMax_, trkNPSStubMin_, trkNStubMin_, trkPtMin_, trkPtTightChi2_, trkZMax_, tTopoToken_, and TTTrack< T >::z0().

115  {
116  std::unique_ptr<TkJetCollection> L1TrackFastJets(new TkJetCollection);
117 
118  // L1 tracks
120  iEvent.getByToken(trackToken_, TTTrackHandle);
121  std::vector<TTTrack<Ref_Phase2TrackerDigi_> >::const_iterator iterL1Track;
122 
123  // Tracker Topology
124  const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
125 
126  edm::Handle<l1t::VertexCollection> TkPrimaryVertexHandle;
127  iEvent.getByToken(pvToken_, TkPrimaryVertexHandle);
128 
129  fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, coneSize_);
130  std::vector<fastjet::PseudoJet> JetInputs;
131 
132  float recoVtx = TkPrimaryVertexHandle->begin()->z0();
133  unsigned int this_l1track = 0;
134  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
135  this_l1track++;
136  float trk_pt = iterL1Track->momentum().perp();
137  float trk_z0 = iterL1Track->z0();
138  float trk_chi2dof = iterL1Track->chi2Red();
139  float trk_bendchi2 = iterL1Track->stubPtConsistency();
140  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_> >, TTStub<Ref_Phase2TrackerDigi_> > >
141  theStubs = iterL1Track->getStubRefs();
142  int trk_nstub = (int)theStubs.size();
143 
144  if (std::abs(trk_z0) > trkZMax_)
145  continue;
146  if (std::abs(iterL1Track->momentum().eta()) > trkEtaMax_)
147  continue;
148  if (trk_pt < trkPtMin_)
149  continue;
150  if (trk_nstub < trkNStubMin_)
151  continue;
152  if (trk_chi2dof > trkChi2dofMax_)
153  continue;
154  if (trk_bendchi2 > trkBendChi2Max_)
155  continue;
156  if (doTightChi2_ && (trk_pt > trkPtTightChi2_ && trk_chi2dof > trkChi2dofTightChi2_))
157  continue;
158 
159  int trk_nPS = 0;
160  for (int istub = 0; istub < trk_nstub; istub++) {
161  DetId detId(theStubs.at(istub)->getDetId());
162  bool tmp_isPS = false;
163  if (detId.det() == DetId::Detector::Tracker) {
164  if (detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3)
165  tmp_isPS = true;
166  else if (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9)
167  tmp_isPS = true;
168  }
169  if (tmp_isPS)
170  trk_nPS++;
171  }
172  if (trk_nPS < trkNPSStubMin_)
173  continue;
174  if (std::abs(recoVtx - trk_z0) > deltaZ0Cut_)
175  continue;
176 
177  fastjet::PseudoJet psuedoJet(iterL1Track->momentum().x(),
178  iterL1Track->momentum().y(),
179  iterL1Track->momentum().z(),
180  iterL1Track->momentum().mag());
181  JetInputs.push_back(psuedoJet); // input tracks for clustering
182  JetInputs.back().set_user_index(this_l1track - 1); // save track index in the collection
183  } // end loop over tracks
184 
185  fastjet::ClusterSequence cs(JetInputs, jet_def); // define the output jet collection
186  std::vector<fastjet::PseudoJet> JetOutputs =
187  fastjet::sorted_by_pt(cs.inclusive_jets(0)); // output jet collection, pT-ordered
188 
189  for (unsigned int ijet = 0; ijet < JetOutputs.size(); ++ijet) {
191  JetOutputs[ijet].px(), JetOutputs[ijet].py(), JetOutputs[ijet].pz(), JetOutputs[ijet].modp());
192  float sumpt = 0;
193  float avgZ = 0;
194  std::vector<edm::Ptr<L1TTTrackType> > L1TrackPtrs;
195  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(cs.constituents(JetOutputs[ijet]));
196 
197  for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
198  auto index = fjConstituents[i].user_index();
199  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, index);
200  L1TrackPtrs.push_back(trkPtr); // L1Tracks in the jet
201  sumpt = sumpt + trkPtr->momentum().perp();
202  avgZ = avgZ + trkPtr->momentum().perp() * trkPtr->z0();
203  }
204  avgZ = avgZ / sumpt;
206  TkJet trkJet(jetP4, jetRef, L1TrackPtrs, avgZ);
207  L1TrackFastJets->push_back(trkJet);
208  } //end loop over Jet Outputs
209 
210  if (displaced_)
211  iEvent.put(std::move(L1TrackFastJets), "L1TrackFastJetsExtended");
212  else
213  iEvent.put(std::move(L1TrackFastJets), "L1TrackFastJets");
214 }
unsigned int tobLayer(const DetId &id) const
edm::EDGetTokenT< VertexCollection > pvToken_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static constexpr auto TOB
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
Definition: DetId.h:17
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
unsigned int tidRing(const DetId &id) const
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ coneSize_

const double L1TrackFastJetProducer::coneSize_
private

Definition at line 75 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ deltaZ0Cut_

const double L1TrackFastJetProducer::deltaZ0Cut_
private

Definition at line 74 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ displaced_

const bool L1TrackFastJetProducer::displaced_
private

Definition at line 79 of file L1TrackFastJetProducer.cc.

Referenced by L1TrackFastJetProducer(), and produce().

◆ doTightChi2_

const bool L1TrackFastJetProducer::doTightChi2_
private

Definition at line 76 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ pvToken_

edm::EDGetTokenT<VertexCollection> L1TrackFastJetProducer::pvToken_
private

Definition at line 82 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trackToken_

const edm::EDGetTokenT<std::vector<TTTrack<Ref_Phase2TrackerDigi_> > > L1TrackFastJetProducer::trackToken_
private

Definition at line 81 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkBendChi2Max_

const double L1TrackFastJetProducer::trkBendChi2Max_
private

Definition at line 69 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkChi2dofMax_

const float L1TrackFastJetProducer::trkChi2dofMax_
private

Definition at line 68 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkChi2dofTightChi2_

const float L1TrackFastJetProducer::trkChi2dofTightChi2_
private

Definition at line 78 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkEtaMax_

const float L1TrackFastJetProducer::trkEtaMax_
private

Definition at line 71 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkNPSStubMin_

const int L1TrackFastJetProducer::trkNPSStubMin_
private

Definition at line 73 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkNStubMin_

const int L1TrackFastJetProducer::trkNStubMin_
private

Definition at line 72 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkPtMin_

const float L1TrackFastJetProducer::trkPtMin_
private

Definition at line 70 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkPtTightChi2_

const float L1TrackFastJetProducer::trkPtTightChi2_
private

Definition at line 77 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ trkZMax_

const float L1TrackFastJetProducer::trkZMax_
private

Definition at line 67 of file L1TrackFastJetProducer.cc.

Referenced by produce().

◆ tTopoToken_

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> L1TrackFastJetProducer::tTopoToken_
private

Definition at line 83 of file L1TrackFastJetProducer.cc.

Referenced by produce().