CMS 3D CMS Logo

L1TCaloRCTToUpgradeConverter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1Trigger/L1TCalorimeter
4 // Class: L1TCaloRCTToUpgradeConverter
5 //
13 //
14 // Original Author: James Brooke
15 // Created: Thu, 05 Dec 2013 17:39:27 GMT
16 //
17 //
28 
32 
34 
37 
39 
40 #include <vector>
41 
42 namespace l1t {
43 
45  public:
47 
48  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
49 
50  private:
51  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
52 
53  // ----------member data ---------------------------
54 
57 
60  };
61 
62 } // namespace l1t
63 
64 using namespace l1t;
65 
67  : rgnToken_{consumes<L1CaloRegionCollection>(ps.getParameter<edm::InputTag>("regionTag"))},
68  emToken_{consumes<L1CaloEmCollection>(ps.getParameter<edm::InputTag>("emTag"))},
69  rgnPutToken_{produces<CaloRegionBxCollection>()},
70  emPutToken_{produces<CaloEmCandBxCollection>()} {}
71 
72 // ------------ method called to produce the data ------------
74  // check status of RCT conditions & renew if needed
75 
76  // store new formats
77  BXVector<CaloEmCand> emcands;
79 
80  // get old formats
81  auto const& ems = iEvent.get(emToken_);
82  auto const& rgns = iEvent.get(rgnToken_);
83 
84  // get the firstBx_ and lastBx_ from the input datatypes (assume bx for em same as rgn)
85  int firstBx = 0;
86  int lastBx = 0;
87  for (auto const& em : ems) {
88  int bx = em.bx();
89  if (bx < firstBx)
90  firstBx = bx;
91  if (bx > lastBx)
92  lastBx = bx;
93  }
94 
95  emcands.setBXRange(firstBx, lastBx);
96  regions.setBXRange(firstBx, lastBx);
97 
98  const ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > p4(0, 0, 0, 0);
99 
100  // loop over EM
101  for (auto const& em : ems) {
102  // get physical units
103  // double pt = 0.;
104  // double eta = 0.;
105  // double phi = 0.;
106  //math::PtEtaPhiMLorentzVector p4( pt+1.e-6, eta, phi, 0. );
107 
108  //CaloStage1Cluster cluster;
109  CaloEmCand EmCand(p4, (int)em.rank(), (int)em.regionId().ieta(), (int)em.regionId().iphi(), (int)em.index());
110 
111  EmCand.setHwIso((int)em.isolated());
112  //std::cout<<"ISO: "<<EmCand.hwIso()<<" "<<em.isolated()<<std::endl;
113 
114  // create new format
115  emcands.push_back(em.bx(), EmCand);
116  }
117 
118  // loop over regions
119  for (auto const& rgn : rgns) {
120  // get physical units
121  // double pt = 0.;
122  // double eta = 0.;
123  // double phi = 0.;
124  //math::PtEtaPhiMLorentzVector p4( pt+1.e-6, eta, phi, 0 );
125 
126  bool tauVeto = rgn.fineGrain(); //equivalent to tauVeto for HB/HE, includes extra info for HF
127  int hwQual = (int)tauVeto;
128 
129  // create new format
130  // several values here are stage 2 only, leave empty
131  CaloRegion region(p4, // LorentzVector& p4,
132  0., // etEm,
133  0., // etHad,
134  (int)rgn.et(), // pt,
135  (int)rgn.id().ieta(), // eta,
136  (int)rgn.id().iphi(), // phi,
137  hwQual, // qual,
138  0, // hwEtEm,
139  0); // hwEtHad
140 
141  // add to output
142  regions.push_back(rgn.bx(), region);
143  }
144 
145  iEvent.emplace(emPutToken_, std::move(emcands));
147 }
148 
149 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
152  desc.add<edm::InputTag>("regionTag");
153  desc.add<edm::InputTag>("emTag");
154  descriptions.addDefault(desc);
155 }
156 
157 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< L1CaloRegionCollection > const rgnToken_
delete x;
Definition: CaloConfig.h:22
L1TCaloRCTToUpgradeConverter(const edm::ParameterSet &ps)
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDPutTokenT< CaloRegionBxCollection > const rgnPutToken_
edm::EDGetTokenT< L1CaloEmCollection > const emToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setHwIso(int iso)
Definition: L1Candidate.h:32
void setBXRange(int bxFirst, int bxLast)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
edm::EDPutTokenT< CaloEmCandBxCollection > const emPutToken_
def move(src, dest)
Definition: eostools.py:511
void push_back(int bx, T object)