CMS 3D CMS Logo

PixelCPEFastESProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 
18 
20 public:
22  std::unique_ptr<PixelClusterParameterEstimator> produce(const TkPixelCPERecord&);
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25 private:
32 
35 };
36 
37 using namespace edm;
38 
40  auto const& myname = p.getParameter<std::string>("ComponentName");
41  auto const& magname = p.getParameter<edm::ESInputTag>("MagneticFieldRecord");
42  useErrorsFromTemplates_ = p.getParameter<bool>("UseErrorsFromTemplates");
43 
44  auto cc = setWhatProduced(this, myname);
45  magfieldToken_ = cc.consumes(magname);
46  pDDToken_ = cc.consumes();
47  hTTToken_ = cc.consumes();
48  lorentzAngleToken_ = cc.consumes(edm::ESInputTag(""));
49  lorentzAngleWidthToken_ = cc.consumes(edm::ESInputTag("", "forWidth"));
51  genErrorDBObjectToken_ = cc.consumes();
52  }
53 }
54 
55 std::unique_ptr<PixelClusterParameterEstimator> PixelCPEFastESProducer::produce(const TkPixelCPERecord& iRecord) {
56  // add the new la width object
57  const SiPixelLorentzAngle* lorentzAngleWidthProduct = nullptr;
58  lorentzAngleWidthProduct = &iRecord.get(lorentzAngleWidthToken_);
59 
60  const SiPixelGenErrorDBObject* genErrorDBObjectProduct = nullptr;
61 
62  // Errors take only from new GenError
63  if (useErrorsFromTemplates_) { // do only when generrors are needed
64  genErrorDBObjectProduct = &iRecord.get(genErrorDBObjectToken_);
65  //} else {
66  //std::cout<<" pass an empty GenError pointer"<<std::endl;
67  }
68  return std::make_unique<PixelCPEFast>(pset_,
69  &iRecord.get(magfieldToken_),
70  iRecord.get(pDDToken_),
71  iRecord.get(hTTToken_),
72  &iRecord.get(lorentzAngleToken_),
73  genErrorDBObjectProduct,
74  lorentzAngleWidthProduct);
75 }
76 
79 
80  // from PixelCPEBase
82 
83  // from PixelCPEFast
85 
86  // used by PixelCPEFast
87  desc.add<double>("EdgeClusterErrorX", 50.0);
88  desc.add<double>("EdgeClusterErrorY", 85.0);
89  desc.add<bool>("UseErrorsFromTemplates", true);
90  desc.add<bool>("TruncatePixelCharge", true);
91 
92  // specific to PixelCPEFastESProducer
93  desc.add<std::string>("ComponentName", "PixelCPEFast");
94  desc.add<edm::ESInputTag>("MagneticFieldRecord", edm::ESInputTag());
95 
96  descriptions.add("PixelCPEFastESProducer", desc);
97 }
98 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleWidthToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< SiPixelGenErrorDBObject, SiPixelGenErrorDBObjectRcd > genErrorDBObjectToken_
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
PixelCPEFastESProducer(const edm::ParameterSet &p)
HLT enums.
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
std::unique_ptr< PixelClusterParameterEstimator > produce(const TkPixelCPERecord &)
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > hTTToken_