Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-04-15 22:41:37

0001 #include "L1Trigger/DTTriggerPhase2/interface/ShowerCandidate.h"
0002 #include "L1Trigger/DTTriggerPhase2/interface/ShowerBuilder.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 
0005 using namespace cmsdt;
0006 
0007 // ============================================================================
0008 // Constructors and destructor
0009 // ============================================================================
0010 ShowerBuilder::ShowerBuilder(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
0011     :  // Unpack information from pset
0012       showerTaggingAlgo_(pset.getParameter<int>("showerTaggingAlgo")),
0013       threshold_for_shower_(pset.getParameter<int>("threshold_for_shower")),
0014       nHits_per_bx_(pset.getParameter<int>("nHits_per_bx")),
0015       obdt_hits_bxpersistence_(pset.getParameter<int>("obdt_hits_bxpersistence")),
0016       obdt_wire_relaxing_time_(pset.getParameter<int>("obdt_wire_relaxing_time")),
0017       bmtl1_hits_bxpersistence_(pset.getParameter<int>("bmtl1_hits_bxpersistence")),
0018       debug_(pset.getUntrackedParameter<bool>("debug")),
0019       scenario_(pset.getParameter<int>("scenario")) {}
0020 
0021 // ============================================================================
0022 // Main methods (initialise, run, finish)
0023 // ============================================================================
0024 void ShowerBuilder::initialise(const edm::EventSetup &iEventSetup) {}
0025 
0026 void ShowerBuilder::run(edm::Event &iEvent,
0027                         const edm::EventSetup &iEventSetup,
0028                         const DTDigiCollection &digis,
0029                         ShowerCandidatePtr &showerCandidate_SL1,
0030                         ShowerCandidatePtr &showerCandidate_SL3) {
0031   // Clear auxiliars
0032   clear();
0033   // Set the incoming hits in the channels
0034   setInChannels(&digis);
0035 
0036   std::map<int, ShowerCandidatePtr> aux_showerCands{// defined as a map to easy acces with SL number
0037                                                     {1, std::make_shared<ShowerCandidate>()},
0038                                                     {3, std::make_shared<ShowerCandidate>()}};
0039 
0040   int nHits = all_hits.size();
0041   if (nHits != 0) {
0042     if (debug_)
0043       LogDebug("ShowerBuilder") << "        - Going to study " << nHits << " hits";
0044     if (showerTaggingAlgo_ == 0) {
0045       // Standalone mode: just save hits and flag if the number of hits is above the threshold
0046       processHits_standAlone(aux_showerCands);
0047     } else if (showerTaggingAlgo_ == 1) {
0048       // Firmware emulation:
0049       // mimics the behavior of sending and receiving hits from the OBDT to the shower algorithm in the BMTL1.
0050       processHitsFirmwareEmulation(aux_showerCands);
0051     }
0052   } else {
0053     if (debug_)
0054       LogDebug("ShowerBuilder") << "        - No hits to study.";
0055   }
0056 
0057   showerCandidate_SL1 = std::move(aux_showerCands[1]);
0058   showerCandidate_SL3 = std::move(aux_showerCands[3]);
0059 }
0060 
0061 // ============================================================================
0062 // Auxiliary methods
0063 // ============================================================================
0064 void ShowerBuilder::clear() {
0065   all_hits.clear();
0066   all_hits_perBx.clear();
0067   showerb::buffer_reset(obdt_buffer);
0068   showerb::buffer_reset(hot_wires_buffer);
0069   showerb::buffer_reset(bmtl1_sl1_buffer);
0070   showerb::buffer_reset(bmtl1_sl3_buffer);
0071 }
0072 
0073 void ShowerBuilder::setInChannels(const DTDigiCollection *digis) {
0074   for (const auto &dtLayerId_It : *digis) {
0075     const DTLayerId id = dtLayerId_It.first;
0076     if (id.superlayer() == 2)
0077       continue;  // Skip SL2 digis
0078     // Now iterate over the digis
0079     for (DTDigiCollection::const_iterator digiIt = (dtLayerId_It.second).first; digiIt != (dtLayerId_It.second).second;
0080          ++digiIt) {
0081       auto dtpAux = DTPrimitive();
0082       dtpAux.setTDCTimeStamp((*digiIt).time());
0083       dtpAux.setChannelId((*digiIt).wire());
0084       dtpAux.setLayerId(id.layer());
0085       dtpAux.setSuperLayerId(id.superlayer());
0086       dtpAux.setCameraId(id.rawId());
0087       all_hits.push_back(dtpAux);
0088     }
0089   }
0090   // sort the hits by time
0091   std::stable_sort(all_hits.begin(), all_hits.end(), showerb::hitTimeSort_shower);
0092 }
0093 
0094 void ShowerBuilder::processHits_standAlone(std::map<int, ShowerCandidatePtr> &showerCands) {
0095   // For each superlayer, fill the buffer and check if nHits >= threshold
0096   for (auto &hit : all_hits) {
0097     showerb::DTPrimPlusBx _hitpbx(-1, hit);
0098     if (hit.superLayerId() == 1)
0099       bmtl1_sl1_buffer.push_back(_hitpbx);
0100     else if (hit.superLayerId() == 3)
0101       bmtl1_sl3_buffer.push_back(_hitpbx);
0102   }
0103 
0104   if (triggerShower(bmtl1_sl1_buffer)) {
0105     if (debug_) {
0106       int nHits_sl1 = bmtl1_sl1_buffer.size();
0107       LogDebug("ShowerBuilder") << "        o Shower found in SL1 with " << nHits_sl1 << " hits";
0108     }
0109     showerCands[1]->flag();
0110     set_shower_properties(showerCands[1], bmtl1_sl1_buffer);
0111   }
0112   if (triggerShower(bmtl1_sl3_buffer)) {
0113     if (debug_) {
0114       int nHits_sl3 = bmtl1_sl3_buffer.size();
0115       LogDebug("ShowerBuilder") << "        o Shower found in SL3 with " << nHits_sl3 << " hits";
0116     }
0117     showerCands[3]->flag();
0118     set_shower_properties(showerCands[3], bmtl1_sl3_buffer);
0119   }
0120 }
0121 
0122 void ShowerBuilder::processHitsFirmwareEmulation(std::map<int, ShowerCandidatePtr> &showerCands) {
0123   // Use a for over BXs to emulate the behavior of the OBDT and BMTL1, this means considering:
0124   // 1. OBDT can only store recived hits during 4 BXs (adjustable with obdt_hits_bxpersistence_)
0125   // 2. OBDT can recive hits from the same wire during 2 BXs (adjustable with obdt_wire_relaxing_time_)
0126   // 3. OBDT can only sends 8 hits per BX to the BMTL1 (adjustable with obdt_hits_per_bx_)
0127   // 4. Shower algorithm in the BMTL1 mantains the hits along 16 BXs and triggers if nHits >= threshold (adjustable with bmtl1_hits_bxpersistence_)
0128   groupHits_byBx();
0129 
0130   // auxiliary variables
0131   int prev_nHits_sl1 = 0;
0132   int nHits_sl1 = 0;
0133   bool shower_sl1_already_set = false;
0134   int prev_nHits_sl3 = 0;
0135   int nHits_sl3 = 0;
0136   bool shower_sl3_already_set = false;
0137   // -------------------
0138 
0139   int min_bx = all_hits_perBx.begin()->first;
0140   int max_bx = all_hits_perBx.rbegin()->first;
0141 
0142   if (debug_)
0143     LogDebug("ShowerBuilder") << "        - bx range: " << min_bx << " - " << max_bx << ", size: " << (max_bx - min_bx);
0144 
0145   for (int bx = min_bx; bx <= max_bx + 17; bx++) {
0146     fill_obdt(bx);
0147 
0148     if (debug_)
0149       LogDebug("ShowerBuilder") << "          ^ " << obdt_buffer.size() << " hits in obdt_buffer at BX " << bx;
0150     fill_bmtl1_buffers();
0151 
0152     nHits_sl1 = bmtl1_sl1_buffer.size();
0153     if (debug_)
0154       LogDebug("ShowerBuilder") << "          ^ " << nHits_sl1 << " hits in bmtl1_sl1_buffer at BX " << bx;
0155     nHits_sl3 = bmtl1_sl3_buffer.size();
0156     if (debug_)
0157       LogDebug("ShowerBuilder") << "          ^ " << nHits_sl3 << " hits in bmtl1_sl3_buffer at BX " << bx;
0158     if (triggerShower(bmtl1_sl1_buffer)) {
0159       if (debug_ && !showerCands[1]->isFlagged()) {
0160         LogDebug("ShowerBuilder") << "        o Shower found in SL1 with " << nHits_sl1 << " hits";
0161       }
0162       showerCands[1]->flag();
0163     }
0164     if (triggerShower(bmtl1_sl3_buffer)) {
0165       if (debug_ && !showerCands[3]->isFlagged()) {
0166         LogDebug("ShowerBuilder") << "        o Shower found in SL3 with " << nHits_sl3 << " hits";
0167       }
0168       showerCands[3]->flag();
0169     }
0170 
0171     if (nHits_sl1 < prev_nHits_sl1 && showerCands[1]->isFlagged() && !shower_sl1_already_set) {
0172       shower_sl1_already_set = true;
0173       set_shower_properties(showerCands[1], bmtl1_sl1_buffer, nHits_sl1, bx);
0174     } else {
0175       prev_nHits_sl1 = nHits_sl1;
0176     }
0177 
0178     if (nHits_sl3 < prev_nHits_sl3 && showerCands[3]->isFlagged() && !shower_sl3_already_set) {
0179       shower_sl3_already_set = true;
0180       set_shower_properties(showerCands[3], bmtl1_sl3_buffer, nHits_sl3, bx);
0181     } else {
0182       prev_nHits_sl3 = nHits_sl3;
0183     }
0184 
0185     bxStep(bx);
0186   }
0187 }
0188 
0189 bool ShowerBuilder::triggerShower(const showerb::ShowerBuffer &buffer) {
0190   int nHits = buffer.size();
0191   if (nHits >= threshold_for_shower_) {
0192     return true;
0193   }
0194   return false;
0195 }
0196 
0197 void ShowerBuilder::set_shower_properties(ShowerCandidatePtr &showerCand,
0198                                           showerb::ShowerBuffer &buffer,
0199                                           int nhits,
0200                                           int bx,
0201                                           int min_wire,
0202                                           int max_wire,
0203                                           float avg_pos,
0204                                           float avg_time) {
0205   DTPrimitives _hits;
0206 
0207   showerb::buffer_get_hits(buffer, _hits);
0208   std::stable_sort(_hits.begin(), _hits.end(), showerb::hitWireSort_shower);
0209 
0210   // change values considering the buffer or the input values
0211   nhits = (nhits == -1) ? _hits.size() : nhits;
0212   min_wire = (min_wire == -1) ? _hits.front().channelId() : min_wire;
0213   max_wire = (max_wire == -1) ? _hits.back().channelId() : max_wire;
0214   avg_pos = (avg_pos == -1) ? showerb::compute_avg_pos(_hits) : avg_pos;
0215   avg_time = (avg_time == -1) ? showerb::compute_avg_time(_hits) : avg_time;
0216 
0217   showerCand->setNhits(nhits);
0218   showerCand->setBX(bx);
0219   showerCand->setMinWire(min_wire);
0220   showerCand->setMaxWire(max_wire);
0221   showerCand->setAvgPos(avg_pos);
0222   showerCand->setAvgTime(avg_time);
0223   showerb::set_wires_profile(showerCand->getWiresProfile(), _hits);
0224 }
0225 
0226 void ShowerBuilder::groupHits_byBx() {
0227   double temp_shift = 0;
0228   if (scenario_ == MC)         //scope for MC
0229     temp_shift = 400;          // value used in standard CMSSW simulation
0230   else if (scenario_ == DATA)  //scope for data
0231     temp_shift = 0;
0232   else if (scenario_ == SLICE_TEST)  //scope for slice test
0233     temp_shift = 400;                // slice test to mimic simulation
0234 
0235   static const double shift_back = temp_shift;
0236 
0237   all_hits_perBx.clear();
0238   // Group hits by BX
0239   for (auto &hit : all_hits) {
0240     // Compute the BX from the TDC time
0241     int bx = hit.tdcTimeStamp() / 25 - shift_back;
0242     all_hits_perBx[bx].push_back(hit);
0243   }
0244 }
0245 
0246 void ShowerBuilder::fill_obdt(const int bx) {
0247   // Fill the OBDT buffer with the hits in the current BX this function ensure that hot wires are not added
0248   if (all_hits_perBx.find(bx) != all_hits_perBx.end()) {
0249     for (auto &hit : all_hits_perBx[bx]) {
0250       if (debug_)
0251         LogDebug("ShowerBuilder") << "          ^ Trying to add hit with wire " << hit.channelId() << " in OBDT at BX "
0252                                   << bx;
0253       if (!showerb::buffer_contains(hot_wires_buffer, hit)) {
0254         if (debug_)
0255           LogDebug("ShowerBuilder") << "          ^ added";
0256         showerb::DTPrimPlusBx _hit(bx, hit);
0257         obdt_buffer.push_back(_hit);
0258         hot_wires_buffer.push_back(_hit);
0259       }
0260     }
0261   }
0262 }
0263 
0264 void ShowerBuilder::fill_bmtl1_buffers() {
0265   // Fill the BMTL1 buffer with the hits in the OBDT buffer only nHits_per_bx_ hits are added
0266   if (obdt_buffer.empty())
0267     return;
0268   if (debug_)
0269     LogDebug("ShowerBuilder") << "          ^ Getting hits from OBDT";
0270   for (int i = 0; i < nHits_per_bx_; i++) {
0271     if (obdt_buffer.empty())
0272       break;
0273     auto _hitpbx = obdt_buffer.front();
0274     if (_hitpbx.second.superLayerId() == 1) {
0275       bmtl1_sl1_buffer.push_back(_hitpbx);
0276     } else if (_hitpbx.second.superLayerId() == 3) {
0277       bmtl1_sl3_buffer.push_back(_hitpbx);
0278     }
0279     obdt_buffer.pop_front();
0280   }
0281 }
0282 
0283 void ShowerBuilder::bxStep(const int _current_bx) {
0284   // Remove old elements from the buffers
0285   showerb::buffer_clear_olds(obdt_buffer, _current_bx, obdt_hits_bxpersistence_);
0286   showerb::buffer_clear_olds(hot_wires_buffer, _current_bx, obdt_hits_bxpersistence_);
0287   showerb::buffer_clear_olds(bmtl1_sl1_buffer, _current_bx, bmtl1_hits_bxpersistence_);
0288   showerb::buffer_clear_olds(bmtl1_sl3_buffer, _current_bx, bmtl1_hits_bxpersistence_);
0289 }