File indexing completed on 2025-01-07 03:06:33
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include <memory>
0012
0013
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/stream/EDProducer.h"
0016 #include "FWCore/Framework/interface/Event.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include "FWCore/Utilities/interface/InputTag.h"
0020 #include "FWCore/Utilities/interface/transform.h"
0021 #include "FWCore/Utilities/interface/StreamID.h"
0022 #include "DataFormats/Math/interface/Point3D.h"
0023 #include "DataFormats/Common/interface/DetSetVector.h"
0024 #include "DataFormats/DetId/interface/DetId.h"
0025 #include "DataFormats/Phase2TrackerDigi/interface/Phase2ITDigiHit.h"
0026 #include "DataFormats/Phase2TrackerDigi/interface/Phase2ITQCore.h"
0027 #include "DataFormats/Phase2TrackerDigi/interface/Phase2ITChip.h"
0028 #include "DataFormats/Phase2TrackerDigi/interface/Phase2ITChipBitStream.h"
0029 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0030 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0031 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0032 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0033 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0034
0035 class Phase2ITQCoreProducer : public edm::stream::EDProducer<> {
0036 public:
0037 Phase2ITQCoreProducer(const edm::ParameterSet&);
0038 ~Phase2ITQCoreProducer() override = default;
0039
0040 private:
0041 void produce(edm::Event&, const edm::EventSetup&) override;
0042
0043 const edm::InputTag src_;
0044 const edm::EDGetTokenT<edm::DetSetVector<PixelDigi>> pixelDigi_token_;
0045 const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoToken_;
0046 typedef math::XYZPointD Point;
0047 typedef std::vector<Point> PointCollection;
0048 };
0049
0050 Phase2ITQCoreProducer::Phase2ITQCoreProducer(const edm::ParameterSet& iConfig)
0051 : src_(iConfig.getParameter<edm::InputTag>("src")),
0052 pixelDigi_token_(consumes(iConfig.getParameter<edm::InputTag>("siPixelDigi"))),
0053 tTopoToken_(esConsumes()) {
0054 produces<edm::DetSetVector<Phase2ITQCore>>();
0055 produces<edm::DetSetVector<Phase2ITChipBitStream>>();
0056 }
0057
0058 namespace {
0059
0060
0061
0062
0063 constexpr int kQCoresInChipRow = (672);
0064 constexpr int kQCoresInChipColumn = (216);
0065 constexpr int kQCoresInChipRowGap = (5);
0066 constexpr int kQCoresInChipColumnGap = (10);
0067 }
0068
0069 void updateHitCoordinatesForLargePixels(Phase2ITDigiHit& hit) {
0070
0071
0072
0073
0074
0075
0076
0077
0078 int row = hit.row();
0079 int col = hit.col();
0080
0081
0082 int updated_row = 0;
0083 int updated_col = 0;
0084
0085
0086 if (row < kQCoresInChipRow) {
0087 updated_row = row;
0088 } else if (row < (kQCoresInChipRow + kQCoresInChipRowGap)) {
0089 updated_row = kQCoresInChipRow - 1;
0090 }
0091 else if (row < (kQCoresInChipRow + 2 * kQCoresInChipRowGap)) {
0092 updated_row = kQCoresInChipRow;
0093 } else {
0094 updated_row = (hit.row() - 2 * kQCoresInChipRowGap);
0095 }
0096
0097
0098 if (col < kQCoresInChipColumn) {
0099 updated_col = col;
0100 } else if (col < kQCoresInChipColumn + kQCoresInChipColumnGap) {
0101 updated_col = kQCoresInChipColumn - kQCoresInChipColumnGap;
0102 } else if (col < (kQCoresInChipColumn + 2 * kQCoresInChipColumn)) {
0103 updated_col = kQCoresInChipColumn;
0104 } else {
0105 updated_col = (hit.col() - 2 * kQCoresInChipColumnGap);
0106 }
0107
0108 hit.set_row(updated_row);
0109 hit.set_col(updated_col);
0110 }
0111
0112 void adjustEdges(std::vector<Phase2ITDigiHit> hitList) {
0113
0114
0115
0116 std::for_each(hitList.begin(), hitList.end(), &updateHitCoordinatesForLargePixels);
0117 }
0118
0119 std::vector<Phase2ITChip> splitByChip(const std::vector<Phase2ITDigiHit>& hitList) {
0120
0121 std::array<std::vector<Phase2ITDigiHit>, 4> hits_per_chip;
0122 for (auto hit : hitList) {
0123 int chip_index = (hit.col() < kQCoresInChipColumn) ? 0 : 1;
0124 if (hit.row() >= kQCoresInChipRow) {
0125 chip_index += 2;
0126 }
0127 hits_per_chip[chip_index].push_back(hit);
0128 }
0129
0130
0131 std::vector<Phase2ITChip> chips;
0132 chips.reserve(4);
0133 for (int chip_index = 0; chip_index < 4; chip_index++) {
0134 chips.push_back(Phase2ITChip(chip_index, hits_per_chip[chip_index]));
0135 }
0136
0137 return chips;
0138 }
0139
0140 std::vector<Phase2ITChip> processHits(std::vector<Phase2ITDigiHit> hitList) {
0141 adjustEdges(hitList);
0142 std::vector<Phase2ITChip> chips = splitByChip(hitList);
0143
0144 return chips;
0145 }
0146
0147
0148 void Phase2ITQCoreProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0149 using namespace edm;
0150 using namespace std;
0151
0152 unique_ptr<edm::DetSetVector<Phase2ITQCore>> aQCoreVector = make_unique<edm::DetSetVector<Phase2ITQCore>>();
0153 unique_ptr<edm::DetSetVector<Phase2ITChipBitStream>> aBitStreamVector =
0154 make_unique<edm::DetSetVector<Phase2ITChipBitStream>>();
0155
0156 auto const& tTopo = iSetup.getData(tTopoToken_);
0157
0158 auto pixelDigiHandle = iEvent.getHandle(pixelDigi_token_);
0159 const edm::DetSetVector<PixelDigi>& pixelDigi = *pixelDigiHandle;
0160 for (const auto& theDigis : pixelDigi) {
0161 DetId tkId = theDigis.id;
0162 std::vector<Phase2ITDigiHit> hitlist;
0163 std::vector<int> id;
0164
0165 if (tkId.subdetId() == PixelSubdetector::PixelBarrel) {
0166 int layer_num = tTopo.pxbLayer(tkId.rawId());
0167 int ladder_num = tTopo.pxbLadder(tkId.rawId());
0168 int module_num = tTopo.pxbModule(tkId.rawId());
0169 id = {tkId.subdetId(), layer_num, ladder_num, module_num};
0170 } else if (tkId.subdetId() == PixelSubdetector::PixelEndcap) {
0171 int module_num = tTopo.pxfModule(tkId());
0172 int disk_num = tTopo.pxfDisk(tkId());
0173 int blade_num = tTopo.pxfBlade(tkId());
0174 int panel_num = tTopo.pxfPanel(tkId());
0175 int side_num = tTopo.pxfSide(tkId());
0176 id = {tkId.subdetId(), module_num, disk_num, blade_num, panel_num, side_num};
0177 }
0178
0179 for (const auto& digi : theDigis) {
0180 hitlist.emplace_back(digi.row(), digi.column(), digi.adc());
0181 }
0182
0183 std::vector<Phase2ITChip> chips = processHits(std::move(hitlist));
0184
0185 DetSet<Phase2ITQCore> DetSetQCores(tkId);
0186 DetSet<Phase2ITChipBitStream> DetSetBitStream(tkId);
0187
0188 for (size_t i = 0; i < chips.size(); i++) {
0189 Phase2ITChip chip = chips[i];
0190 std::vector<Phase2ITQCore> qcores = chip.get_organized_QCores();
0191 for (auto& qcore : qcores) {
0192 DetSetQCores.push_back(qcore);
0193 }
0194
0195 Phase2ITChipBitStream aChipBitStream(i, chip.get_chip_code());
0196 DetSetBitStream.push_back(aChipBitStream);
0197 }
0198
0199 aBitStreamVector->insert(DetSetBitStream);
0200 aQCoreVector->insert(DetSetQCores);
0201 }
0202
0203 iEvent.put(std::move(aQCoreVector));
0204 iEvent.put(std::move(aBitStreamVector));
0205 }
0206
0207 DEFINE_FWK_MODULE(Phase2ITQCoreProducer);