CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Types | Private Member Functions | Private Attributes
SiPixelFakeLorentzAngleESSource Class Reference

#include <CalibTracker/SiPixelESProducer/src/SiPixelFakeLorentzAngleESSource.cc>

Inheritance diagram for SiPixelFakeLorentzAngleESSource:
edm::ESProducer edm::EventSetupRecordIntervalFinder edm::ESProxyFactoryProducer edm::eventsetup::DataProxyProvider

Public Member Functions

virtual std::unique_ptr< SiPixelLorentzAngleproduce (const SiPixelLorentzAngleRcd &)
 
 SiPixelFakeLorentzAngleESSource (const edm::ParameterSet &)
 
 ~SiPixelFakeLorentzAngleESSource () override=default
 
- Public Member Functions inherited from edm::ESProducer
 ESProducer ()
 
 ESProducer (const ESProducer &)=delete
 
ESProxyIndex const * getTokenIndices (unsigned int iIndex) const
 
ESRecordIndex const * getTokenRecordIndices (unsigned int iIndex) const
 
bool hasMayConsumes () const noexcept
 
size_t numberOfTokenIndices (unsigned int iIndex) const
 
ESProducer const & operator= (const ESProducer &)=delete
 
SerialTaskQueueChainqueue ()
 
template<typename Record >
std::optional< std::vector< ESProxyIndex > > updateFromMayConsumes (unsigned int iIndex, const Record &iRecord) const
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &) final
 
 ~ESProducer () noexcept(false) override
 
- Public Member Functions inherited from edm::ESProxyFactoryProducer
 ESProxyFactoryProducer ()
 
 ESProxyFactoryProducer (const ESProxyFactoryProducer &)=delete
 
const ESProxyFactoryProduceroperator= (const ESProxyFactoryProducer &)=delete
 
 ~ESProxyFactoryProducer () noexcept(false) override
 
- Public Member Functions inherited from edm::eventsetup::DataProxyProvider
void createKeyedProxies (EventSetupRecordKey const &key, unsigned int nConcurrentIOVs)
 
 DataProxyProvider ()
 
 DataProxyProvider (const DataProxyProvider &)=delete
 
const ComponentDescriptiondescription () const
 
void fillRecordsNotAllowingConcurrentIOVs (std::set< EventSetupRecordKey > &recordsNotAllowingConcurrentIOVs) const
 
virtual void initConcurrentIOVs (EventSetupRecordKey const &key, unsigned int nConcurrentIOVs)
 
bool isUsingRecord (const EventSetupRecordKey &key) const
 
KeyedProxieskeyedProxies (const EventSetupRecordKey &iRecordKey, unsigned int iovIndex=0)
 
const DataProxyProvideroperator= (const DataProxyProvider &)=delete
 
void setAppendToDataLabel (const edm::ParameterSet &)
 
void setDescription (const ComponentDescription &iDescription)
 
std::set< EventSetupRecordKeyusingRecords () const
 
virtual ~DataProxyProvider () noexcept(false)
 
- Public Member Functions inherited from edm::EventSetupRecordIntervalFinder
bool concurrentFinder () const
 
const eventsetup::ComponentDescriptiondescriptionForFinder () const
 
 EventSetupRecordIntervalFinder ()
 
 EventSetupRecordIntervalFinder (const EventSetupRecordIntervalFinder &)=delete
 
std::set< eventsetup::EventSetupRecordKeyfindingForRecords () const
 
const ValidityIntervalfindIntervalFor (const eventsetup::EventSetupRecordKey &, const IOVSyncValue &)
 
bool nonconcurrentAndIOVNeedsUpdate (const eventsetup::EventSetupRecordKey &key, const IOVSyncValue &syncValue) const
 
const EventSetupRecordIntervalFinderoperator= (const EventSetupRecordIntervalFinder &)=delete
 
void resetInterval (const eventsetup::EventSetupRecordKey &)
 
void setDescriptionForFinder (const eventsetup::ComponentDescription &iDescription)
 
