Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:11:11

0001 //-------------------------------------------------
0002 //
0003 /**  \class DTTrigProd
0004  *     Main EDProducer for the DTTPG
0005  *
0006  *
0007  *
0008  *   \author C. Battilana
0009  *
0010  */
0011 //
0012 //--------------------------------------------------
0013 
0014 // Framework related classes
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016 #include "FWCore/Framework/interface/stream/EDProducer.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/Frameworkfwd.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 #include "FWCore/Framework/interface/EventSetup.h"
0021 #include "FWCore/Framework/interface/ESHandle.h"
0022 #include "FWCore/Utilities/interface/ESGetToken.h"
0023 #include "FWCore/Framework/interface/Run.h"
0024 #include "FWCore/Utilities/interface/EDPutToken.h"
0025 
0026 #include "L1TriggerConfig/DTTPGConfig/interface/DTConfigManager.h"
0027 #include "L1TriggerConfig/DTTPGConfig/interface/DTConfigManagerRcd.h"
0028 #include "L1Trigger/DTTrigger/interface/DTTrig.h"
0029 
0030 // Data Formats classes
0031 #include "L1Trigger/DTSectorCollector/interface/DTSectCollPhSegm.h"
0032 #include "L1Trigger/DTSectorCollector/interface/DTSectCollThSegm.h"
0033 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
0034 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhDigi.h"
0035 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
0036 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThDigi.h"
0037 
0038 // Collaborating classes
0039 #include <iostream>
0040 
0041 using namespace edm;
0042 using namespace std;
0043 
0044 // DataFormats interface
0045 typedef vector<DTSectCollPhSegm> SectCollPhiColl;
0046 typedef SectCollPhiColl::const_iterator SectCollPhiColl_iterator;
0047 typedef vector<DTSectCollThSegm> SectCollThetaColl;
0048 typedef SectCollThetaColl::const_iterator SectCollThetaColl_iterator;
0049 
0050 class DTTrigProd : public edm::stream::EDProducer<> {
0051 public:
0052   //! Constructor
0053   DTTrigProd(const edm::ParameterSet& pset);
0054 
0055   //! Create Trigger Units before starting event processing
0056   void beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) override;
0057 
0058   //! Producer: process every event and generates trigger data
0059   void produce(edm::Event& iEvent, const edm::EventSetup& iEventSetup) override;
0060 
0061 private:
0062   // Trigger istance
0063   DTTrig my_trig;
0064 
0065   edm::EDPutTokenT<L1MuDTChambPhContainer> phToken_;
0066   edm::EDPutTokenT<L1MuDTChambThContainer> thToken_;
0067   edm::ESGetToken<DTConfigManager, DTConfigManagerRcd> dtConfigToken_;
0068 
0069   // Trigger Configuration Manager CCB validity flag
0070   bool my_CCBValid = false;
0071 
0072   // Sector Format Flag true=[0-11] false=[1-12]
0073   const bool my_DTTFnum;
0074 
0075   // Debug Flag
0076   const bool my_debug;
0077 
0078   // Lut dump file parameters
0079   const bool my_lut_dump_flag;
0080   const short int my_lut_btic;
0081 };
0082 
0083 DTTrigProd::DTTrigProd(const ParameterSet& pset)
0084     : my_trig(pset, consumesCollector()),
0085       phToken_{produces<L1MuDTChambPhContainer>()},
0086       thToken_{produces<L1MuDTChambThContainer>()},
0087       dtConfigToken_(esConsumes<DTConfigManager, DTConfigManagerRcd, edm::Transition::BeginRun>()),
0088       my_DTTFnum{pset.getParameter<bool>("DTTFSectorNumbering")},
0089       my_debug{pset.getUntrackedParameter<bool>("debug")},
0090       my_lut_dump_flag{pset.getUntrackedParameter<bool>("lutDumpFlag")},
0091       my_lut_btic{static_cast<short int>(pset.getUntrackedParameter<int>("lutBtic"))} {}
0092 
0093 void DTTrigProd::beginRun(edm::Run const& iRun, const edm::EventSetup& iEventSetup) {
0094   if (my_debug)
0095     cout << "DTTrigProd::beginRun  " << iRun.id().run() << endl;
0096 
0097   ESHandle<DTConfigManager> dtConfig = iEventSetup.getHandle(dtConfigToken_);
0098 
0099   my_CCBValid = dtConfig->CCBConfigValidity();
0100 
0101   my_trig.createTUs(iEventSetup);
0102   if (my_debug)
0103     cout << "[DTTrigProd] TU's Created" << endl;
0104 
0105   if (my_lut_dump_flag) {
0106     cout << "Dumping luts...." << endl;
0107     my_trig.dumpLuts(my_lut_btic, dtConfig.product());
0108   }
0109 }
0110 
0111 void DTTrigProd::produce(Event& iEvent, const EventSetup& iEventSetup) {
0112   vector<L1MuDTChambPhDigi> outPhi;
0113   vector<L1MuDTChambThDigi> outTheta;
0114 
0115   // SV check if CCB configuration is valid, otherwise just produce empty collections
0116   if (!my_CCBValid) {
0117     if (my_debug)
0118       cout << "[DTTrigProd] CCB configuration is not valid for this run, empty collection will be produced " << endl;
0119   } else {
0120     my_trig.triggerReco(iEvent, iEventSetup);
0121     // BX offset used to correct DTTPG output
0122     int bx_offset = my_trig.getBXOffset();
0123 
0124     if (my_debug)
0125       cout << "[DTTrigProd] Trigger algorithm run for " << iEvent.id() << endl;
0126 
0127     // Convert Phi Segments
0128     SectCollPhiColl myPhiSegments;
0129     myPhiSegments = my_trig.SCPhTrigs();
0130 
0131     SectCollPhiColl_iterator SCPCend = myPhiSegments.end();
0132     for (SectCollPhiColl_iterator it = myPhiSegments.begin(); it != SCPCend; ++it) {
0133       int step = (*it).step() - bx_offset;  // Shift correct BX to 0 (needed for DTTF data processing)
0134       int sc_sector = (*it).SCId().sector();
0135       if (my_DTTFnum == true)
0136         sc_sector--;  // Modified for DTTF numbering [0-11]
0137       outPhi.push_back(L1MuDTChambPhDigi(step,
0138                                          (*it).ChamberId().wheel(),
0139                                          sc_sector,
0140                                          (*it).ChamberId().station(),
0141                                          (*it).phi(),
0142                                          (*it).phiB(),
0143                                          (*it).code(),
0144                                          !(*it).isFirst(),
0145                                          0));
0146     }
0147 
0148     // Convert Theta Segments
0149     SectCollThetaColl myThetaSegments;
0150     myThetaSegments = my_trig.SCThTrigs();
0151 
0152     SectCollThetaColl_iterator SCTCend = myThetaSegments.end();
0153     for (SectCollThetaColl_iterator it = myThetaSegments.begin(); it != SCTCend; ++it) {
0154       int pos[7], qual[7];
0155       for (int i = 0; i < 7; i++) {
0156         pos[i] = (*it).position(i);
0157         qual[i] = (*it).quality(i);
0158       }
0159       int step = (*it).step() - bx_offset;  // Shift correct BX to 0 (needed for DTTF data processing)
0160       int sc_sector = (*it).SCId().sector();
0161       if (my_DTTFnum == true)
0162         sc_sector--;  // Modified for DTTF numbering [0-11]
0163       outTheta.push_back(
0164           L1MuDTChambThDigi(step, (*it).ChamberId().wheel(), sc_sector, (*it).ChamberId().station(), pos, qual));
0165     }
0166   }
0167 
0168   // Write everything into the event (CB write empty collection as default actions if emulator does not run)
0169   iEvent.emplace(phToken_, std::move(outPhi));
0170   iEvent.emplace(thToken_, std::move(outTheta));
0171 }
0172 
0173 DEFINE_FWK_MODULE(DTTrigProd);