CMS 3D CMS Logo

CTPPSInterpolatedOpticalFunctionsESSource.cc
Go to the documentation of this file.
1 // authors: Jan Kaspar (jan.kaspar@gmail.com)
2 
8 
10 
12 
15 
19 
21 public:
24 
25  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
26 
27  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> produce(const CTPPSInterpolatedOpticsRcd &);
28 
29 private:
34  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> currentData_;
37  const bool useNewLHCInfo_;
38 };
39 
40 //----------------------------------------------------------------------------------------------------
41 //----------------------------------------------------------------------------------------------------
42 
44  : currentCrossingAngle_(LHCInfoCombined::crossingAngleInvalid),
45  currentDataValid_(false),
46  useNewLHCInfo_(iConfig.getParameter<bool>("useNewLHCInfo")) {
47  auto cc = setWhatProduced(this, iConfig.getParameter<std::string>("opticsLabel"));
48  opticsToken_ = cc.consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("opticsLabel")));
49  lhcInfoToken_ = cc.consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")));
50  lhcInfoPerLSToken_ = cc.consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("lhcInfoPerLSLabel")));
51  lhcInfoPerFillToken_ = cc.consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("lhcInfoPerFillLabel")));
52 }
53 
54 //----------------------------------------------------------------------------------------------------
55 
58 
59  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
60  desc.add<std::string>("lhcInfoPerFillLabel", "")->setComment("label of the LHCInfoPerFill record");
61  desc.add<std::string>("lhcInfoPerLSLabel", "")->setComment("label of the LHCInfoPerLS record");
62  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics records");
63  desc.add<bool>("useNewLHCInfo", false)->setComment("flag whether to use new LHCInfoPer* records or old LHCInfo");
64 
65  descriptions.add("ctppsInterpolatedOpticalFunctionsESSourceDefault", desc);
66 }
67 
68 //----------------------------------------------------------------------------------------------------
69 
70 std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> CTPPSInterpolatedOpticalFunctionsESSource::produce(
71  const CTPPSInterpolatedOpticsRcd &iRecord) {
72  // get the input data
73  LHCOpticalFunctionsSetCollection const &ofColl = iRecord.get(opticsToken_);
74 
75  auto lhcInfoCombined = LHCInfoCombined::createLHCInfoCombined<
79 
80  // is there anything to do?
81  if (currentDataValid_ && lhcInfoCombined.crossingAngle() == currentCrossingAngle_)
82  return currentData_;
83 
84  // is crossing angle reasonable (LHCInfo is correctly filled in DB)?
85  if (lhcInfoCombined.isCrossingAngleInvalid()) {
86  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource")
87  << "Invalid crossing angle, no optical functions produced.";
88 
89  currentDataValid_ = false;
91  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
92 
93  return currentData_;
94  }
95 
96  // set new crossing angle
97  currentCrossingAngle_ = lhcInfoCombined.crossingAngle();
98  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource")
99  << "Crossing angle has changed to " << currentCrossingAngle_ << ".";
100 
101  // is input optics available ?
102  if (ofColl.empty()) {
103  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource")
104  << "No input optics available, no optical functions produced.";
105 
106  currentDataValid_ = false;
107  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
108 
109  return currentData_;
110  }
111 
112  // regular case with single-xangle input
113  if (ofColl.size() == 1) {
114  const auto &it = ofColl.begin();
115 
116  // does the input xangle correspond to the actual one?
117  if (fabs(currentCrossingAngle_ - it->first) > 1e-6)
118  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource")
119  << "Cannot interpolate: input given only for xangle " << it->first << " while interpolation requested for "
120  << currentCrossingAngle_ << ".";
121 
122  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
123  for (const auto &rp_p : it->second) {
124  const auto rpId = rp_p.first;
125  LHCInterpolatedOpticalFunctionsSet iof(rp_p.second);
126  iof.initializeSplines();
127  currentData_->emplace(rpId, std::move(iof));
128  }
129 
130  currentDataValid_ = true;
131  }
132 
133  // regular case with multi-xangle input
134  if (ofColl.size() > 1) {
135  // find the closest xangle points for interpolation
136  auto it1 = ofColl.begin();
137  auto it2 = std::next(it1);
138 
139  if (currentCrossingAngle_ > it1->first) {
140  for (; it1 != ofColl.end(); ++it1) {
141  it2 = std::next(it1);
142 
143  if (it2 == ofColl.end()) {
144  it2 = it1;
145  it1 = std::prev(it1);
146  break;
147  }
148 
149  if (it1->first <= currentCrossingAngle_ && currentCrossingAngle_ < it2->first)
150  break;
151  }
152  }
153 
154  const auto &xangle1 = it1->first;
155  const auto &xangle2 = it2->first;
156 
157  const auto &ofs1 = it1->second;
158  const auto &ofs2 = it2->second;
159 
160  // do the interpoaltion RP by RP
161  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
162  for (const auto &rp_p : ofs1) {
163  const auto rpId = rp_p.first;
164  const auto &rp_it2 = ofs2.find(rpId);
165  if (rp_it2 == ofs2.end())
166  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "RP mismatch between ofs1 and ofs2.";
167 
168  const auto &of1 = rp_p.second;
169  const auto &of2 = rp_it2->second;
170 
171  const size_t num_xi_vals1 = of1.getXiValues().size();
172  const size_t num_xi_vals2 = of2.getXiValues().size();
173 
174  if (num_xi_vals1 != num_xi_vals2)
175  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Size mismatch between ofs1 and ofs2.";
176 
177  const size_t num_xi_vals = num_xi_vals1;
178 
180  iof.m_z = of1.getScoringPlaneZ();
182  iof.m_xi_values.resize(num_xi_vals);
183 
184  for (size_t fi = 0; fi < of1.getFcnValues().size(); ++fi) {
185  iof.m_fcn_values[fi].resize(num_xi_vals);
186 
187  for (size_t pi = 0; pi < num_xi_vals; ++pi) {
188  double xi = of1.getXiValues()[pi];
189  double xi_control = of2.getXiValues()[pi];
190 
191  if (fabs(xi - xi_control) > 1e-6)
192  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Xi mismatch between ofs1 and ofs2.";
193 
194  iof.m_xi_values[pi] = xi;
195 
196  double v1 = of1.getFcnValues()[fi][pi];
197  double v2 = of2.getFcnValues()[fi][pi];
198  iof.m_fcn_values[fi][pi] = v1 + (v2 - v1) / (xangle2 - xangle1) * (currentCrossingAngle_ - xangle1);
199  }
200  }
201 
202  iof.initializeSplines();
203 
204  currentDataValid_ = true;
205  currentData_->emplace(rpId, std::move(iof));
206  }
207  }
208 
209  return currentData_;
210 }
211 
212 //----------------------------------------------------------------------------------------------------
213 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:166
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ESGetToken< LHCInfoPerLS, LHCInfoPerLSRcd > lhcInfoPerLSToken_
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > currentData_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
static constexpr float crossingAngleInvalid
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > produce(const CTPPSInterpolatedOpticsRcd &)
const Double_t pi
std::vector< double > m_xi_values
void initializeSplines()
builds splines from m_*_values fields
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
Log< level::Info, false > LogInfo
edm::ESGetToken< LHCOpticalFunctionsSetCollection, CTPPSOpticsRcd > opticsToken_
std::vector< std::vector< double > > m_fcn_values
length unit cm
Set of optical functions corresponding to one scoring plane along LHC, including splines for interpol...
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< LHCInfoPerFill, LHCInfoPerFillRcd > lhcInfoPerFillToken_
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
Collection of optical functions for two crossing angle values and various scoring planes...
double m_z
position of the scoring plane, in LHC/TOTEM convention, cm
def move(src, dest)
Definition: eostools.py:511
static LHCInfoCombined createLHCInfoCombined(const edm::eventsetup::DependentRecordImplementation< RecordT, ListT > &iRecord, const edm::ESGetToken< LHCInfoPerLS, LHCInfoPerLSRcd > &tokenInfoPerLS, const edm::ESGetToken< LHCInfoPerFill, LHCInfoPerFillRcd > &tokenInfoPerFill, const edm::ESGetToken< LHCInfo, LHCInfoRcd > &tokenInfo, bool useNewLHCInfo)