CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PixelCPEGenericESProducer.cc
Go to the documentation of this file.
11 
16 
17 // new record
19 
20 #include <string>
21 #include <memory>
22 
24 public:
26  std::unique_ptr<PixelClusterParameterEstimator> produce(const TkPixelCPERecord&);
27  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
28 
29 private:
36 
40  std::string CPEgenericMode_; // user's choice of CPE generic
41 };
42 
43 using namespace edm;
44 
46  CPEgenericMode_ = p.getParameter<std::string>("ComponentName");
47  // Use LA-width from DB. If both (upper and this) are false LA-width is calcuated from LA-offset
48  useLAWidthFromDB_ = p.getParameter<bool>("useLAWidthFromDB");
49  // Use Alignment LA-offset
50  const bool doLorentzFromAlignment = p.getParameter<bool>("doLorentzFromAlignment");
51  char const* laLabel = ""; // standard LA, from calibration, label=""
52  if (doLorentzFromAlignment) {
53  laLabel = "fromAlignment";
54  }
55 
56  auto magname = p.getParameter<edm::ESInputTag>("MagneticFieldRecord");
57  UseErrorsFromTemplates_ = p.getParameter<bool>("UseErrorsFromTemplates");
58 
59  pset_ = p;
60  auto c = setWhatProduced(this, CPEgenericMode_);
61  magfieldToken_ = c.consumes(magname);
62  pDDToken_ = c.consumes();
63  hTTToken_ = c.consumes();
64  lorentzAngleToken_ = c.consumes(edm::ESInputTag("", laLabel));
65  if (useLAWidthFromDB_) {
66  lorentzAngleWidthToken_ = c.consumes(edm::ESInputTag("", "forWidth"));
67  }
68  if (UseErrorsFromTemplates_) {
69  genErrorDBObjectToken_ = c.consumes();
70  }
71 }
72 
73 std::unique_ptr<PixelClusterParameterEstimator> PixelCPEGenericESProducer::produce(const TkPixelCPERecord& iRecord) {
74  // add the new la width object
75  const SiPixelLorentzAngle* lorentzAngleWidthProduct = nullptr;
76  if (useLAWidthFromDB_) { // use the width LA
77  lorentzAngleWidthProduct = &iRecord.get(lorentzAngleWidthToken_);
78  }
79  //std::cout<<" la width "<<lorentzAngleWidthProduct<<std::endl; //dk
80 
81  const SiPixelGenErrorDBObject* genErrorDBObjectProduct = nullptr;
82 
83  // Errors take only from new GenError
84  if (UseErrorsFromTemplates_) { // do only when generrors are needed
85  genErrorDBObjectProduct = &iRecord.get(genErrorDBObjectToken_);
86  //} else {
87  //std::cout<<" pass an empty GenError pointer"<<std::endl;
88  }
89 
90  return (CPEgenericMode_ == "PixelCPEGenericForBricked")
91  ? std::make_unique<PixelCPEGenericForBricked>(pset_,
92  &iRecord.get(magfieldToken_),
93  iRecord.get(pDDToken_),
94  iRecord.get(hTTToken_),
95  &iRecord.get(lorentzAngleToken_),
96  genErrorDBObjectProduct,
97  lorentzAngleWidthProduct)
98  : std::make_unique<PixelCPEGeneric>(pset_,
99  &iRecord.get(magfieldToken_),
100  iRecord.get(pDDToken_),
101  iRecord.get(hTTToken_),
102  &iRecord.get(lorentzAngleToken_),
103  genErrorDBObjectProduct,
104  lorentzAngleWidthProduct);
105 }
106 
109 
110  // from PixelCPEBase
112 
113  // from PixelCPEGeneric
115 
116  // specific to PixelCPEGenericESProducer
117  desc.add<std::string>("ComponentName", "PixelCPEGeneric");
118  desc.add<edm::ESInputTag>("MagneticFieldRecord", edm::ESInputTag(""));
119  descriptions.add("_generic_default", desc);
120 }
121 
const edm::EventSetup & c
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleToken_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > hTTToken_
tuple doLorentzFromAlignment
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > lorentzAngleWidthToken_
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::unique_ptr< PixelClusterParameterEstimator > produce(const TkPixelCPERecord &)
edm::ESGetToken< SiPixelGenErrorDBObject, SiPixelGenErrorDBObjectRcd > genErrorDBObjectToken_
static void fillPSetDescription(edm::ParameterSetDescription &desc)
PixelCPEGenericESProducer(const edm::ParameterSet &p)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)