CMS 3D CMS Logo

PixelCPEClusterRepairESProducer.cc
Go to the documentation of this file.
12 
17 
18 #include <string>
19 #include <memory>
20 
22 public:
25  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
26  std::unique_ptr<PixelClusterParameterEstimator> produce(const TkPixelCPERecord&);
27 
28 private:
35 
37  bool DoLorentz_;
38 };
39 
40 using namespace edm;
41 
43  std::string myname = p.getParameter<std::string>("ComponentName");
44 
45  //DoLorentz_ = p.getParameter<bool>("DoLorentz"); // True when LA from alignment is used
46  DoLorentz_ = p.getParameter<bool>("DoLorentz");
47 
48  pset_ = p;
49  auto c = setWhatProduced(this, myname);
50  c.setConsumes(magfieldToken_)
51  .setConsumes(pDDToken_)
52  .setConsumes(hTTToken_)
53  .setConsumes(templateDBobjectToken_)
54  .setConsumes(templateDBobject2DToken_);
55  if (DoLorentz_) {
56  c.setConsumes(lorentzAngleToken_, edm::ESInputTag("", "fromAlignment"));
57  }
58 
59  //std::cout<<" from ES Producer Templates "<<myname<<" "<<DoLorentz_<<std::endl; //dk
60 }
61 
63 
65  // templates2
67  desc.add<std::string>("ComponentName", "PixelCPEClusterRepair");
68 
69  // from PixelCPEBase
71 
72  // from PixelCPEClusterRepair
74 
75  // specific to PixelCPEClusterRepairESProducer
76  desc.add<bool>("DoLorentz", true);
77  descriptions.add("_templates2_default", desc);
78 }
79 
80 std::unique_ptr<PixelClusterParameterEstimator> PixelCPEClusterRepairESProducer::produce(
81  const TkPixelCPERecord& iRecord) {
82  // Normal, default LA actually is NOT needed
83  // null is ok becuse LA is not use by templates in this mode
84  const SiPixelLorentzAngle* lorentzAngleProduct = nullptr;
85  if (DoLorentz_) { // LA correction from alignment
86  lorentzAngleProduct = &iRecord.get(lorentzAngleToken_);
87  }
88 
89  return std::make_unique<PixelCPEClusterRepair>(pset_,
90  &iRecord.get(magfieldToken_),
91  iRecord.get(pDDToken_),
92  iRecord.get(hTTToken_),
93  lorentzAngleProduct,
94  &iRecord.get(templateDBobjectToken_),
95  &iRecord.get(templateDBobject2DToken_));
96 }
97 
edm::ESGetToken< SiPixelTemplateDBObject, SiPixelTemplateDBObjectESProducerRcd > templateDBobjectToken_
T getParameter(std::string const &) const
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:138
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > hTTToken_
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleToken_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
PixelCPEClusterRepairESProducer(const edm::ParameterSet &p)
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectESProducerRcd > templateDBobject2DToken_
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
std::unique_ptr< PixelClusterParameterEstimator > produce(const TkPixelCPERecord &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.