CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Phase2StripCPEESProducer.cc
Go to the documentation of this file.
5 
11 
16 
18 
19 #include <memory>
20 #include <map>
21 
23 public:
25  std::unique_ptr<ClusterParameterEstimator<Phase2TrackerCluster1D> > produce(const TkPhase2OTCPERecord& iRecord);
26 
27 private:
28  enum CPE_t { DEFAULT, GEOMETRIC };
29 
35 };
36 
38  std::string name = p.getParameter<std::string>("ComponentType");
39 
40  std::map<std::string, CPE_t> enumMap;
41  enumMap[std::string("Phase2StripCPE")] = DEFAULT;
42  enumMap[std::string("Phase2StripCPEGeometric")] = GEOMETRIC;
43  if (enumMap.find(name) == enumMap.end())
44  throw cms::Exception("Unknown StripCPE type") << name;
45 
46  cpeNum_ = enumMap[name];
47  pset_ = p.getParameter<edm::ParameterSet>("parameters");
48  auto c = setWhatProduced(this, name);
49  if (cpeNum_ != GEOMETRIC) {
50  magfieldToken_ = c.consumes();
51  pDDToken_ = c.consumes();
52  lorentzAngleToken_ = c.consumes();
53  }
54 }
55 
56 std::unique_ptr<ClusterParameterEstimator<Phase2TrackerCluster1D> > Phase2StripCPEESProducer::produce(
57  const TkPhase2OTCPERecord& iRecord) {
58  std::unique_ptr<ClusterParameterEstimator<Phase2TrackerCluster1D> > cpe_;
59  switch (cpeNum_) {
60  case DEFAULT:
61  cpe_ = std::make_unique<Phase2StripCPE>(
62  pset_, iRecord.get(magfieldToken_), iRecord.get(pDDToken_), iRecord.get(lorentzAngleToken_));
63  break;
64 
65  case GEOMETRIC:
66  cpe_ = std::make_unique<Phase2StripCPEGeometric>(pset_);
67  break;
68  }
69  return cpe_;
70 }
71 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
const edm::EventSetup & c
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
Phase2StripCPEESProducer(const edm::ParameterSet &)
std::unique_ptr< ClusterParameterEstimator< Phase2TrackerCluster1D > > produce(const TkPhase2OTCPERecord &iRecord)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
edm::ESGetToken< SiPhase2OuterTrackerLorentzAngle, SiPhase2OuterTrackerLorentzAngleRcd > lorentzAngleToken_