CMS 3D CMS Logo

PhysCandPacker.cc
Go to the documentation of this file.
3 
4 #include "CaloTokens.h"
5 #include "PhysCandPacker.h"
6 
7 template<typename T, typename F>
9 process(unsigned int id1, unsigned int id2, const BXVector<T>& coll, F filter)
10 {
11  std::vector<uint32_t> load[2];
12 
13  for (int i = coll.getFirstBX(); i <= coll.getLastBX(); ++i) {
14  uint16_t jetbit[4] = {0, 0, 0, 0};
15  int n = 0;
16  for (auto j = coll.begin(i); j != coll.end(i) && n < 4; ++j) {
17  if (!filter(*j))
18  continue;
19  //std::cout << j->hwPt() << " @ " << j->hwEta() << ", " << j->hwPhi() << " > " << j->hwQual() << " > " << j->hwIso() << std::endl;
20  jetbit[n++] = std::min(j->hwPt(), 0x3F) |
21  (abs(j->hwEta()) & 0x7) << 6 |
22  ((j->hwEta() >> 3) & 0x1) << 9 |
23  (j->hwPhi() & 0x1F) << 10;
24  }
25  uint32_t word0=(jetbit[0] & 0xFFFF);
26  uint32_t word1=(jetbit[1] & 0xFFFF);
27  uint32_t word2=(jetbit[2] & 0xFFFF);
28  uint32_t word3=(jetbit[3] & 0xFFFF);
29 
30  load[0].push_back(word0);
31  load[0].push_back(word2);
32 
33  load[1].push_back(word1);
34  load[1].push_back(word3);
35 
36  }
37 
38  return {l1t::Block(id1, load[0]),l1t::Block(id2, load[1])};
39 }
40 
41 namespace l1t {
42  namespace stage1 {
43  Blocks
45  {
47  event.getByToken(static_cast<const CaloTokens*>(toks)->getEGammaToken(), egammas);
48 
49  return process(85,87, *egammas, [](const l1t::EGamma& eg) -> bool { return eg.hwIso() == 1; });
50  }
51 
52  Blocks
54  {
56  event.getByToken(static_cast<const CaloTokens*>(toks)->getEGammaToken(), egammas);
57 
58  return process(89,91, *egammas, [](const l1t::EGamma& eg) -> bool { return eg.hwIso() == 0; });
59  }
60 
61  Blocks
63  {
65  event.getByToken(static_cast<const CaloTokens*>(toks)->getJetToken(), jets);
66 
67  return process(77,79, *jets, [](const l1t::Jet& jet) -> bool { return !(jet.hwQual() & 2); });
68  }
69 
70  Blocks
72  {
74  event.getByToken(static_cast<const CaloTokens*>(toks)->getJetToken(), jets);
75 
76  return process(81,83, *jets, [](const l1t::Jet& jet) -> bool { return jet.hwQual() & 2; });
77  }
78 
79  Blocks
81  {
83  event.getByToken(static_cast<const CaloTokens*>(toks)->getTauToken(), taus);
84 
85  return process(101,103, *taus, [](const l1t::Tau& tau) -> bool { return true; });
86  }
87 
88  Blocks
90  {
92  event.getByToken(static_cast<const CaloTokens*>(toks)->getIsoTauToken(), taus);
93 
94  return process(105,107, *taus, [](const l1t::Tau& tau) -> bool { return true; });
95  }
96  }
97 }
98 
const_iterator end(int bx) const
Definition: Tau.h:16
Blocks pack(const edm::Event &, const PackerTokens *) override
delete x;
Definition: CaloConfig.h:22
int hwIso() const
Definition: L1Candidate.h:52
Definition: Jet.h:16
Blocks pack(const edm::Event &, const PackerTokens *) override
vector< PseudoJet > jets
std::vector< Block > Blocks
Definition: Block.h:72
Blocks pack(const edm::Event &, const PackerTokens *) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
int hwQual() const
Definition: L1Candidate.h:51
int getFirstBX() const
JetCorrectorParametersCollection coll
Definition: classes.h:10
#define DEFINE_L1T_PACKER(type)
Definition: PackerFactory.h:22
def load(fileName)
Definition: svgfig.py:546
Blocks pack(const edm::Event &, const PackerTokens *) override
Blocks pack(const edm::Event &, const PackerTokens *) override
int getLastBX() const
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
l1t::Blocks process(unsigned int id1, unsigned int id2, const BXVector< T > &coll, F filter)
const_iterator begin(int bx) const
Definition: event.py:1
Blocks pack(const edm::Event &, const PackerTokens *) override