CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
L1TMuonBarrelKalmanTrackProducer.cc
Go to the documentation of this file.
1 #include <memory>
4 
7 
10 
15 
16 //
17 // class declaration
18 //
19 
21 public:
24 
25  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
26 
27 private:
28  void beginStream(edm::StreamID) override;
29  void produce(edm::Event&, const edm::EventSetup&) override;
30  void endStream() override;
32  std::vector<int> bx_;
35 };
37  : src_(consumes<std::vector<L1MuKBMTCombinedStub> >(iConfig.getParameter<edm::InputTag>("src"))),
38  bx_(iConfig.getParameter<std::vector<int> >("bx")),
39  algo_(new L1TMuonBarrelKalmanAlgo(iConfig.getParameter<edm::ParameterSet>("algoSettings"))),
40  trackFinder_(new L1TMuonBarrelKalmanTrackFinder(iConfig.getParameter<edm::ParameterSet>("trackFinderSettings"))) {
41  produces<L1MuKBMTrackBxCollection>();
42  produces<l1t::RegionalMuonCandBxCollection>("BMTF");
43 }
44 
46  if (algo_ != nullptr)
47  delete algo_;
48 
49  if (trackFinder_ != nullptr)
50  delete trackFinder_;
51 
52  // do anything here that needs to be done at destruction time
53  // (e.g. close files, deallocate resources etc.)
54 }
55 
56 //
57 // member functions
58 //
59 
60 // ------------ method called to produce the data ------------
62  using namespace edm;
64  iEvent.getByToken(src_, stubHandle);
65 
67  for (uint i = 0; i < stubHandle->size(); ++i) {
68  L1MuKBMTCombinedStubRef r(stubHandle, i);
69  stubs.push_back(r);
70  }
71 
72  std::unique_ptr<l1t::RegionalMuonCandBxCollection> outBMTF(new l1t::RegionalMuonCandBxCollection());
73  std::unique_ptr<L1MuKBMTrackBxCollection> out(new L1MuKBMTrackBxCollection());
74  outBMTF->setBXRange(bx_.front(), bx_.back());
75  out->setBXRange(bx_.front(), bx_.back());
76  for (const auto& bx : bx_) {
78  for (const auto& track : tmp) {
79  out->push_back(bx, track);
80  algo_->addBMTFMuon(bx, track, outBMTF);
81  }
82  }
83  iEvent.put(std::move(outBMTF), "BMTF");
84  iEvent.put(std::move(out));
85 }
86 
87 // ------------ method called once each stream before processing any runs, lumis or events ------------
89 
90 // ------------ method called once each stream after processing all runs, lumis and events ------------
92 
94  //The following says we do not know what parameters are allowed so do no validation
95  // Please change this to state exactly what you do use, even if it is no parameters
97  desc.setUnknown();
98  descriptions.addDefault(desc);
99 }
100 
101 //define this as a plug-in
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
L1MuKBMTrackCollection process(L1TMuonBarrelKalmanAlgo *, const L1MuKBMTCombinedStubRefVector &stubs, int bx)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
L1TMuonBarrelKalmanTrackFinder * trackFinder_
std::vector< edm::Ref< L1MuKBMTCombinedStubCollection > > L1MuKBMTCombinedStubRefVector
L1TMuonBarrelKalmanTrackProducer(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< L1MuKBMTCombinedStub > > src_
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
def move
Definition: eostools.py:511
BXVector< L1MuKBMTrack > L1MuKBMTrackBxCollection
Definition: L1MuKBMTrack.h:17
void addBMTFMuon(int, const L1MuKBMTrack &, std::unique_ptr< l1t::RegionalMuonCandBxCollection > &)
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
tmp
align.sh
Definition: createJobs.py:716