CMS 3D CMS Logo

InputMakerPhase2.cc
Go to the documentation of this file.
1 /*
2  * InputMakerPhase2.cpp
3  *
4  * Created on: May 20, 2020
5  * Author: kbunkow
6  */
7 
12 
13 #include <iostream>
14 
17  event.getByToken(inputTokenDtPh, dtPhDigis);
18  event.getByToken(inputTokenDtTh, dtThDigis);
19 }
20 
22  unsigned int iProcessor,
23  l1t::tftype procTyp,
24  int bxFrom,
25  int bxTo,
26  std::vector<std::unique_ptr<IOMTFEmulationObserver> >& observers) {
27  if (!dtPhDigis)
28  return;
29 
30  boost::property_tree::ptree procDataTree;
31 
32  for (const auto& digiIt : *dtPhDigis->getContainer()) {
33  DTChamberId detid(digiIt.whNum(), digiIt.stNum(), digiIt.scNum() + 1);
34 
36  if (!acceptDigi(detid, iProcessor, procTyp))
37  continue;
38 
39  // HACK for Phase-2 (DT TPs are centered in bX=20)
40  if (digiIt.bxNum() - 20 >= bxFrom && digiIt.bxNum() - 20 <= bxTo) {
41  addDTphiDigi(muonStubsInLayers, digiIt, dtThDigis.product(), iProcessor, procTyp);
42 
43  auto& dtP2Digi = procDataTree.add_child("dtP2Digi", boost::property_tree::ptree());
44  dtP2Digi.add("<xmlattr>.whNum", digiIt.whNum());
45  dtP2Digi.add("<xmlattr>.scNum", digiIt.scNum());
46  dtP2Digi.add("<xmlattr>.stNum", digiIt.stNum());
47  dtP2Digi.add("<xmlattr>.slNum", digiIt.slNum());
48  dtP2Digi.add("<xmlattr>.quality", digiIt.quality());
49  dtP2Digi.add("<xmlattr>.rpcFlag", digiIt.rpcFlag());
50  dtP2Digi.add("<xmlattr>.phi", digiIt.phi());
51  dtP2Digi.add("<xmlattr>.phiBend", digiIt.phiBend());
52  }
53  }
54 
55  if (!mergePhiAndTheta) {
56  for (auto& thetaDigi : (*(dtThDigis->getContainer()))) {
57  if (thetaDigi.bxNum() >= bxFrom && thetaDigi.bxNum() <= bxTo) {
58  addDTetaStubs(muonStubsInLayers, thetaDigi, iProcessor, procTyp);
59  }
60  }
61  }
62 
63  for (auto& obs : observers)
64  obs->addProcesorData("linkData", procDataTree);
65 }
66 
67 //dtThDigis is provided as argument, because in the OMTF implementation the phi and eta digis are merged (even thought it is artificial)
69  const L1Phase2MuDTPhDigi& digi,
70  const L1MuDTChambThContainer* dtThDigis,
71  unsigned int iProcessor,
72  l1t::tftype procTyp) {
73  DTChamberId detid(digi.whNum(), digi.stNum(), digi.scNum() + 1);
74 
75  MuonStub stub;
76 
77  //converting the quality to the same encoding as in phase-1, as it is important for extrapolation
78  if (digi.quality() >= 6)
79  stub.qualityHw = digi.quality() - 2;
80  else if (digi.quality() >= 3) {
81  if (digi.slNum() == 3)
82  stub.qualityHw = 3;
83  else if (digi.slNum() == 1)
84  stub.qualityHw = 2;
85  } else {
86  if (digi.slNum() == 3)
87  stub.qualityHw = 1;
88  else if (digi.slNum() == 1)
89  stub.qualityHw = 0;
90  }
91 
92  if (stub.qualityHw < config.getMinDtPhiQuality())
93  return;
94 
95  unsigned int hwNumber = config.getLayerNumber(detid.rawId());
96  if (config.getHwToLogicLayer().find(hwNumber) == config.getHwToLogicLayer().end())
97  return;
98 
99  auto iter = config.getHwToLogicLayer().find(hwNumber);
100  unsigned int iLayer = iter->second;
101  unsigned int iInput = OMTFinputMaker::getInputNumber(&config, detid.rawId(), iProcessor, procTyp);
102  //MuonStub& stub = muonStubsInLayers[iLayer][iInput];
103 
104  stub.type = MuonStub::DT_PHI_ETA;
105 
106  stub.phiHw = angleConverter.getProcessorPhi(
107  OMTFinputMaker::getProcessorPhiZero(&config, iProcessor), procTyp, digi.scNum(), digi.phi());
108 
109  //TODO the dtThDigis are not good yet,so passing an empty container to the angleConverter
110  //then it should return middle of chambers
111  //remove when the dtThDigis are fixed on the DT side
112  L1MuDTChambThContainer dtThDigisEmpty;
113  stub.etaHw = angleConverter.getGlobalEta(detid, &dtThDigisEmpty, digi.bxNum() - 20);
114  //in phase2, the phiB is 13 bits, and range is [-2, 2 rad] so 4 rad, 2^13 units/(4 rad) = 1^11/rad.
115  //need to convert them to 512units==1rad (to use OLD PATTERNS...)
116  stub.phiBHw = digi.phiBend() * config.dtPhiBUnitsRad() / 2048;
117  //the cut if (stub.qualityHw >= config.getMinDtPhiBQuality()) is done in the ProcessorBase<GoldenPatternType>::restrictInput
118  //as is is done like that in the firmware
119 
120  // need to shift 20-BX to roll-back the shift introduced by the DT TPs
121  stub.bx = digi.bxNum() - 20;
122  //stub.timing = digi.getTiming(); //TODO what about sub-bx timing, is is available?
123 
124  stub.logicLayer = iLayer;
125  stub.detId = detid;
126 
127  OmtfName board(iProcessor, &config);
128  LogTrace("l1tOmtfEventPrint") << board.name() << " L1Phase2MuDTPhDigi: detid " << detid << " digi "
129  << " whNum " << digi.whNum() << " scNum " << digi.scNum() << " stNum " << digi.stNum()
130  << " slNum " << digi.slNum() << " quality " << digi.quality() << " rpcFlag "
131  << digi.rpcFlag() << " phi " << digi.phi() << " phiBend " << digi.phiBend()
132  << std::endl;
133  OMTFinputMaker::addStub(&config, muonStubsInLayers, iLayer, iInput, stub);
134 }
135 
137  const L1MuDTChambThDigi& thetaDigi,
138  unsigned int iProcessor,
139  l1t::tftype procTyp) {
140  //in the Phase1 omtf the theta stubs are merged with the phi in the addDTphiDigi
141  //TODO implement if needed
142 }
143 
145  unsigned int iProcessor,
146  l1t::tftype procType) {
147  return OMTFinputMaker::acceptDtDigi(&config, dTChamberId, iProcessor, procType);
148 }
149 
151  MuStubsInputTokens& muStubsInputTokens,
152  edm::EDGetTokenT<L1Phase2MuDTPhContainer> inputTokenDTPhPhase2,
153  const OMTFConfiguration* config,
154  std::unique_ptr<OmtfAngleConverter> angleConverter)
155  : OMTFinputMaker(edmParameterSet, muStubsInputTokens, config, std::move(angleConverter)) {
156  edm::LogImportant("OMTFReconstruction") << "constructing InputMakerPhase2" << std::endl;
157 
158  if (edmParameterSet.exists("usePhase2DTPrimitives") && edmParameterSet.getParameter<bool>("usePhase2DTPrimitives")) {
159  if (edmParameterSet.getParameter<bool>("dropDTPrimitives") != true)
160  throw cms::Exception(
161  "L1TMuonOverlapPhase2 InputMakerPhase2::InputMakerPhase2 usePhase2DTPrimitives is true, but dropDTPrimitives "
162  "is not true");
163  //if the Phase2DTPrimitives are used, then the phase1 DT primitives should be dropped
164  edm::LogImportant("OMTFReconstruction") << " using Phase2 DT trigger primitives" << std::endl;
165  digiToStubsConverters.emplace_back(std::make_unique<DtPhase2DigiToStubsConverterOmtf>(
166  config, this->angleConverter.get(), inputTokenDTPhPhase2, muStubsInputTokens.inputTokenDtTh));
167  }
168 }
virtual bool acceptDigi(const DTChamberId &dTChamberId, unsigned int iProcessor, l1t::tftype procType)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool acceptDigi(const DTChamberId &dTChamberId, unsigned int iProcessor, l1t::tftype procType) override
void makeStubs(MuonStubPtrs2D &muonStubsInLayers, unsigned int iProcessor, l1t::tftype procTyp, int bxFrom, int bxTo, std::vector< std::unique_ptr< IOMTFEmulationObserver > > &observers) override
virtual int getGlobalEta(const DTChamberId dTChamberId, const L1MuDTChambThContainer *dtThDigis, int bxNum) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
T const * product() const
Definition: Handle.h:70
Definition: config.py:1
std::vector< MuonStubPtrs1D > MuonStubPtrs2D
Definition: MuonStub.h:69
int qualityHw
Definition: MuonStub.h:47
#define LogTrace(id)
const OmtfAngleConverter & angleConverter
InputMakerPhase2(const edm::ParameterSet &edmParameterSet, MuStubsInputTokens &muStubsInputTokens, edm::EDGetTokenT< L1Phase2MuDTPhContainer > inputTokenDTPhPhase2, const OMTFConfiguration *config, std::unique_ptr< OmtfAngleConverter > angleConverter)
The_Container const * getContainer() const
edm::Handle< L1Phase2MuDTPhContainer > dtPhDigis
static int getProcessorPhiZero(const OMTFConfiguration *config, unsigned int iProcessor)
edm::EDGetTokenT< L1MuDTChambThContainer > inputTokenDtTh
Log< level::Error, true > LogImportant
virtual int getProcessorPhi(int phiZero, l1t::tftype part, int dtScNum, int dtPhi) const
void addDTphiDigi(MuonStubPtrs2D &muonStubsInLayers, const L1Phase2MuDTPhDigi &digi, const L1MuDTChambThContainer *dtThDigis, unsigned int iProcessor, l1t::tftype procTyp) override
std::vector< std::unique_ptr< DigiToStubsConverterBase > > digiToStubsConverters
virtual void addDTphiDigi(MuonStubPtrs2D &muonStubsInLayers, const L1Phase2MuDTPhDigi &digi, const L1MuDTChambThContainer *dtThDigis, unsigned int iProcessor, l1t::tftype procTyp)=0
edm::Handle< L1MuDTChambThContainer > dtThDigis
edm::EDGetTokenT< L1MuDTChambThContainer > inputTokenDtTh
unsigned int getInputNumber(unsigned int rawId, unsigned int iProcessor, l1t::tftype type)
void loadDigis(const edm::Event &event) override
Segment_Container const * getContainer() const
std::string name() const
Definition: OmtfName.cc:51
edm::EDGetTokenT< L1Phase2MuDTPhContainer > inputTokenDtPh
static bool acceptDtDigi(const OMTFConfiguration *config, const DTChamberId &dTChamberId, unsigned int iProcessor, l1t::tftype procType)
static void addStub(const OMTFConfiguration *config, MuonStubPtrs2D &muonStubsInLayers, unsigned int iLayer, unsigned int iInput, MuonStub &stub)
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
void addDTetaStubs(MuonStubPtrs2D &muonStubsInLayers, const L1MuDTChambThDigi &thetaDigi, unsigned int iProcessor, l1t::tftype procTyp) override
virtual void addDTetaStubs(MuonStubPtrs2D &muonStubsInLayers, const L1MuDTChambThDigi &thetaDigi, unsigned int iProcessor, l1t::tftype procTyp)=0
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1