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 
20 public:
23 
24  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
25 
26  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> produce(const CTPPSInterpolatedOpticsRcd &);
27 
28 private:
31  std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> currentData_;
34 };
35 
36 //----------------------------------------------------------------------------------------------------
37 //----------------------------------------------------------------------------------------------------
38 
40  : currentCrossingAngle_(-1.), currentDataValid_(false) {
41  auto cc = setWhatProduced(this, iConfig.getParameter<std::string>("opticsLabel"));
42  opticsToken_ = cc.consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("opticsLabel")));
43  lhcInfoToken_ = cc.consumes(edm::ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")));
44 }
45 
46 //----------------------------------------------------------------------------------------------------
47 
50 
51  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
52  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics records");
53 
54  descriptions.add("ctppsInterpolatedOpticalFunctionsESSource", desc);
55 }
56 
57 //----------------------------------------------------------------------------------------------------
58 
59 std::shared_ptr<LHCInterpolatedOpticalFunctionsSetCollection> CTPPSInterpolatedOpticalFunctionsESSource::produce(
60  const CTPPSInterpolatedOpticsRcd &iRecord) {
61  // get the input data
62  LHCOpticalFunctionsSetCollection const &ofColl = iRecord.get(opticsToken_);
63 
64  LHCInfo const &lhcInfo = iRecord.get(lhcInfoToken_);
65 
66  // is there anything to do?
68  return currentData_;
69 
70  // is crossing angle reasonable (LHCInfo is correctly filled in DB)?
71  if (lhcInfo.crossingAngle() == 0.) {
72  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource")
73  << "Invalid crossing angle, no optical functions produced.";
74 
75  currentDataValid_ = false;
77  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
78 
79  return currentData_;
80  }
81 
82  // set new crossing angle
84  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource")
85  << "Crossing angle has changed to " << currentCrossingAngle_ << ".";
86 
87  // is input optics available ?
88  if (ofColl.empty()) {
89  edm::LogInfo("CTPPSInterpolatedOpticalFunctionsESSource")
90  << "No input optics available, no optical functions produced.";
91 
92  currentDataValid_ = false;
93  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
94 
95  return currentData_;
96  }
97 
98  // regular case with single-xangle input
99  if (ofColl.size() == 1) {
100  const auto &it = ofColl.begin();
101 
102  // does the input xangle correspond to the actual one?
103  if (fabs(currentCrossingAngle_ - it->first) > 1e-6)
104  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource")
105  << "Cannot interpolate: input given only for xangle " << it->first << " while interpolation requested for "
106  << currentCrossingAngle_ << ".";
107 
108  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
109  for (const auto &rp_p : it->second) {
110  const auto rpId = rp_p.first;
111  LHCInterpolatedOpticalFunctionsSet iof(rp_p.second);
112  iof.initializeSplines();
113  currentData_->emplace(rpId, std::move(iof));
114  }
115 
116  currentDataValid_ = true;
117  }
118 
119  // regular case with multi-xangle input
120  if (ofColl.size() > 1) {
121  // find the closest xangle points for interpolation
122  auto it1 = ofColl.begin();
123  auto it2 = std::next(it1);
124 
125  if (currentCrossingAngle_ > it1->first) {
126  for (; it1 != ofColl.end(); ++it1) {
127  it2 = std::next(it1);
128 
129  if (it2 == ofColl.end()) {
130  it2 = it1;
131  it1 = std::prev(it1);
132  break;
133  }
134 
135  if (it1->first <= currentCrossingAngle_ && currentCrossingAngle_ < it2->first)
136  break;
137  }
138  }
139 
140  const auto &xangle1 = it1->first;
141  const auto &xangle2 = it2->first;
142 
143  const auto &ofs1 = it1->second;
144  const auto &ofs2 = it2->second;
145 
146  // do the interpoaltion RP by RP
147  currentData_ = std::make_shared<LHCInterpolatedOpticalFunctionsSetCollection>();
148  for (const auto &rp_p : ofs1) {
149  const auto rpId = rp_p.first;
150  const auto &rp_it2 = ofs2.find(rpId);
151  if (rp_it2 == ofs2.end())
152  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "RP mismatch between ofs1 and ofs2.";
153 
154  const auto &of1 = rp_p.second;
155  const auto &of2 = rp_it2->second;
156 
157  const size_t num_xi_vals1 = of1.getXiValues().size();
158  const size_t num_xi_vals2 = of2.getXiValues().size();
159 
160  if (num_xi_vals1 != num_xi_vals2)
161  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Size mismatch between ofs1 and ofs2.";
162 
163  const size_t num_xi_vals = num_xi_vals1;
164 
166  iof.m_z = of1.getScoringPlaneZ();
168  iof.m_xi_values.resize(num_xi_vals);
169 
170  for (size_t fi = 0; fi < of1.getFcnValues().size(); ++fi) {
171  iof.m_fcn_values[fi].resize(num_xi_vals);
172 
173  for (size_t pi = 0; pi < num_xi_vals; ++pi) {
174  double xi = of1.getXiValues()[pi];
175  double xi_control = of2.getXiValues()[pi];
176 
177  if (fabs(xi - xi_control) > 1e-6)
178  throw cms::Exception("CTPPSInterpolatedOpticalFunctionsESSource") << "Xi mismatch between ofs1 and ofs2.";
179 
180  iof.m_xi_values[pi] = xi;
181 
182  double v1 = of1.getFcnValues()[fi][pi];
183  double v2 = of2.getFcnValues()[fi][pi];
184  iof.m_fcn_values[fi][pi] = v1 + (v2 - v1) / (xangle2 - xangle1) * (currentCrossingAngle_ - xangle1);
185  }
186  }
187 
188  iof.initializeSplines();
189 
190  currentDataValid_ = true;
191  currentData_->emplace(rpId, std::move(iof));
192  }
193  }
194 
195  return currentData_;
196 }
197 
198 //----------------------------------------------------------------------------------------------------
199 
edm::ESInputTag
Definition: ESInputTag.h:87
MessageLogger.h
funct::false
false
Definition: Factorize.h:29
CTPPSInterpolatedOpticalFunctionsESSource::currentData_
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > currentData_
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:31
ESHandle.h
LHCOpticalFunctionsSet::m_fcn_values
std::vector< std::vector< double > > m_fcn_values
length unit cm
Definition: LHCOpticalFunctionsSet.h:36
edm::ESProducer::setWhatProduced
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
LHCInfo.h
ESProducer.h
LHCInfo
Definition: LHCInfo.h:12
LHCOpticalFunctionsSet::nFunctions
Definition: LHCOpticalFunctionsSet.h:16
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
CTPPSOpticsRcd.h
CTPPSInterpolatedOpticalFunctionsESSource::lhcInfoToken_
edm::ESGetToken< LHCInfo, LHCInfoRcd > lhcInfoToken_
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:30
LHCInfo::crossingAngle
const float crossingAngle() const
Definition: LHCInfo.cc:182
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
CTPPSInterpolatedOpticalFunctionsESSource::produce
std::shared_ptr< LHCInterpolatedOpticalFunctionsSetCollection > produce(const CTPPSInterpolatedOpticsRcd &)
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:59
LHCInterpolatedOpticalFunctionsSetCollection.h
CTPPSInterpolatedOpticalFunctionsESSource::~CTPPSInterpolatedOpticalFunctionsESSource
~CTPPSInterpolatedOpticalFunctionsESSource() override
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:22
CTPPSInterpolatedOpticalFunctionsESSource::currentCrossingAngle_
float currentCrossingAngle_
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:32
edm::eventsetup::DependentRecordImplementation::get
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
Definition: DependentRecordImplementation.h:109
LHCOpticalFunctionsSet::m_z
double m_z
position of the scoring plane, in LHC/TOTEM convention, cm
Definition: LHCOpticalFunctionsSet.h:33
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:112
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
LHCInterpolatedOpticalFunctionsSet
Set of optical functions corresponding to one scoring plane along LHC, including splines for interpol...
Definition: LHCInterpolatedOpticalFunctionsSet.h:14
CTPPSInterpolatedOpticalFunctionsESSource::currentDataValid_
bool currentDataValid_
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:33
edm::ParameterSet
Definition: ParameterSet.h:47
CTPPSInterpolatedOpticsRcd.h
cc
CTPPSInterpolatedOpticalFunctionsESSource::opticsToken_
edm::ESGetToken< LHCOpticalFunctionsSetCollection, CTPPSOpticsRcd > opticsToken_
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:29
edm::ESGetToken< LHCOpticalFunctionsSetCollection, CTPPSOpticsRcd >
profile_2016_postTS2_cff.rpId
rpId
Definition: profile_2016_postTS2_cff.py:21
protons_cff.xi
xi
Definition: protons_cff.py:35
ModuleFactory.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
DEFINE_FWK_EVENTSETUP_MODULE
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
Exception
Definition: hltDiff.cc:245
CTPPSInterpolatedOpticalFunctionsESSource
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:19
LHCInterpolatedOpticalFunctionsSet::initializeSplines
void initializeSplines()
builds splines from m_*_values fields
Definition: LHCInterpolatedOpticalFunctionsSet.cc:7
LHCOpticalFunctionsSetCollection.h
CTPPSInterpolatedOpticalFunctionsESSource::CTPPSInterpolatedOpticalFunctionsESSource
CTPPSInterpolatedOpticalFunctionsESSource(const edm::ParameterSet &)
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:39
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESProducer
Definition: ESProducer.h:104
pi
const Double_t pi
Definition: trackSplitPlot.h:36
cms::Exception
Definition: Exception.h:70
LHCOpticalFunctionsSetCollection
Collection of optical functions for two crossing angle values and various scoring planes....
Definition: LHCOpticalFunctionsSetCollection.h:16
ParameterSet.h
CTPPSInterpolatedOpticalFunctionsESSource::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: CTPPSInterpolatedOpticalFunctionsESSource.cc:48
LHCInfoRcd.h
LHCOpticalFunctionsSet::m_xi_values
std::vector< double > m_xi_values
Definition: LHCOpticalFunctionsSet.h:35
GetRecoTauVFromDQM_MC_cff.next
next
Definition: GetRecoTauVFromDQM_MC_cff.py:31
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
CTPPSInterpolatedOpticsRcd
Definition: CTPPSInterpolatedOpticsRcd.h:13