File indexing completed on 2023-03-31 23:02:21
0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/global/EDProducer.h"
0003
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011
0012 #include "RecoTracker/PixelTrackFitting/interface/PixelFitter.h"
0013 #include "RecoTracker/PixelTrackFitting/interface/PixelFitterByHelixProjections.h"
0014
0015 #include "MagneticField/Engine/interface/MagneticField.h"
0016 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0017 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0018 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0019
0020 class PixelFitterByHelixProjectionsProducer : public edm::global::EDProducer<> {
0021 public:
0022 explicit PixelFitterByHelixProjectionsProducer(const edm::ParameterSet& iConfig)
0023 : theTopoToken(esConsumes()),
0024 theFieldToken(esConsumes()),
0025 thePutToken(produces<PixelFitter>()),
0026 thescaleErrorsForBPix1(iConfig.getParameter<bool>("scaleErrorsForBPix1")),
0027 thescaleFactor(iConfig.getParameter<double>("scaleFactor")) {}
0028 ~PixelFitterByHelixProjectionsProducer() override {}
0029
0030 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0031 edm::ParameterSetDescription desc;
0032 desc.add<bool>("scaleErrorsForBPix1", false);
0033 desc.add<double>("scaleFactor", 0.65)->setComment("The default value was derived for phase1 pixel");
0034 descriptions.add("pixelFitterByHelixProjectionsDefault", desc);
0035 }
0036
0037 private:
0038 void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0039
0040 const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> theTopoToken;
0041 const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> theFieldToken;
0042 const edm::EDPutTokenT<PixelFitter> thePutToken;
0043 const bool thescaleErrorsForBPix1;
0044 const float thescaleFactor;
0045 };
0046
0047 void PixelFitterByHelixProjectionsProducer::produce(edm::StreamID,
0048 edm::Event& iEvent,
0049 const edm::EventSetup& iSetup) const {
0050 iEvent.emplace(
0051 thePutToken,
0052 std::make_unique<PixelFitterByHelixProjections>(
0053 &iSetup.getData(theTopoToken), &iSetup.getData(theFieldToken), thescaleErrorsForBPix1, thescaleFactor));
0054 }
0055
0056 DEFINE_FWK_MODULE(PixelFitterByHelixProjectionsProducer);