CMS 3D CMS Logo

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