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:
31 
35 };
36 
37 using namespace edm;
38 
40  std::string myname = p.getParameter<std::string>("ComponentName");
41 
42  useLAFromDB_ = p.getParameter<bool>("useLAFromDB");
43  doLorentzFromAlignment_ = p.getParameter<bool>("doLorentzFromAlignment");
44 
45  pset_ = p;
46  auto c = setWhatProduced(this, myname);
47  magfieldToken_ = c.consumes();
48  pDDToken_ = c.consumes();
49  hTTToken_ = c.consumes();
50  templateDBobjectToken_ = c.consumes();
51  if (useLAFromDB_ || doLorentzFromAlignment_) {
52  char const* laLabel = doLorentzFromAlignment_ ? "fromAlignment" : "";
53  lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));
54  }
55 }
56 
57 std::unique_ptr<PixelClusterParameterEstimator> PixelCPETemplateRecoESProducer::produce(
58  const TkPixelCPERecord& iRecord) {
59  // Normal, default LA is used in case of template failure, load it unless
60  // turned off
61  // if turned off, null is ok, becomes zero
62  const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
63  if (useLAFromDB_ || doLorentzFromAlignment_) {
64  lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
65  }
66 
67  return std::make_unique<PixelCPETemplateReco>(pset_,
68  &iRecord.get(magfieldToken_),
69  iRecord.get(pDDToken_),
70  iRecord.get(hTTToken_),
71  lorentzAngleProduct,
72  &iRecord.get(templateDBobjectToken_));
73 }
74 
77 
78  // from PixelCPEBase
80 
81  // from PixelCPETemplateReco
83 
84  // specific to PixelCPETemplateRecoESProducer
85  desc.add<std::string>("ComponentName", "PixelCPETemplateReco");
86  descriptions.add("_templates_default", desc);
87 }
88 
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > templateDBobjectToken_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
std::unique_ptr< PixelClusterParameterEstimator > produce(const TkPixelCPERecord &)
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
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_