CMS 3D CMS Logo

SiPhase2OuterTrackerFakeLorentzAngleESSource.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPhase2OuterTrackerFakeLorentzAngleESSource
4 // Class: SiPhase2OuterTrackerFakeLorentzAngleESSource
5 //
11 //
12 // Original Author: Marco Musich
13 // Created: Jul 31st, 2020
14 //
15 //
16 
17 // user include files
26 //
27 // constructors and destructor
28 //
30  const edm::ParameterSet& conf_)
31  : LAvalue_(conf_.getParameter<double>("LAValue")), recordName_(conf_.getParameter<std::string>("recordName")) {
32  edm::LogInfo("SiPhase2OuterTrackerFakeLorentzAngleESSource::SiPhase2OuterTrackerFakeLorentzAngleESSource");
33  // the following line is needed to tell the framework what
34  // data is being produced
35  if (recordName_ == "LorentzAngle") {
37  m_tTopoToken = cc.consumes();
38  m_geomDetToken = cc.consumes();
39  findingRecord<SiPhase2OuterTrackerLorentzAngleRcd>();
40  } else if (recordName_ == "SimLorentzAngle") {
42  m_tTopoToken = cc.consumes();
43  m_geomDetToken = cc.consumes();
44  findingRecord<SiPhase2OuterTrackerLorentzAngleSimRcd>();
45  }
46 }
47 
49 
50 namespace fakeOTLA {
51  template <class T>
52  std::unique_ptr<T> produceRecord(const float value, const GeometricDet& geomDet) {
53  using namespace edm::es;
54  T* obj = new T();
55  for (const auto detId : TrackerGeometryUtils::getOuterTrackerDetIds(geomDet)) {
56  const DetId detectorId = DetId(detId);
57  const int subDet = detectorId.subdetId();
58  if (detectorId.det() == DetId::Detector::Tracker) {
59  if (subDet == StripSubdetector::TOB || subDet == StripSubdetector::TID) {
60  if (!obj->putLorentzAngle(detId, value))
61  edm::LogError("SiPhase2OuterTrackerFakeLorentzAngleESSource")
62  << "[SiPhase2OuterTrackerFakeLorentzAngleESSource::produce] detid already exists" << std::endl;
63  } // if it's a OT DetId
64  } // check if Tracker
65  } // loop on DetIds
66  return std::unique_ptr<T>(obj);
67  }
68 } // namespace fakeOTLA
69 
70 std::unique_ptr<SiPhase2OuterTrackerLorentzAngle> SiPhase2OuterTrackerFakeLorentzAngleESSource::produceOTLA(
72  const auto& geomDet = rcd.getRecord<TrackerTopologyRcd>().get(m_geomDetToken);
73  return fakeOTLA::produceRecord<SiPhase2OuterTrackerLorentzAngle>(LAvalue_, geomDet);
74 }
75 
76 std::unique_ptr<SiPhase2OuterTrackerLorentzAngle> SiPhase2OuterTrackerFakeLorentzAngleESSource::produceOTSimLA(
78  const auto& geomDet = rcd.getRecord<TrackerTopologyRcd>().get(m_geomDetToken);
79  return fakeOTLA::produceRecord<SiPhase2OuterTrackerLorentzAngle>(LAvalue_, geomDet);
80 }
81 
83  const edm::IOVSyncValue& iosv,
84  edm::ValidityInterval& oValidity) {
86  oValidity = infinity;
87 }
88 
91  desc.add<double>("LAValue", 0.07);
92  desc.add<std::string>("recordName", "LorentzAngle");
93  descriptions.add("siPhase2OTFakeLorentzAngleESSource", desc);
94 }
95 
96 //define this as a plug-in
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
edm::ESGetToken< GeometricDet, IdealGeometryRecord > m_geomDetToken
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
static const IOVSyncValue & endOfTime()
Definition: IOVSyncValue.cc:82
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > m_tTopoToken
static const IOVSyncValue & beginOfTime()
Definition: IOVSyncValue.cc:88
virtual std::unique_ptr< SiPhase2OuterTrackerLorentzAngle > produceOTLA(const SiPhase2OuterTrackerLorentzAngleRcd &)
const double infinity
void setIntervalFor(const edm::eventsetup::EventSetupRecordKey &, const edm::IOVSyncValue &, edm::ValidityInterval &) override
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: value.py:1
static constexpr auto TOB
Log< level::Info, false > LogInfo
#define DEFINE_FWK_EVENTSETUP_SOURCE(type)
Definition: SourceFactory.h:92
Definition: DetId.h:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< uint32_t > getOuterTrackerDetIds(const GeometricDet &geomDet)
Definition: utils.cc:20
long double T
std::unique_ptr< T > produceRecord(const float value, const GeometricDet &geomDet)
static constexpr auto TID
virtual std::unique_ptr< SiPhase2OuterTrackerLorentzAngle > produceOTSimLA(const SiPhase2OuterTrackerLorentzAngleSimRcd &)