CMS 3D CMS Logo

CTPPSOpticalFunctionsESSource.cc
Go to the documentation of this file.
1 // Original Author: Jan Kašpar
2 
9 
12 
13 //----------------------------------------------------------------------------------------------------
14 
19 public:
21  ~CTPPSOpticalFunctionsESSource() override = default;
22 
23  std::unique_ptr<LHCOpticalFunctionsSetCollection> produce(const CTPPSOpticsRcd &);
25 
26 private:
28  const edm::IOVSyncValue &,
29  edm::ValidityInterval &) override;
30 
32 
33  struct FileInfo {
34  double m_xangle;
36  };
37 
38  struct RPInfo {
41  };
42 
43  struct Entry {
45  std::vector<FileInfo> m_fileInfo;
46  std::unordered_map<unsigned int, RPInfo> m_rpInfo;
47  };
48 
49  std::vector<Entry> m_entries;
50 
52  unsigned int m_currentEntry;
53 };
54 
55 //----------------------------------------------------------------------------------------------------
56 //----------------------------------------------------------------------------------------------------
57 
59  : m_label(conf.getParameter<std::string>("label")), m_currentEntryValid(false), m_currentEntry(0) {
60  for (const auto &entry_pset : conf.getParameter<std::vector<edm::ParameterSet>>("configuration")) {
61  edm::EventRange validityRange = entry_pset.getParameter<edm::EventRange>("validityRange");
62 
63  std::vector<FileInfo> fileInfo;
64  for (const auto &pset : entry_pset.getParameter<std::vector<edm::ParameterSet>>("opticalFunctions")) {
65  const double &xangle = pset.getParameter<double>("xangle");
66  const std::string &fileName = pset.getParameter<edm::FileInPath>("fileName").fullPath();
67  fileInfo.push_back({xangle, fileName});
68  }
69 
70  std::unordered_map<unsigned int, RPInfo> rpInfo;
71  for (const auto &pset : entry_pset.getParameter<std::vector<edm::ParameterSet>>("scoringPlanes")) {
72  const unsigned int rpId = pset.getParameter<unsigned int>("rpId");
73  const std::string dirName = pset.getParameter<std::string>("dirName");
74  const double z = pset.getParameter<double>("z");
75  const RPInfo entry = {dirName, z};
76  rpInfo.emplace(rpId, entry);
77  }
78 
79  m_entries.push_back({validityRange, fileInfo, rpInfo});
80  }
81 
82  setWhatProduced(this, m_label);
83  findingRecord<CTPPSOpticsRcd>();
84 }
85 
86 //----------------------------------------------------------------------------------------------------
87 
89  const edm::IOVSyncValue &iosv,
90  edm::ValidityInterval &oValidity) {
91  for (unsigned int idx = 0; idx < m_entries.size(); ++idx) {
92  const auto &entry = m_entries[idx];
93 
94  // is within an entry ?
95  if (edm::contains(entry.m_validityRange, iosv.eventID())) {
96  m_currentEntryValid = true;
98  oValidity = edm::ValidityInterval(edm::IOVSyncValue(entry.m_validityRange.startEventID()),
99  edm::IOVSyncValue(entry.m_validityRange.endEventID()));
100  return;
101  }
102  }
103 
104  // not within any entry
105  m_currentEntryValid = false;
106  m_currentEntry = 0;
107 
108  edm::LogInfo("") << "No configuration entry found for event " << iosv.eventID()
109  << ", no optical functions will be available.";
110 
111  const edm::EventID start(iosv.eventID().run(), iosv.eventID().luminosityBlock(), iosv.eventID().event());
112  const edm::EventID end(iosv.eventID().run(), iosv.eventID().luminosityBlock(), iosv.eventID().event());
114 }
115 
116 //----------------------------------------------------------------------------------------------------
117 
118 std::unique_ptr<LHCOpticalFunctionsSetCollection> CTPPSOpticalFunctionsESSource::produce(const CTPPSOpticsRcd &) {
119  // prepare output, empty by default
120  auto output = std::make_unique<LHCOpticalFunctionsSetCollection>();
121 
122  // fill the output
123  if (m_currentEntryValid) {
124  const auto &entry = m_entries[m_currentEntry];
125 
126  for (const auto &fi : entry.m_fileInfo) {
127  std::unordered_map<unsigned int, LHCOpticalFunctionsSet> xa_data;
128 
129  for (const auto &rpi : entry.m_rpInfo) {
130  LHCOpticalFunctionsSet fcn(fi.m_fileName, rpi.second.m_dirName, rpi.second.m_scoringPlaneZ);
131  xa_data.emplace(rpi.first, std::move(fcn));
132  }
133 
134  output->emplace(fi.m_xangle, xa_data);
135  }
136  }
137 
138  // commit the output
139  return output;
140 }
141 
142 //----------------------------------------------------------------------------------------------------
143 
146 
147  desc.add<std::string>("label", "")->setComment("label of the optics record");
148 
149  edm::ParameterSetDescription config_desc;
150 
151  config_desc.add<edm::EventRange>("validityRange", edm::EventRange())->setComment("interval of validity");
152 
154  of_desc.add<double>("xangle")->setComment("half crossing angle value in urad");
155  of_desc.add<edm::FileInPath>("fileName")->setComment("ROOT file with optical functions");
156  std::vector<edm::ParameterSet> of;
157  config_desc.addVPSet("opticalFunctions", of_desc, of)
158  ->setComment("list of optical functions at different crossing angles");
159 
161  sp_desc.add<unsigned int>("rpId")->setComment("associated detector DetId");
162  sp_desc.add<std::string>("dirName")->setComment("associated path to the optical functions file");
163  sp_desc.add<double>("z")->setComment("longitudinal position at scoring plane/detector");
164  std::vector<edm::ParameterSet> sp;
165  config_desc.addVPSet("scoringPlanes", sp_desc, sp)->setComment("list of sensitive planes/detectors stations");
166 
167  std::vector<edm::ParameterSet> config;
168  desc.addVPSet("configuration", config_desc, sp)->setComment("list of configuration blocks");
169 
170  descriptions.add("ctppsOpticalFunctionsESSource", desc);
171 }
172 
173 //----------------------------------------------------------------------------------------------------
174 
Definition: start.py:1
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:165
void setComment(std::string const &value)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
bool contains(EventRange const &lh, EventID const &rh)
Definition: EventRange.cc:37
Loads optical functions from ROOT files.
dictionary config
Read in AllInOne config in JSON format.
Definition: DMR_cfg.py:21
std::pair< Time_t, Time_t > ValidityInterval
Definition: Time.h:17
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:39
~CTPPSOpticalFunctionsESSource() override=default
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setIntervalFor(const edm::eventsetup::EventSetupRecordKey &, const edm::IOVSyncValue &, edm::ValidityInterval &) override
Log< level::Info, false > LogInfo
RunNumber_t run() const
Definition: EventID.h:38
#define DEFINE_FWK_EVENTSETUP_SOURCE(type)
Definition: SourceFactory.h:92
Set of optical functions corresponding to one scoring plane along LHC.
void fcn(int &, double *, double &, double *, int)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const EventID & eventID() const
Definition: IOVSyncValue.h:40
std::unordered_map< unsigned int, RPInfo > m_rpInfo
std::unique_ptr< LHCOpticalFunctionsSetCollection > produce(const CTPPSOpticsRcd &)
CTPPSOpticalFunctionsESSource(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &)
def move(src, dest)
Definition: eostools.py:511
EventNumber_t event() const
Definition: EventID.h:40