Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:47

0001 #include <memory>
0002 
0003 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
0004 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
0005 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
0006 #include "L1Trigger/L1TMuon/interface/deprecate/DTCollector.h"
0007 
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 
0013 using namespace L1TMuon;
0014 
0015 DTCollector::DTCollector(const edm::ParameterSet& ps)
0016     : SubsystemCollector(ps), bx_min(ps.getParameter<int>("BX_min")), bx_max(ps.getParameter<int>("BX_max")) {
0017   if (ps.getParameter<bool>("runBunchCrossingCleaner")) {
0018     const edm::ParameterSet& bxccfg = ps.getParameterSet("bxCleanerCfg");
0019     _bxc = std::make_unique<DTBunchCrossingCleaner>(bxccfg);
0020   } else {
0021     _bxc.reset(nullptr);
0022   }
0023 }
0024 
0025 void DTCollector::extractPrimitives(const edm::Event& ev,
0026                                     const edm::EventSetup& es,
0027                                     TriggerPrimitiveCollection& out) const {
0028   TriggerPrimitiveCollection cleaned, temp, chamb_list;
0029   edm::Handle<L1MuDTChambPhContainer> phiDigis;
0030   edm::Handle<L1MuDTChambThContainer> thetaDigis;
0031   ev.getByLabel(_src, phiDigis);
0032   ev.getByLabel(_src, thetaDigis);
0033   for (int wheel = -2; wheel <= 2; ++wheel) {
0034     for (int station = 1; station <= 4; ++station) {
0035       for (int sector = 0; sector <= 11; ++sector) {
0036         chamb_list.clear();
0037         for (int bx = bx_min; bx <= bx_max; ++bx) {
0038           std::unique_ptr<const L1MuDTChambPhDigi> phi_segm_1(phiDigis->chPhiSegm1(wheel, station, sector, bx));
0039           std::unique_ptr<const L1MuDTChambPhDigi> phi_segm_2(phiDigis->chPhiSegm2(wheel, station, sector, bx));
0040           std::unique_ptr<const L1MuDTChambThDigi> theta_segm(thetaDigis->chThetaSegm(wheel, station, sector, bx));
0041 
0042           int bti_group_1 = -1, bti_group_2 = -1;
0043 
0044           if (theta_segm) {
0045             bti_group_1 = findBTIGroupForThetaDigi(*theta_segm, 1);
0046             bti_group_2 = findBTIGroupForThetaDigi(*theta_segm, 2);
0047           }
0048 
0049           if (phi_segm_1 && bti_group_1 != -1) {
0050             chamb_list.push_back(processDigis(*phi_segm_1, *theta_segm, bti_group_1));
0051           } else if (phi_segm_1 && bti_group_1 == -1) {
0052             chamb_list.push_back(processDigis(*phi_segm_1, 1));
0053           } else if (!phi_segm_1 && bti_group_1 != -1) {
0054             chamb_list.push_back(processDigis(*theta_segm, bti_group_1));
0055           }
0056 
0057           if (phi_segm_2 && bti_group_2 != -1) {
0058             chamb_list.push_back(processDigis(*phi_segm_2, *theta_segm, bti_group_2));
0059           } else if (phi_segm_2 && bti_group_2 == -1) {
0060             chamb_list.push_back(processDigis(*phi_segm_2, 2));
0061           } else if (!phi_segm_2 && bti_group_2 != -1) {
0062             chamb_list.push_back(processDigis(*phi_segm_2, bti_group_2));
0063           }
0064 
0065           phi_segm_1.release();
0066           phi_segm_2.release();
0067           theta_segm.release();
0068         }
0069         if (_bxc) {
0070           temp = _bxc->clean(chamb_list);
0071           cleaned.insert(cleaned.end(), temp.begin(), temp.end());
0072         } else {
0073           cleaned.insert(cleaned.end(), chamb_list.begin(), chamb_list.end());
0074         }
0075       }
0076     }
0077   }
0078   out.insert(out.end(), cleaned.begin(), cleaned.end());
0079 }
0080 
0081 TriggerPrimitive DTCollector::processDigis(const L1MuDTChambPhDigi& digi, const int& segment_number) const {
0082   DTChamberId detid(digi.whNum(), digi.stNum(), digi.scNum() + 1);
0083   return TriggerPrimitive(detid, digi, segment_number);
0084 }
0085 
0086 TriggerPrimitive DTCollector::processDigis(const L1MuDTChambThDigi& digi_th, const int bti_group) const {
0087   DTChamberId detid(digi_th.whNum(), digi_th.stNum(), digi_th.scNum() + 1);
0088   return TriggerPrimitive(detid, digi_th, bti_group);
0089 }
0090 
0091 TriggerPrimitive DTCollector::processDigis(const L1MuDTChambPhDigi& digi_phi,
0092                                            const L1MuDTChambThDigi& digi_theta,
0093                                            const int bti_group) const {
0094   DTChamberId detid(digi_phi.whNum(), digi_phi.stNum(), digi_phi.scNum() + 1);
0095   return TriggerPrimitive(detid, digi_phi, digi_theta, bti_group);
0096 }
0097 
0098 int DTCollector::findBTIGroupForThetaDigi(const L1MuDTChambThDigi& digi, const int pos) const {
0099   //if( digi.stNum() == 4 ) return -1; // there is no theta layer there
0100   int result = -1;
0101   for (int i = 0; i < 7; ++i) {
0102     if (digi.position(i) == pos)
0103       result = i;
0104   }
0105   return result;
0106 }
0107 
0108 #include "L1Trigger/L1TMuon/interface/deprecate/SubsystemCollectorFactory.h"
0109 DEFINE_EDM_PLUGIN(SubsystemCollectorFactory, DTCollector, "DTCollector");