virtual ~EventSetupRecordIntervalFinder () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::eventsetup::DataProxyProvider
static void prevalidate (ConfigurationDescriptions &)
 

Protected Member Functions

void setIntervalFor (const edm::eventsetup::EventSetupRecordKey &, const edm::IOVSyncValue &, edm::ValidityInterval &) override
 
- Protected Member Functions inherited from edm::ESProducer
template<typename T >
auto setWhatProduced (T *iThis, const es::Label &iLabel={})
 
template<typename T >
auto setWhatProduced (T *iThis, const char *iLabel)
 
template<typename T >
auto setWhatProduced (T *iThis, const std::string &iLabel)
 
template<typename T , typename TDecorator >
auto setWhatProduced (T *iThis, const TDecorator &iDec, const es::Label &iLabel={})
 
template<typename T , typename TReturn , typename TRecord >
auto setWhatProduced (T *iThis, TReturn(T ::*iMethod)(const TRecord &), const es::Label &iLabel={})
 
template<typename T , typename TReturn , typename TRecord , typename TArg >
auto setWhatProduced (T *iThis, TReturn(T ::*iMethod)(const TRecord &), const TArg &iDec, const es::Label &iLabel={})
 
template<typename TFunc >
auto setWhatProduced (TFunc &&func, const es::Label &iLabel={})
 
template<typename TReturn , typename TRecord , typename TFunc , typename TDecorator >
ESConsumesCollectorT< TRecord > setWhatProduced (TFunc &&func, TDecorator &&iDec, const es::Label &iLabel={})
 
void usesResources (std::vector< std::string > const &)
 
- Protected Member Functions inherited from edm::ESProxyFactoryProducer
template<class TFactory >
void registerFactory (std::unique_ptr< TFactory > iFactory, const std::string &iLabel=std::string())
 
virtual void registerFactoryWithKey (const EventSetupRecordKey &iRecord, std::unique_ptr< eventsetup::ProxyFactoryBase > iFactory, const std::string &iLabel=std::string())
 
KeyedProxiesVector registerProxies (const EventSetupRecordKey &, unsigned int iovIndex) override
 
- Protected Member Functions inherited from edm::eventsetup::DataProxyProvider
template<class T >
void usingRecord ()
 
void usingRecordWithKey (const EventSetupRecordKey &key)
 
- Protected Member Functions inherited from edm::EventSetupRecordIntervalFinder
template<class T >
void findingRecord ()
 
void findingRecordWithKey (const eventsetup::EventSetupRecordKey &)
 

Private Types

typedef std::vector< edm::ParameterSetParameters
 

Private Member Functions

int HVgroup (int panel, int module)
 

Private Attributes

float bPixLorentzAnglePerTesla_
 
Parameters BPixParameters_
 
const edm::FileInPath fp_
 
float fPixLorentzAnglePerTesla_
 
Parameters FPixParameters_
 
Parameters ModuleParameters_
 
const std::string myLabel_
 
const edm::FileInPath t_topo_fp_
 

Additional Inherited Members

- Protected Types inherited from edm::ESProxyFactoryProducer
using EventSetupRecordKey = eventsetup::EventSetupRecordKey
 
- Protected Types inherited from edm::eventsetup::DataProxyProvider
using KeyedProxiesVector = std::vector< std::pair< DataKey, std::shared_ptr< DataProxy > >>
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 35 of file SiPixelFakeLorentzAngleESSource.h.

Member Typedef Documentation

◆ Parameters

Definition at line 55 of file SiPixelFakeLorentzAngleESSource.h.

Constructor & Destructor Documentation

◆ SiPixelFakeLorentzAngleESSource()

SiPixelFakeLorentzAngleESSource::SiPixelFakeLorentzAngleESSource ( const edm::ParameterSet conf_)

Definition at line 35 of file SiPixelFakeLorentzAngleESSource.cc.

References edm::ESProducer::setWhatProduced().

