CMS 3D CMS Logo

PixelCPETemplateRecoESProducer.cc
Go to the documentation of this file.
10 
15 
16 #include <string>
17 #include <memory>
18 
20 public:
22  std::unique_ptr<PixelClusterParameterEstimator> produce(const TkPixelCPERecord&);
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25 private:
32 
36 };
37 
38 using namespace edm;
39 
41  std::string myname = p.getParameter<std::string>("ComponentName");
42 
43  useLAFromDB_ = p.getParameter<bool>("useLAFromDB");
44  doLorentzFromAlignment_ = p.getParameter<bool>("doLorentzFromAlignment");
45 
46  pset_ = p;
47  auto c = setWhatProduced(this, myname);
48  magfieldToken_ = c.consumes();
49  pDDToken_ = c.consumes();
50  hTTToken_ = c.consumes();
51  templateDBobjectToken_ = c.consumes();
52  templateStoreToken_ = c.consumes();
53  if (useLAFromDB_ || doLorentzFromAlignment_) {
54  char const* laLabel = doLorentzFromAlignment_ ? "fromAlignment" : "";
55  lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));
56  }
57 }
58 
59 std::unique_ptr<PixelClusterParameterEstimator> PixelCPETemplateRecoESProducer::produce(
60  const TkPixelCPERecord& iRecord) {
61  // Normal, default LA is used in case of template failure, load it unless
62  // turned off
63  // if turned off, null is ok, becomes zero
64  const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
65  if (useLAFromDB_ || doLorentzFromAlignment_) {
66  lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
67  }
68 
69  return std::make_unique<PixelCPETemplateReco>(pset_,
70  &iRecord.get(magfieldToken_),
71  iRecord.get(pDDToken_),
72  iRecord.get(hTTToken_),
73  lorentzAngleProduct,
74  &iRecord.get(templateStoreToken_),
75  &iRecord.get(templateDBobjectToken_));
76 }
77 
80 
81  // from PixelCPEBase
83 
84  // from PixelCPETemplateReco
86 
87  // specific to PixelCPETemplateRecoESProducer
88  desc.add<std::string>("ComponentName", "PixelCPETemplateReco");
89  descriptions.add("_templates_default", desc);
90 }
91 
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > templateDBobjectToken_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
std::unique_ptr< PixelClusterParameterEstimator > produce(const TkPixelCPERecord &)
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
edm::ESGetToken< std::vector< SiPixelTemplateStore >, SiPixelTemplateDBObjectESProducerRcd > templateStoreToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
PixelCPETemplateRecoESProducer(const edm::ParameterSet &p)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > hTTToken_
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleToken_