File indexing completed on 2024-09-07 04:37:02
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include <memory>
0012
0013
0014 #include "DataFormats/Common/interface/RefVector.h"
0015 #include "FWCore/Framework/interface/Frameworkfwd.h"
0016 #include "FWCore/Framework/interface/stream/EDProducer.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/Framework/interface/EventSetup.h"
0020 #include "DataFormats/Common/interface/Handle.h"
0021 #include "FWCore/Utilities/interface/InputTag.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0024 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0025 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0026 #include "DataFormats/Math/interface/LorentzVector.h"
0027
0028
0029 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0030 #include "DataFormats/L1TCorrelator/interface/TkJet.h"
0031 #include "DataFormats/L1TCorrelator/interface/TkJetFwd.h"
0032 #include "DataFormats/L1Trigger/interface/Vertex.h"
0033
0034
0035 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0036 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0037
0038 #include <fastjet/JetDefinition.hh>
0039
0040 #include <string>
0041 #include "TMath.h"
0042 #include "TH1.h"
0043
0044 using namespace l1t;
0045 using namespace edm;
0046 using namespace std;
0047
0048
0049
0050
0051
0052
0053
0054 class L1TrackFastJetProducer : public edm::stream::EDProducer<> {
0055 public:
0056 typedef TTTrack<Ref_Phase2TrackerDigi_> L1TTTrackType;
0057 typedef std::vector<L1TTTrackType> L1TTTrackCollectionType;
0058 typedef edm::RefVector<L1TTTrackCollectionType> L1TTTrackRefCollectionType;
0059
0060 explicit L1TrackFastJetProducer(const edm::ParameterSet&);
0061 ~L1TrackFastJetProducer() override;
0062
0063 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0064
0065 private:
0066 void produce(edm::Event&, const edm::EventSetup&) override;
0067
0068
0069 const double coneSize_;
0070 const bool displaced_;
0071
0072 const EDGetTokenT<L1TTTrackRefCollectionType> trackToken_;
0073 };
0074
0075
0076 L1TrackFastJetProducer::L1TrackFastJetProducer(const edm::ParameterSet& iConfig)
0077 : coneSize_((float)iConfig.getParameter<double>("coneSize")),
0078 displaced_(iConfig.getParameter<bool>("displaced")),
0079 trackToken_(consumes<L1TTTrackRefCollectionType>(iConfig.getParameter<InputTag>("L1TrackInputTag"))) {
0080 if (displaced_)
0081 produces<TkJetCollection>("L1TrackFastJetsExtended");
0082 else
0083 produces<TkJetCollection>("L1TrackFastJets");
0084 }
0085
0086
0087 L1TrackFastJetProducer::~L1TrackFastJetProducer() {}
0088
0089
0090 void L1TrackFastJetProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0091 std::unique_ptr<TkJetCollection> L1TrackFastJets(new TkJetCollection);
0092
0093
0094 edm::Handle<L1TTTrackRefCollectionType> TTTrackHandle;
0095 iEvent.getByToken(trackToken_, TTTrackHandle);
0096
0097 fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, coneSize_);
0098 std::vector<fastjet::PseudoJet> JetInputs;
0099
0100 for (unsigned int this_l1track = 0; this_l1track < TTTrackHandle->size(); this_l1track++) {
0101 edm::Ptr<L1TTTrackType> iterL1Track(TTTrackHandle, this_l1track);
0102
0103 fastjet::PseudoJet psuedoJet(iterL1Track->momentum().x(),
0104 iterL1Track->momentum().y(),
0105 iterL1Track->momentum().z(),
0106 iterL1Track->momentum().mag());
0107 JetInputs.push_back(psuedoJet);
0108 JetInputs.back().set_user_index(this_l1track);
0109 }
0110
0111 fastjet::ClusterSequence cs(JetInputs, jet_def);
0112 std::vector<fastjet::PseudoJet> JetOutputs =
0113 fastjet::sorted_by_pt(cs.inclusive_jets(0));
0114
0115 for (unsigned int ijet = 0; ijet < JetOutputs.size(); ++ijet) {
0116 math::XYZTLorentzVector jetP4(
0117 JetOutputs[ijet].px(), JetOutputs[ijet].py(), JetOutputs[ijet].pz(), JetOutputs[ijet].modp());
0118 float sumpt = 0;
0119 float avgZ = 0;
0120 std::vector<edm::Ptr<L1TTTrackType> > L1TrackPtrs;
0121 std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(cs.constituents(JetOutputs[ijet]));
0122
0123 for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
0124 auto index = fjConstituents[i].user_index();
0125 edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, index);
0126 L1TrackPtrs.push_back(trkPtr);
0127 sumpt = sumpt + trkPtr->momentum().perp();
0128 avgZ = avgZ + trkPtr->momentum().perp() * trkPtr->z0();
0129 }
0130 avgZ = avgZ / sumpt;
0131 edm::Ref<JetBxCollection> jetRef;
0132 TkJet trkJet(jetP4, jetRef, L1TrackPtrs, avgZ);
0133 L1TrackFastJets->push_back(trkJet);
0134 }
0135
0136 if (displaced_)
0137 iEvent.put(std::move(L1TrackFastJets), "L1TrackFastJetsExtended");
0138 else
0139 iEvent.put(std::move(L1TrackFastJets), "L1TrackFastJets");
0140 }
0141
0142 void L1TrackFastJetProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0143
0144 edm::ParameterSetDescription desc;
0145 desc.add<edm::InputTag>(
0146 "L1PVertexInputTag",
0147 edm::InputTag("l1tTrackVertexAssociationProducerForJets", "Level1TTTracksSelectedAssociated"));
0148 desc.add<double>("coneSize", 0.5);
0149 desc.add<bool>("displaced", false);
0150 }
0151
0152
0153 DEFINE_FWK_MODULE(L1TrackFastJetProducer);