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 
33  bool DoLorentz_;
34 };
35 
36 using namespace edm;
37 
39  std::string myname = p.getParameter<std::string>("ComponentName");
40 
41  //DoLorentz_ = p.getParameter<bool>("DoLorentz"); // True when LA from alignment is used
42  DoLorentz_ = p.getParameter<bool>("DoLorentz");
43 
44  pset_ = p;
45  auto c = setWhatProduced(this, myname);
46  c.setConsumes(magfieldToken_).setConsumes(pDDToken_).setConsumes(hTTToken_).setConsumes(templateDBobjectToken_);
47  if (DoLorentz_) {
48  c.setConsumes(lorentzAngleToken_, edm::ESInputTag("", "fromAlignment"));
49  }
50  //std::cout<<" from ES Producer Templates "<<myname<<" "<<DoLorentz_<<std::endl; //dk
51 }
52 
53 std::unique_ptr<PixelClusterParameterEstimator> PixelCPETemplateRecoESProducer::produce(
54  const TkPixelCPERecord& iRecord) {
55  // Normal, deafult LA actually is NOT needed
56  // null is ok becuse LA is not use by templates in this mode
57  const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
58  if (DoLorentz_) { // LA correction from alignment
59  lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
60  }
61 
62  return std::make_unique<PixelCPETemplateReco>(pset_,
63  &iRecord.get(magfieldToken_),
64  iRecord.get(pDDToken_),
65  iRecord.get(hTTToken_),
66  lorentzAngleProduct,
67  &iRecord.get(templateDBobjectToken_));
68 }
69 
72 
73  // from PixelCPEBase
75 
76  // from PixelCPETemplateReco
78 
79  // specific to PixelCPETemplateRecoESProducer
80  desc.add<std::string>("ComponentName", "PixelCPETemplateReco");
81  desc.add<bool>("DoLorentz", true);
82  descriptions.add("_templates_default", desc);
83 }
84 
T getParameter(std::string const &) const
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:138
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > templateDBobjectToken_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
std::unique_ptr< PixelClusterParameterEstimator > produce(const TkPixelCPERecord &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleToken_