CMS 3D CMS Logo

L1TMuonBarrelKalmanTrackProducer.cc
Go to the documentation of this file.
1 #include <memory>
4 
7 
10 
11 
16 
17 //
18 // class declaration
19 //
20 
22  public:
25 
26  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
27 
28  private:
29  void beginStream(edm::StreamID) override;
30  void produce(edm::Event&, const edm::EventSetup&) override;
31  void endStream() override;
33  std::vector<int> bx_;
34  std::unique_ptr<L1TMuonBarrelKalmanAlgo> algo_;
35  std::unique_ptr<L1TMuonBarrelKalmanTrackFinder> trackFinder_;
36 
37 
38 
39 
40 
41 };
43  src_(consumes<std::vector<L1MuKBMTCombinedStub> >(iConfig.getParameter<edm::InputTag>("src"))),
44  bx_(iConfig.getParameter<std::vector<int> >("bx")),
45  algo_(new L1TMuonBarrelKalmanAlgo(iConfig.getParameter<edm::ParameterSet>("algoSettings"))),
46  trackFinder_(new L1TMuonBarrelKalmanTrackFinder(iConfig.getParameter<edm::ParameterSet>("trackFinderSettings")))
47 {
48  produces <std::vector<L1MuKBMTrack> >();
49  produces <l1t::RegionalMuonCandBxCollection>("BMTF");
50 
51 }
52 
53 
55 {
56 
57  // do anything here that needs to be done at destruction time
58  // (e.g. close files, deallocate resources etc.)
59 
60 }
61 
62 
63 
64 
65 //
66 // member functions
67 //
68 
69 // ------------ method called to produce the data ------------
70 void
72 {
73  using namespace edm;
75  std::vector<L1MuKBMTrack> outAll;
76  iEvent.getByToken(src_,stubHandle);
77 
79  for (uint i=0;i<stubHandle->size();++i) {
80  L1MuKBMTCombinedStubRef r(stubHandle,i);
81  stubs.push_back(r);
82  }
83 
84 
85  std::unique_ptr<l1t::RegionalMuonCandBxCollection> outBMTF(new l1t::RegionalMuonCandBxCollection());
86  outBMTF->setBXRange(-3,3);
88  for (const auto& bx : bx_) {
89  L1MuKBMTrackCollection tmp = trackFinder_->process(algo_.get(),stubs,bx);
90  for (const auto& track :tmp) {
91  algo_->addBMTFMuon(bx,track,outBMTF);
92  }
93  if (!tmp.empty())
94  out.insert(out.end(),tmp.begin(),tmp.end());
95  }
96 
97 
98 
99 
100 
101  iEvent.put(std::move(outBMTF),"BMTF");
102  std::unique_ptr<L1MuKBMTrackCollection > out1 = std::make_unique<L1MuKBMTrackCollection >(out);
103  iEvent.put(std::move(out1));
104 
105 }
106 
107 // ------------ method called once each stream before processing any runs, lumis or events ------------
108 void
110 {
111 }
112 
113 // ------------ method called once each stream after processing all runs, lumis and events ------------
114 void
116 }
117 
118 void
120  //The following says we do not know what parameters are allowed so do no validation
121  // Please change this to state exactly what you do use, even if it is no parameters
123  desc.setUnknown();
124  descriptions.addDefault(desc);
125 }
126 
127 //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:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void produce(edm::Event &, const edm::EventSetup &) override
L1TMuonBarrelKalmanTrackProducer(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< L1MuKBMTCombinedStub > > src_
std::vector< edm::Ref< L1MuKBMTCombinedStubCollection > > L1MuKBMTCombinedStubRefVector
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
std::unique_ptr< L1TMuonBarrelKalmanTrackFinder > trackFinder_
std::unique_ptr< L1TMuonBarrelKalmanAlgo > algo_
def uint(string)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
HLT enums.
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:14
def move(src, dest)
Definition: eostools.py:510