CMS 3D CMS Logo

CTPPSInterpolatedOpticalFunctionsESSource.cc
Go to the documentation of this file.
1 // authors: Jan Kaspar (jan.kaspar@gmail.com)
2 
8 
10 
14 
18 
19 
21 {
22  public:
25 
26  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
27 
28  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> produce(const CTPPSInterpolatedOpticsRcd&);
29 
30  private:
33 
36  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> currentData_;
37 };
38 
39 //----------------------------------------------------------------------------------------------------
40 //----------------------------------------------------------------------------------------------------
41 
43  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
44  opticsLabel_(iConfig.getParameter<std::string>("opticsLabel")),
47 {
49 }
50 
51 //----------------------------------------------------------------------------------------------------
52 
54 {
56 
57  desc.add<std::string>("lhcInfoLabel", "")
58  ->setComment("label of the LHCInfo record");
59 
60  desc.add<std::string>("opticsLabel", "")
61  ->setComment("label of the optics record");
62 
63  descriptions.add("ctppsInterpolatedOpticalFunctionsESSource", desc);
64 }
65 
66 //----------------------------------------------------------------------------------------------------
67 
68 std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> CTPPSInterpolatedOpticalFunctionsESSource::produce(const CTPPSInterpolatedOpticsRcd& iRecord)
69 {
70  // get the input data
72  iRecord.getRecord<CTPPSOpticsRcd>().get(opticsLabel_, hOFColl);
73 
74  edm::ESHandle<LHCInfo> hLHCInfo;
75  iRecord.getRecord<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
76 
77  // is there anything to do?
79  return currentData_;
80 
81  // is crossing angle reasonable (LHCInfo is correctly filled in DB)?
82  if (hLHCInfo->crossingAngle() == 0.)
83  {
84  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource") << "Invalid crossing angle, no optical functions produced.";
85 
86  currentDataValid_ = false;
88  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
89 
90  return currentData_;
91  }
92 
93  // set new crossing angle
95  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource") << "Crossing angle has changed to " << currentCrossingAngle_ << ".";
96 
97  // is input optics available ?
98  if (hOFColl->empty())
99  {
100  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource") << "No input optics available, no optical functions produced.";
101 
102  currentDataValid_ = false;
103  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
104 
105  return currentData_;
106  }
107 
108  // regular case with single-xangle input
109  if (hOFColl->size() == 1)
110  {
111  const auto &it = hOFColl->begin();
112 
113  // does the input xangle correspond to the actual one?
114  if (fabs(currentCrossingAngle_ - it->first) > 1e-6)
115  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Cannot interpolate: input given only for xangle "
116  << it->first << " while interpolation requested for " << currentCrossingAngle_ << ".";
117 
118  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
119  for (const auto &rp_p : it->second)
120  {
121  const auto rpId = rp_p.first;
122  LHCInterpolatedOpticalFunctionsSet iof(rp_p.second);
123  iof.initializeSplines();
124  currentData_->emplace(rpId, std::move(iof));
125  }
126 
127  currentDataValid_ = true;
128  }
129 
130  // regular case with multi-xangle input
131  if (hOFColl->size() > 1)
132  {
133  // find the closest xangle points for interpolation
134  auto it1 = hOFColl->begin();
135  auto it2 = std::next(it1);
136 
137  if (currentCrossingAngle_ > it1->first)
138  {
139  for (; it1 != hOFColl->end(); ++it1)
140  {
141  it2 = std::next(it1);
142 
143  if (it2 == hOFColl->end())
144  {
145  it2 = it1;
146  it1 = std::prev(it1);
147  break;
148  }
149 
150  if (it1->first <= currentCrossingAngle_ && currentCrossingAngle_ < it2->first)
151  break;
152  }
153  }
154 
155  const auto &xangle1 = it1->first;
156  const auto &xangle2 = it2->first;
157 
158  const auto &ofs1 = it1->second;
159  const auto &ofs2 = it2->second;
160 
161  // do the interpoaltion RP by RP
162  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
163  for (const auto &rp_p : ofs1)
164  {
165  const auto rpId = rp_p.first;
166  const auto &rp_it2 = ofs2.find(rpId);
167  if (rp_it2 == ofs2.end())
168  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "RP mismatch between ofs1 and ofs2.";
169 
170  const auto &of1 = rp_p.second;
171  const auto &of2 = rp_it2->second;
172 
173  const size_t num_xi_vals1 = of1.getXiValues().size();
174  const size_t num_xi_vals2 = of2.getXiValues().size();
175 
176  if (num_xi_vals1 != num_xi_vals2)
177  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Size mismatch between ofs1 and ofs2.";
178 
179  const size_t num_xi_vals = num_xi_vals1;
180 
182  iof.m_z = of1.getScoringPlaneZ();
184  iof.m_xi_values.resize(num_xi_vals);
185 
186  for (size_t fi = 0; fi < of1.getFcnValues().size(); ++fi)
187  {
188  iof.m_fcn_values[fi].resize(num_xi_vals);
189 
190  for (size_t pi = 0; pi < num_xi_vals; ++pi)
191  {
192  double xi = of1.getXiValues()[pi];
193  double xi_control = of2.getXiValues()[pi];
194 
195  if (fabs(xi - xi_control) > 1e-6)
196  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Xi mismatch between ofs1 and ofs2.";
197 
198  iof.m_xi_values[pi] = xi;
199 
200  double v1 = of1.getFcnValues()[fi][pi];
201  double v2 = of2.getFcnValues()[fi][pi];
202  iof.m_fcn_values[fi][pi] = v1 + (v2 - v1) / (xangle2 - xangle1) * (currentCrossingAngle_ - xangle1);
203  }
204  }
205 
206  iof.initializeSplines();
207 
208  currentDataValid_ = true;
209  currentData_->emplace(rpId, std::move(iof));
210  }
211  }
212 
213  return currentData_;
214 }
215 
216 //----------------------------------------------------------------------------------------------------
217 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:124
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > currentData_
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > produce(const CTPPSInterpolatedOpticsRcd &)
const Double_t pi
float const crossingAngle() const
Definition: LHCInfo.cc:176
std::vector< double > m_xi_values
void initializeSplines()
builds splines from m_*_values fields
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
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)
double m_z
position of the scoring plane, in LHC/TOTEM convention, cm
def move(src, dest)
Definition: eostools.py:511