36  : fp_(conf_.getParameter<edm::FileInPath>("file")),
37  t_topo_fp_(conf_.getParameter<edm::FileInPath>("topologyInput")),
38  myLabel_(conf_.getParameter<std::string>("appendToDataLabel")),
39  BPixParameters_(conf_.getParameter<Parameters>("BPixParameters")),
40  FPixParameters_(conf_.getParameter<Parameters>("FPixParameters")),
41  ModuleParameters_(conf_.getParameter<Parameters>("ModuleParameters")),
42  bPixLorentzAnglePerTesla_((float)conf_.getUntrackedParameter<double>("bPixLorentzAnglePerTesla", -9999.)),
43  fPixLorentzAnglePerTesla_((float)conf_.getUntrackedParameter<double>("fPixLorentzAnglePerTesla", -9999.)) {
44  edm::LogInfo("SiPixelFakeLorentzAngleESSource::SiPixelFakeLorentzAngleESSource");
45  // the following line is needed to tell the framework what data is being produced
46  setWhatProduced(this);
47  findingRecord<SiPixelLorentzAngleRcd>();
48 }
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T getUntrackedParameter(std::string const &, T const &) const
vector< ParameterSet > Parameters
Log< level::Info, false > LogInfo

◆ ~SiPixelFakeLorentzAngleESSource()

SiPixelFakeLorentzAngleESSource::~SiPixelFakeLorentzAngleESSource ( )
overridedefault

Member Function Documentation

◆ fillDescriptions()

void SiPixelFakeLorentzAngleESSource::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 239 of file SiPixelFakeLorentzAngleESSource.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addOptional(), edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and AlCaHLTBitMon_QueryRunRegistry::string.

