@@ -38,9 +38,7 @@ class Phase2ITQCoreProducer : public edm::stream::EDProducer<> {
3838 ~Phase2ITQCoreProducer () override = default ;
3939
4040private:
41- virtual void beginJob (const edm::EventSetup&);
4241 virtual void produce (edm::Event&, const edm::EventSetup&);
43- virtual void endJob ();
4442
4543 const edm::InputTag src_;
4644 const edm::EDGetTokenT<edm::DetSetVector<PixelDigi>> pixelDigi_token_;
@@ -51,8 +49,8 @@ class Phase2ITQCoreProducer : public edm::stream::EDProducer<> {
5149
5250Phase2ITQCoreProducer::Phase2ITQCoreProducer (const edm::ParameterSet& iConfig)
5351 : src_(iConfig.getParameter<edm::InputTag>(" src" )),
54- pixelDigi_token_(consumes<edm::DetSetVector<PixelDigi>> (iConfig.getParameter<edm::InputTag>(" siPixelDigi" ))),
55- tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd> ()) {
52+ pixelDigi_token_(consumes(iConfig.getParameter<edm::InputTag>(" siPixelDigi" ))),
53+ tTopoToken_(esConsumes()) {
5654 produces<edm::DetSetVector<Phase2ITQCore>>();
5755 produces<edm::DetSetVector<Phase2ITChipBitStream>>();
5856}
@@ -62,13 +60,13 @@ namespace {
6260 // Dimension for 4 chips module = 1354 X 434 = (672 + 5 + 672 + 5) X (216 + 1 + 216 + 1)
6361 // Spacing 1 in column and 5 in row is introduced for each chip in between
6462 // if neighboring chip exists
65- const int kQCoresInChipRow = (672 );
66- const int kQCoresInChipColumn = (216 );
67- const int kQCoresInChipRowGap = (5 );
68- const int kQCoresInChipColumnGap = (10 );
63+ constexpr int kQCoresInChipRow = (672 );
64+ constexpr int kQCoresInChipColumn = (216 );
65+ constexpr int kQCoresInChipRowGap = (5 );
66+ constexpr int kQCoresInChipColumnGap = (10 );
6967} // namespace
7068
71- Phase2ITDigiHit updateHitCoordinatesForLargePixels (Phase2ITDigiHit& hit) {
69+ void updateHitCoordinatesForLargePixels (Phase2ITDigiHit& hit) {
7270 /*
7371 In-place modification of Hit coordinates to take into account large pixels
7472 Hits corresponding to large pixels are remapped so they lie on the boundary of the chip
@@ -109,21 +107,18 @@ Phase2ITDigiHit updateHitCoordinatesForLargePixels(Phase2ITDigiHit& hit) {
109107
110108 hit.set_row (updated_row);
111109 hit.set_col (updated_col);
112-
113- return hit;
114110}
115111
116- std::vector<Phase2ITDigiHit> adjustEdges (std::vector<Phase2ITDigiHit> hitList) {
112+ void adjustEdges (std::vector<Phase2ITDigiHit> hitList) {
117113 /*
118114 In-place modification of Hit coordinates to take into account large pixels
119115 */
120116 std::for_each (hitList.begin (), hitList.end (), &updateHitCoordinatesForLargePixels);
121- return hitList;
122117}
123118
124- std::vector<Phase2ITChip> splitByChip (std::vector<Phase2ITDigiHit> hitList) {
119+ std::vector<Phase2ITChip> splitByChip (const std::vector<Phase2ITDigiHit>& hitList) {
125120 // Split the hit list by read out chip
126- std::vector <std::vector<Phase2ITDigiHit>> hits_per_chip ( 4 ) ;
121+ std::array <std::vector<Phase2ITDigiHit>, 4 > hits_per_chip;
127122 for (auto hit : hitList) {
128123 int chip_index = (hit.col () < kQCoresInChipColumn ) ? 0 : 1 ;
129124 if (hit.row () >= kQCoresInChipRow ) {
@@ -134,6 +129,7 @@ std::vector<Phase2ITChip> splitByChip(std::vector<Phase2ITDigiHit> hitList) {
134129
135130 // Generate Phase2ITChip objects from the hit lists
136131 std::vector<Phase2ITChip> chips;
132+ chips.reserve (4 );
137133 for (int chip_index = 0 ; chip_index < 4 ; chip_index++) {
138134 chips.push_back (Phase2ITChip (chip_index, hits_per_chip[chip_index]));
139135 }
@@ -142,16 +138,8 @@ std::vector<Phase2ITChip> splitByChip(std::vector<Phase2ITDigiHit> hitList) {
142138}
143139
144140std::vector<Phase2ITChip> processHits (std::vector<Phase2ITDigiHit> hitList) {
145- std::vector<Phase2ITDigiHit> newHitList;
146- newHitList = adjustEdges (hitList);
147-
148- std::vector<Phase2ITChip> chips = splitByChip (newHitList);
149- std::vector<bool > code;
150-
151- for (size_t i = 0 ; i < chips.size (); i++) {
152- Phase2ITChip chip = chips[i];
153- code = chip.get_chip_code ();
154- }
141+ adjustEdges (hitList);
142+ std::vector<Phase2ITChip> chips = splitByChip (hitList);
155143
156144 return chips;
157145}
@@ -167,13 +155,10 @@ void Phase2ITQCoreProducer::produce(edm::Event& iEvent, const edm::EventSetup& i
167155
168156 auto const & tTopo = iSetup.getData (tTopoToken_);
169157
170- edm::Handle<edm::DetSetVector<PixelDigi>> pixelDigiHandle;
171- iEvent.getByToken (pixelDigi_token_, pixelDigiHandle);
172-
173- edm::DetSetVector<PixelDigi>::const_iterator iterDet;
174- for (iterDet = pixelDigiHandle->begin (); iterDet != pixelDigiHandle->end (); iterDet++) {
175- DetId tkId = iterDet->id ;
176- edm::DetSet<PixelDigi> theDigis = (*pixelDigiHandle)[tkId];
158+ auto pixelDigiHandle = iEvent.getHandle (pixelDigi_token_);
159+ const edm::DetSetVector<PixelDigi>& pixelDigi = *pixelDigiHandle;
160+ for (const auto & theDigis : pixelDigi) {
161+ DetId tkId = theDigis.id ;
177162 std::vector<Phase2ITDigiHit> hitlist;
178163 std::vector<int > id;
179164
@@ -191,11 +176,11 @@ void Phase2ITQCoreProducer::produce(edm::Event& iEvent, const edm::EventSetup& i
191176 id = {tkId.subdetId (), module_num, disk_num, blade_num, panel_num, side_num};
192177 }
193178
194- for (auto iterDigi = theDigis. begin (); iterDigi != theDigis. end (); ++iterDigi ) {
195- hitlist.emplace_back (Phase2ITDigiHit (iterDigi-> row (), iterDigi-> column (), iterDigi-> adc () ));
179+ for (const auto & digi : theDigis) {
180+ hitlist.emplace_back (digi. row (), digi. column (), digi. adc ());
196181 }
197182
198- std::vector<Phase2ITChip> chips = processHits (hitlist);
183+ std::vector<Phase2ITChip> chips = processHits (std::move ( hitlist) );
199184
200185 DetSet<Phase2ITQCore> DetSetQCores (tkId);
201186 DetSet<Phase2ITChipBitStream> DetSetBitStream (tkId);
@@ -219,8 +204,4 @@ void Phase2ITQCoreProducer::produce(edm::Event& iEvent, const edm::EventSetup& i
219204 iEvent.put (std::move (aBitStreamVector));
220205}
221206
222- void Phase2ITQCoreProducer::beginJob (edm::EventSetup const &) {}
223-
224- void Phase2ITQCoreProducer::endJob () {}
225-
226207DEFINE_FWK_MODULE (Phase2ITQCoreProducer);
0 commit comments