CMS 3D CMS Logo

L1TMuonBarrelKalmanStubProducer.cc
Go to the documentation of this file.
1 #include <memory>
4 
7 
10 
15 
17 
18 //For masks
19 
25 
26 //
27 // class declaration
28 //
29 
31 public:
34 
35  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
36 
37 private:
38  void beginStream(edm::StreamID) override;
39  void produce(edm::Event&, const edm::EventSetup&) override;
40  void endStream() override;
43  std::unique_ptr<L1TMuonBarrelKalmanStubProcessor> proc_;
44  const int verbose_;
47 };
48 
49 //
50 // constants, enums and typedefs
51 //
52 
53 //
54 // static data member definitions
55 //
56 
57 //
58 // constructors and destructor
59 //
61  : srcPhi_(consumes<L1MuDTChambPhContainer>(iConfig.getParameter<edm::InputTag>("srcPhi"))),
62  srcTheta_(consumes<L1MuDTChambThContainer>(iConfig.getParameter<edm::InputTag>("srcTheta"))),
63  proc_(std::make_unique<L1TMuonBarrelKalmanStubProcessor>(iConfig)),
64  verbose_(iConfig.getParameter<int>("verbose")),
65  bmtfParamsToken_(esConsumes()),
66  putToken_(produces<L1MuKBMTCombinedStubCollection>()) {}
67 
69 
70 //
71 // member functions
72 //
73 
74 // ------------ method called to produce the data ------------
76  using namespace edm;
78  iEvent.getByToken(srcPhi_, phiIn);
79 
81  iEvent.getByToken(srcTheta_, thetaIn);
82 
83  //Get parameters
84 
85  const L1TMuonBarrelParams& bmtfParams = iSetup.getData(bmtfParamsToken_);
86 
87  L1MuKBMTCombinedStubCollection stubs = proc_->makeStubs(phiIn.product(), thetaIn.product(), bmtfParams);
88  if (verbose_ == 1)
89  for (const auto& stub : stubs) {
90  printf("Stub: wheel=%d sector=%d station =%d tag=%d eta1=%d qeta1=%d eta2=%d qeta2=%d\n",
91  stub.whNum(),
92  stub.scNum(),
93  stub.stNum(),
94  stub.tag(),
95  stub.eta1(),
96  stub.qeta1(),
97  stub.eta2(),
98  stub.qeta2());
99  }
100 
101  if (verbose_ == 2) {
102  std::cout << "NEW" << std::endl;
103  for (uint sector = 0; sector < 12; ++sector)
104  proc_->makeInputPattern(phiIn.product(), thetaIn.product(), sector);
105  }
106 
107  iEvent.emplace(putToken_, std::move(stubs));
108 }
109 
110 // ------------ method called once each stream before processing any runs, lumis or events ------------
112 
113 // ------------ method called once each stream after processing all runs, lumis and events ------------
115 
117  //The following says we do not know what parameters are allowed so do no validation
118  // Please change this to state exactly what you do use, even if it is no parameters
120  desc.setUnknown();
121  descriptions.addDefault(desc);
122 }
123 
124 //define this as a plug-in
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T const * product() const
Definition: Handle.h:70
L1TMuonBarrelKalmanStubProducer(const edm::ParameterSet &)
const edm::EDPutTokenT< L1MuKBMTCombinedStubCollection > putToken_
void produce(edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
std::unique_ptr< L1TMuonBarrelKalmanStubProcessor > proc_
const edm::EDGetTokenT< L1MuDTChambThContainer > srcTheta_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
HLT enums.
const edm::ESGetToken< L1TMuonBarrelParams, L1TMuonBarrelParamsRcd > bmtfParamsToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511
std::vector< L1MuKBMTCombinedStub > L1MuKBMTCombinedStubCollection
const edm::EDGetTokenT< L1MuDTChambPhContainer > srcPhi_