239  {
241  desc.setComment("ESSource to supply per-module SiPixelLorentzAngle payloads in the EventSetup");
242 
243  desc.add<edm::FileInPath>(
244  "file", edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt"))
245  ->setComment("Tracker skimmed geometry");
246  desc.add<edm::FileInPath>("topologyInput",
247  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml"))
248  ->setComment("Tracker Topology");
249 
250  desc.add<std::string>("appendToDataLabel", "")->setComment("label to which write the data");
251  desc.addUntracked<double>("bPixLorentzAnglePerTesla", -9999.)->setComment("LA value for all BPix");
252  desc.addUntracked<double>("fPixLorentzAnglePerTesla", -9999.)->setComment("LA value for all FPix");
253 
254  edm::ParameterSetDescription desc_BPixParameters;
255  desc_BPixParameters.addOptional<int>("layer");
256  desc_BPixParameters.addOptional<int>("ladder");
257  desc_BPixParameters.addOptional<int>("module");
258  desc_BPixParameters.addOptional<int>("side");
259  desc_BPixParameters.add<double>("angle");
260  std::vector<edm::ParameterSet> default_BPixParametersToAdd;
261  desc.addVPSet("BPixParameters", desc_BPixParameters, default_BPixParametersToAdd)
262  ->setComment("LA values for given BPix regions");
263 
264  edm::ParameterSetDescription desc_FPixParameters;
265  desc_FPixParameters.addOptional<int>("side");
266  desc_FPixParameters.addOptional<int>("disk");
267  desc_FPixParameters.addOptional<int>("ring");
268  desc_FPixParameters.addOptional<int>("blade");
269  desc_FPixParameters.addOptional<int>("panel");
270  desc_FPixParameters.addOptional<int>("HVgroup");
271  desc_FPixParameters.add<double>("angle");
272  std::vector<edm::ParameterSet> default_FPixParametersToAdd;
273  desc.addVPSet("FPixParameters", desc_FPixParameters, default_FPixParametersToAdd)
274  ->setComment("LA values for given FPix regions");
275 
276  edm::ParameterSetDescription desc_ModuleParameters;
277  desc_ModuleParameters.add<unsigned int>("rawid");
278  desc_ModuleParameters.add<double>("angle");
279  std::vector<edm::ParameterSet> default_ModuleParametersToAdd;
280  desc.addVPSet("ModuleParameters", desc_ModuleParameters, default_ModuleParametersToAdd)
281  ->setComment("LA values for given modules");
282 
283  descriptions.addWithDefaultLabel(desc);
284 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ HVgroup()

int SiPixelFakeLorentzAngleESSource::HVgroup ( int  panel,
int  module 
)
private

Definition at line 222 of file SiPixelFakeLorentzAngleESSource.cc.

References callgraph::module.

Referenced by produce().

222  {
223  if (1 == panel && (1 == module || 2 == module)) {
224  return 1;
225  } else if (1 == panel && (3 == module || 4 == module)) {
226  return 2;
227  } else if (2 == panel && 1 == module) {
228  return 1;
229  } else if (2 == panel && (2 == module || 3 == module)) {
230  return 2;
231  } else {
232  edm::LogError("SiPixelFakeLorentzAngleESSource")
233  << " ERROR! in SiPixelFakeLorentzAngleESSource::HVgroup(...), panel = " << panel << ", module = " << module
234  << std::endl;
235  return 0;
236  }
237 }
Log< level::Error, false > LogError

◆ produce()

std::unique_ptr< SiPixelLorentzAngle > SiPixelFakeLorentzAngleESSource::produce ( const SiPixelLorentzAngleRcd )
virtual

Definition at line 50 of file SiPixelFakeLorentzAngleESSource.cc.

References PixelEndcapName::bladeName(), bPixLorentzAnglePerTesla_, BPixParameters_, PixelEndcapName::diskName(), dqmMemoryStats::float, newFWLiteAna::found, fp_, fPixLorentzAnglePerTesla_, FPixParameters_, StandaloneTrackerTopology::fromTrackerParametersXMLFile(), edm::FileInPath::fullPath(), HVgroup(), PVValHelper::ladder, phase1PixelTopology::layer, LogDebug, callgraph::module, ModuleParameters_, getGTfromDQMFile::obj, PixelEndcapName::pannelName(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, TrackerTopology::pxbLadder(), TrackerTopology::pxbLayer(), TrackerTopology::pxbModule(), DetId::rawId(), DQM::reader, relativeConstraints::ring, PixelEndcapName::ringName(), TrackerTopology::side(), DetId::subdetId(), and t_topo_fp_.

50  {
51  using namespace edm::es;
52  unsigned int nmodules = 0;
55  const std::vector<uint32_t>& DetIds = reader.getAllDetIds();
56 
58 
59  // Loop over detectors
60  for (const auto& detit : DetIds) {
61  nmodules++;
62  const DetId detid(detit);
63  auto rawId = detid.rawId();
64  int found = 0;
65  int side = tTopo.side(detid); // 1:-z 2:+z for fpix, for bpix gives 0
66 
67  // fill bpix values for LA
68  if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
69  int layer = tTopo.pxbLayer(detid);
70  // Barrel ladder id 1-20,32,44.
71  int ladder = tTopo.pxbLadder(detid);
72  // Barrel Z-index=1,8
73  int module = tTopo.pxbModule(detid);
74  if (module < 5) {
75  side = 1;
76  } else {
77  side = 2;
78  }
79 
80  LogDebug("SiPixelFakeLorentzAngleESSource") << " pixel barrel:"
81  << " layer=" << layer << " ladder=" << ladder << " module=" << module
82  << " rawId=" << rawId << " side " << side;
83 
84  // use a commmon value (e.g. for MC)
85  if (bPixLorentzAnglePerTesla_ != -9999.) { // use common value for all
86  edm::LogInfo("SiPixelFakeLorentzAngleESSource")
87  << " LA = " << bPixLorentzAnglePerTesla_ << " common for all bpix" << std::endl;
88  if (!obj->putLorentzAngle(detid.rawId(), bPixLorentzAnglePerTesla_))
89  edm::LogError("SiPixelFakeLorentzAngleESSource")
90  << "ERROR!: detid " << rawId << " already exists" << std::endl;
91  } else {
92  //first individuals are put
93  for (const auto& it : ModuleParameters_) {
94  if (it.getParameter<unsigned int>("rawid") == detid.rawId()) {
95  float lorentzangle = (float)it.getParameter<double>("angle");
96  if (!found) {
97  obj->putLorentzAngle(detid.rawId(), lorentzangle);
98  edm::LogInfo("SiPixelFakeLorentzAngleESSource")
99  << " LA= " << lorentzangle << " individual value " << detid.rawId() << std::endl;
100  found = 1;
101  } else
102  edm::LogError("SiPixelFakeLorentzAngleESSource") << "ERROR!: detid already exists" << std::endl;
103  }
104  }
105 
106  //modules already put are automatically skipped
107  for (const auto& it : BPixParameters_) {
108  if (it.exists("layer"))
109  if (it.getParameter<int>("layer") != layer)
110  continue;
111  if (it.exists("ladder"))
112  if (it.getParameter<int>("ladder") != ladder)
113  continue;
114  if (it.exists("module"))
115  if (it.getParameter<int>("module") != module)
116  continue;
117  if (it.exists("side"))
118  if (it.getParameter<int>("side") != side)
119  continue;
120  if (!found) {
121  float lorentzangle = (float)it.getParameter<double>("angle");
122  obj->putLorentzAngle(detid.rawId(), lorentzangle);
123  edm::LogInfo("SiPixelFakeLorentzAngleESSource") << " LA= " << lorentzangle << std::endl;
124  found = 2;
125  } else if (found == 1) {
126  edm::LogWarning("SiPixelFakeLorentzAngleESSource")
127  << "detid already given in ModuleParameters, skipping ..." << std::endl;
128  } else
129  edm::LogError("SiPixelFakeLorentzAngleESSource")
130  << " ERROR!: detid " << rawId << " already exists" << std::endl;
131  }
132  }
133  } else if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
134  // fill fpix values for LA (for phase2 fpix & epix)
135 
136  // For fpix we also need to find the ring number which is not available from topology
137  // Convert to online
138  PixelEndcapName pen(detid, &tTopo, true); // use det-id phaseq
139 
140  //PixelEndcapName::HalfCylinder sh = pen.halfCylinder(); //enum
141  //string nameF = pen.name();
142  //int plaquetteName = pen.plaquetteName();
143  int disk = pen.diskName();
144  int blade = pen.bladeName();
145  int panel = pen.pannelName();
146  int ring = pen.ringName();
147 
148  LogDebug("SiPixelFakeLorentzAngleESSource") << " pixel endcap:"
149  << " side=" << side << " disk=" << disk << " blade =" << blade
150  << " pannel=" << panel << " ring=" << ring << " rawId=" << rawId;
151 
152  // use a commmon value (e.g. for MC)
153  if (fPixLorentzAnglePerTesla_ != -9999.) { // use common value for all
154  edm::LogInfo("SiPixelFakeLorentzAngleESSource")
155  << " LA = " << fPixLorentzAnglePerTesla_ << " common for all fpix" << std::endl;
156  if (!obj->putLorentzAngle(detid.rawId(), fPixLorentzAnglePerTesla_))
157  edm::LogError("SiPixelFakeLorentzAngleESSource") << " ERROR! detid already exists" << std::endl;
158  } else {
159  //first individuals are put
160  for (const auto& it : ModuleParameters_) {
161  if (it.getParameter<unsigned int>("rawid") == detid.rawId()) {
162  float lorentzangle = (float)it.getParameter<double>("angle");
163  if (!found) {
164  obj->putLorentzAngle(detid.rawId(), lorentzangle);
165  edm::LogInfo("SiPixelFakeLorentzAngleESSource")
166  << " LA= " << lorentzangle << " individual value " << detid.rawId() << std::endl;
167  found = 1;
168  } else
169  edm::LogError("SiPixelFakeLorentzAngleESSource")
170  << "ERROR!: detid " << rawId << " already exists" << std::endl;
171  } // if
172  } // for
173 
174  //modules already put are automatically skipped
175  for (const auto& it : FPixParameters_) {
176  if (it.exists("side"))
177  if (it.getParameter<int>("side") != side)
178  continue;
179  if (it.exists("disk"))
180  if (it.getParameter<int>("disk") != disk)
181  continue;
182  if (it.exists("ring"))
183  if (it.getParameter<int>("ring") != ring)
184  continue;
185  if (it.exists("blade"))
186  if (it.getParameter<int>("blade") != blade)
187  continue;
188  if (it.exists("panel"))
189  if (it.getParameter<int>("panel") != panel)
190  continue;
191  if (it.exists("HVgroup"))
192  if (it.getParameter<int>("HVgroup") != HVgroup(panel, ring))
193  continue;
194  if (!found) {
195  float lorentzangle = (float)it.getParameter<double>("angle");
196  obj->putLorentzAngle(detid.rawId(), lorentzangle);
197  edm::LogInfo("SiPixelFakeLorentzAngleESSource") << " LA= " << lorentzangle << std::endl;
198  found = 2;
199  } else if (found == 1) {
200  edm::LogWarning("SiPixelFakeLorentzAngleESSource")
201  << " detid " << rawId << " already given in ModuleParameters, skipping ..." << std::endl;
202  } else
203  edm::LogError("SiPixelFakeLorentzAngleESSource")
204  << "ERROR!: detid" << rawId << "already exists" << std::endl;
205  } // for
206  } // if
207  } // bpix/fpix
208  } // iterate on detids
209 
210  edm::LogInfo("SiPixelFakeLorentzAngleESSource") << "Modules = " << nmodules << std::endl;
211 
212  return std::unique_ptr<SiPixelLorentzAngle>(obj);
213 }
unsigned int pxbLayer(const DetId &id) const
std::string fullPath() const
Definition: FileInPath.cc:161
reader
Definition: DQM.py:105
unsigned int pxbLadder(const DetId &id) const
Log< level::Error, false > LogError
unsigned int side(const DetId &id) const
constexpr std::array< uint8_t, layerIndexSize > layer
Log< level::Info, false > LogInfo
Definition: DetId.h:17
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
unsigned int pxbModule(const DetId &id) const
Log< level::Warning, false > LogWarning
#define LogDebug(id)

◆ setIntervalFor()

void SiPixelFakeLorentzAngleESSource::setIntervalFor ( const edm::eventsetup::EventSetupRecordKey ,
const edm::IOVSyncValue iosv,
edm::ValidityInterval oValidity 
)
overrideprotectedvirtual

Implements edm::EventSetupRecordIntervalFinder.

Definition at line 215 of file SiPixelFakeLorentzAngleESSource.cc.

References edm::IOVSyncValue::beginOfTime(), edm::IOVSyncValue::endOfTime(), and infinity.

217  {
219  oValidity = infinity;
220 }
static const IOVSyncValue & endOfTime()
Definition: IOVSyncValue.cc:82
static const IOVSyncValue & beginOfTime()
Definition: IOVSyncValue.cc:88
const double infinity

Member Data Documentation

◆ bPixLorentzAnglePerTesla_

float SiPixelFakeLorentzAngleESSource::bPixLorentzAnglePerTesla_
private

Definition at line 60 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().

◆ BPixParameters_

Parameters SiPixelFakeLorentzAngleESSource::BPixParameters_
private

Definition at line 56 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().

◆ fp_

const edm::FileInPath SiPixelFakeLorentzAngleESSource::fp_
private

Definition at line 52 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().

◆ fPixLorentzAnglePerTesla_

float SiPixelFakeLorentzAngleESSource::fPixLorentzAnglePerTesla_
private

Definition at line 61 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().

◆ FPixParameters_

Parameters SiPixelFakeLorentzAngleESSource::FPixParameters_
private

Definition at line 57 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().

◆ ModuleParameters_

Parameters SiPixelFakeLorentzAngleESSource::ModuleParameters_
private

Definition at line 58 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().

◆ myLabel_

const std::string SiPixelFakeLorentzAngleESSource::myLabel_
private

Definition at line 54 of file SiPixelFakeLorentzAngleESSource.h.

◆ t_topo_fp_

const edm::FileInPath SiPixelFakeLorentzAngleESSource::t_topo_fp_
private

Definition at line 53 of file SiPixelFakeLorentzAngleESSource.h.

Referenced by produce().