CMS 3D CMS Logo

ProtonReconstructionAlgorithm.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  * Laurent Forthomme
5  ****************************************************************************/
6 
8 
10 
14 
15 #include "TMinuitMinimizer.h"
16 
17 using namespace std;
18 using namespace edm;
19 
20 //----------------------------------------------------------------------------------------------------
21 
23  bool improved_estimate,
24  const std::string &multiRPAlgorithm,
25  unsigned int verbosity)
26  : verbosity_(verbosity),
27  fitVtxY_(fit_vtx_y),
28  useImprovedInitialEstimate_(improved_estimate),
29  multi_rp_algorithm_(mrpaUndefined),
30  initialized_(false),
31  fitter_(new ROOT::Fit::Fitter),
32  chiSquareCalculator_(new ChiSquareCalculator) {
33  // needed for thread safety
34  TMinuitMinimizer::UseStaticMinuit(false);
35 
36  // check and set multi-RP algorithm
37  if (multiRPAlgorithm == "chi2")
39  else if (multiRPAlgorithm == "newton")
41  else if (multiRPAlgorithm == "anal-iter")
43 
45  throw cms::Exception("ProtonReconstructionAlgorithm") << "Algorithm '" << multiRPAlgorithm << "' not understood.";
46 
47  // initialise fitter
48  double pStart[] = {0, 0, 0, 0};
49  fitter_->SetFCN(4, *chiSquareCalculator_, pStart, 0, true);
50 }
51 
52 //----------------------------------------------------------------------------------------------------
53 
55  // reset cache
56  release();
57 
58  // build optics data for each object
59  for (const auto &p : opticalFunctions) {
60  const LHCInterpolatedOpticalFunctionsSet &ofs = p.second;
61 
62  // make record
63  RPOpticsData rpod;
64  rpod.optics = &p.second;
70 
71  vector<double> xiValues =
72  ofs.getXiValues(); // local copy made since the TSpline constructor needs non-const parameters
73  vector<double> xDValues = ofs.getFcnValues()[LHCOpticalFunctionsSet::exd];
74  rpod.s_xi_vs_x_d = make_shared<TSpline3>("", xDValues.data(), xiValues.data(), xiValues.size());
75 
76  // calculate auxiliary data
77  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in = {0., 0., 0., 0., 0.};
79  rpod.optics->transport(k_in, k_out);
80  rpod.x0 = k_out.x;
81  rpod.y0 = k_out.y;
82 
84  rpod.ch0 -= rpod.x0;
85 
87 
88  // insert record
89  const CTPPSDetId rpId(p.first);
90  m_rp_optics_.emplace(rpId, std::move(rpod));
91  }
92 
93  // update settings
94  initialized_ = true;
95 }
96 
97 //----------------------------------------------------------------------------------------------------
98 
99 void ProtonReconstructionAlgorithm::doLinearFit(const std::vector<double> &vx,
100  const std::vector<double> &vy,
101  double &b,
102  double &a) {
103  double s_1 = 0., s_x = 0., s_xx = 0., s_y = 0., s_xy = 0.;
104  for (unsigned int i = 0; i < vx.size(); ++i) {
105  s_1 += 1.;
106  s_x += vx[i];
107  s_xx += vx[i] * vx[i];
108  s_y += vy[i];
109  s_xy += vx[i] * vy[i];
110  }
111 
112  const double d = s_xx * s_1 - s_x * s_x;
113  a = (s_1 * s_xy - s_x * s_y) / d;
114  b = (-s_x * s_xy + s_xx * s_y) / d;
115 }
116 
117 //----------------------------------------------------------------------------------------------------
118 
120  initialized_ = false;
121 
122  m_rp_optics_.clear();
123 }
124 
125 //----------------------------------------------------------------------------------------------------
126 
128  // extract proton parameters
130  0., parameters[1], parameters[3], parameters[2], parameters[0]};
131 
132  // calculate chi^2 by looping over hits
133  double s2 = 0.;
134 
135  for (const auto &track : *tracks) {
136  const CTPPSDetId rpId(track->rpId());
137 
138  // transport proton to the RP
139  auto oit = m_rp_optics->find(rpId);
141  oit->second.optics->transport(k_in, k_out);
142 
143  // proton position wrt. beam
144  const double x = k_out.x - oit->second.x0;
145  const double y = k_out.y - oit->second.y0;
146 
147  // calculate chi^2 contributions, convert track data mm --> cm
148  const double x_diff_norm = (x - track->x() * 1E-1) / (track->xUnc() * 1E-1);
149  const double y_diff_norm = (y - track->y() * 1E-1) / (track->yUnc() * 1E-1);
150 
151  // increase chi^2
152  s2 += x_diff_norm * x_diff_norm + y_diff_norm * y_diff_norm;
153  }
154 
155  return s2;
156 }
157 
158 //----------------------------------------------------------------------------------------------------
159 
161  double x_N,
162  double x_F,
165  return (x_N - i_N.s_x_d_vs_xi->Eval(xi)) * i_F.s_L_x_vs_xi->Eval(xi) -
166  (x_F - i_F.s_x_d_vs_xi->Eval(xi)) * i_N.s_L_x_vs_xi->Eval(xi);
167 }
168 
169 //----------------------------------------------------------------------------------------------------
170 
172  const LHCInfo &lhcInfo,
173  std::ostream &os) const {
174  // make sure optics is available for all tracks
175  for (const auto &it : tracks) {
176  auto oit = m_rp_optics_.find(it->rpId());
177  if (oit == m_rp_optics_.end())
178  throw cms::Exception("ProtonReconstructionAlgorithm")
179  << "Optics data not available for RP " << it->rpId() << ", i.e. " << CTPPSDetId(it->rpId()) << ".";
180  }
181 
182  // initial estimate of xi and th_x
183  double xi_init = 0., th_x_init = 0.;
184 
186  double x_N = tracks[0]->x() * 1E-1, // conversion: mm --> cm
187  x_F = tracks[1]->x() * 1E-1;
188 
189  const RPOpticsData &i_N = m_rp_optics_.find(tracks[0]->rpId())->second,
190  &i_F = m_rp_optics_.find(tracks[1]->rpId())->second;
191 
192  const double a = i_F.ch1 * i_N.la1 - i_N.ch1 * i_F.la1;
193  const double b =
194  i_F.ch0 * i_N.la1 - i_N.ch0 * i_F.la1 + i_F.ch1 * i_N.la0 - i_N.ch1 * i_F.la0 + x_N * i_F.la1 - x_F * i_N.la1;
195  const double c = x_N * i_F.la0 - x_F * i_N.la0 + i_F.ch0 * i_N.la0 - i_N.ch0 * i_F.la0;
196  const double D = b * b - 4. * a * c;
197  const double sqrt_D = (D >= 0.) ? sqrt(D) : 0.;
198 
199  xi_init = (-b + sqrt_D) / 2. / a;
200  th_x_init = (x_N - i_N.ch0 - i_N.ch1 * xi_init) / (i_N.la0 + i_N.la1 * xi_init);
201  } else {
202  double s_xi0 = 0., s_1 = 0.;
203  for (const auto &track : tracks) {
204  auto oit = m_rp_optics_.find(track->rpId());
205  double xi = oit->second.s_xi_vs_x_d->Eval(track->x() * 1E-1 + oit->second.x0); // conversion: mm --> cm
206 
207  s_1 += 1.;
208  s_xi0 += xi;
209  }
210 
211  xi_init = s_xi0 / s_1;
212  }
213 
214  if (!std::isfinite(xi_init))
215  xi_init = 0.;
216  if (!std::isfinite(th_x_init))
217  th_x_init = 0.;
218 
219  // initial estimate of th_y and vtx_y
220  double y[2] = {0}, v_y[2] = {0}, L_y[2] = {0};
221  unsigned int y_idx = 0;
222  for (const auto &track : tracks) {
223  if (y_idx >= 2)
224  break;
225 
226  auto oit = m_rp_optics_.find(track->rpId());
227 
228  y[y_idx] = track->y() * 1E-1 - oit->second.s_y_d_vs_xi->Eval(xi_init); // track y: mm --> cm
229  v_y[y_idx] = oit->second.s_v_y_vs_xi->Eval(xi_init);
230  L_y[y_idx] = oit->second.s_L_y_vs_xi->Eval(xi_init);
231 
232  y_idx++;
233  }
234 
235  double vtx_y_init = 0.;
236  double th_y_init = 0.;
237 
238  if (fitVtxY_) {
239  const double det_y = v_y[0] * L_y[1] - L_y[0] * v_y[1];
240  vtx_y_init = (det_y != 0.) ? (L_y[1] * y[0] - L_y[0] * y[1]) / det_y : 0.;
241  th_y_init = (det_y != 0.) ? (v_y[0] * y[1] - v_y[1] * y[0]) / det_y : 0.;
242  } else {
243  vtx_y_init = 0.;
244  th_y_init = (y[1] / L_y[1] + y[0] / L_y[0]) / 2.;
245  }
246 
247  if (!std::isfinite(vtx_y_init))
248  vtx_y_init = 0.;
249  if (!std::isfinite(th_y_init))
250  th_y_init = 0.;
251 
252  unsigned int armId = CTPPSDetId((*tracks.begin())->rpId()).arm();
253 
254  if (verbosity_)
255  os << "ProtonReconstructionAlgorithm::reconstructFromMultiRP(" << armId << ")" << std::endl
256  << " initial estimate: xi_init = " << xi_init << ", th_x_init = " << th_x_init
257  << ", th_y_init = " << th_y_init << ", vtx_y_init = " << vtx_y_init << "." << std::endl;
258 
259  // prepare result containers
260  bool valid = false;
261  double chi2 = 0.;
262  double xi = 0., th_x = 0., th_y = 0., vtx_y = 0.;
263 
264  using FP = reco::ForwardProton;
266 
267  if (multi_rp_algorithm_ == mrpaChi2) {
268  // minimisation
269  fitter_->Config().ParSettings(0).Set("xi", xi_init, 0.005);
270  fitter_->Config().ParSettings(1).Set("th_x", th_x_init, 2E-6);
271  fitter_->Config().ParSettings(2).Set("th_y", th_y_init, 2E-6);
272  fitter_->Config().ParSettings(3).Set("vtx_y", vtx_y_init, 10E-6);
273 
274  if (!fitVtxY_)
275  fitter_->Config().ParSettings(3).Fix();
276 
277  chiSquareCalculator_->tracks = &tracks;
278  chiSquareCalculator_->m_rp_optics = &m_rp_optics_;
279 
280  fitter_->FitFCN();
281  fitter_->FitFCN(); // second minimisation in case the first one had troubles
282 
283  // extract proton parameters
284  const ROOT::Fit::FitResult &result = fitter_->Result();
285  const double *params = result.GetParams();
286 
287  valid = result.IsValid();
288  chi2 = result.Chi2();
289  xi = params[0];
290  th_x = params[1];
291  th_y = params[2];
292  vtx_y = params[3];
293 
294  map<unsigned int, signed int> index_map = {
295  {(unsigned int)FP::Index::xi, 0},
296  {(unsigned int)FP::Index::th_x, 1},
297  {(unsigned int)FP::Index::th_y, 2},
298  {(unsigned int)FP::Index::vtx_y, ((fitVtxY_) ? 3 : -1)},
299  {(unsigned int)FP::Index::vtx_x, -1},
300  };
301 
302  for (unsigned int i = 0; i < (unsigned int)FP::Index::num_indices; ++i) {
303  signed int fit_i = index_map[i];
304 
305  for (unsigned int j = 0; j < (unsigned int)FP::Index::num_indices; ++j) {
306  signed int fit_j = index_map[j];
307 
308  cm(i, j) = (fit_i >= 0 && fit_j >= 0) ? result.CovMatrix(fit_i, fit_j) : 0.;
309  }
310  }
311  }
312 
314  // settings
315  const unsigned int maxIterations = 100;
316  const double maxXiDiff = 1E-6;
317  const double xi_ep = 1E-5;
318 
319  // collect input
320  double x_N = tracks[0]->x() * 1E-1, // conversion: mm --> cm
321  x_F = tracks[1]->x() * 1E-1;
322 
323  double y_N = tracks[0]->y() * 1E-1, // conversion: mm --> cm
324  y_F = tracks[1]->y() * 1E-1;
325 
326  const RPOpticsData &i_N = m_rp_optics_.find(tracks[0]->rpId())->second,
327  &i_F = m_rp_optics_.find(tracks[1]->rpId())->second;
328 
329  // horizontal reconstruction - run iterations
330  valid = false;
331  double xi_prev = xi_init;
332  for (unsigned int it = 0; it < maxIterations; ++it) {
334  const double g = newtonGoalFcn(xi_prev, x_N, x_F, i_N, i_F);
335  const double gp = (newtonGoalFcn(xi_prev + xi_ep, x_N, x_F, i_N, i_F) - g) / xi_ep;
336 
337  xi = xi_prev - g / gp;
338  }
339 
341  const double d_x_eff_N = i_N.s_x_d_vs_xi->Eval(xi_prev) / xi_prev;
342  const double d_x_eff_F = i_F.s_x_d_vs_xi->Eval(xi_prev) / xi_prev;
343 
344  const double l_x_N = i_N.s_L_x_vs_xi->Eval(xi_prev);
345  const double l_x_F = i_F.s_L_x_vs_xi->Eval(xi_prev);
346 
347  xi = (l_x_N * x_F - l_x_F * x_N) / (l_x_N * d_x_eff_F - l_x_F * d_x_eff_N);
348  }
349 
350  if (abs(xi - xi_prev) < maxXiDiff) {
351  valid = true;
352  break;
353  }
354 
355  xi_prev = xi;
356  }
357 
358  th_x = (x_F - i_F.s_x_d_vs_xi->Eval(xi)) / i_F.s_L_x_vs_xi->Eval(xi);
359 
360  // vertical reconstruction
361  const double y_eff_N = y_N - i_N.s_y_d_vs_xi->Eval(xi);
362  const double y_eff_F = y_F - i_F.s_y_d_vs_xi->Eval(xi);
363 
364  const double l_y_N = i_N.s_L_y_vs_xi->Eval(xi);
365  const double l_y_F = i_F.s_L_y_vs_xi->Eval(xi);
366 
367  const double v_y_N = i_N.s_v_y_vs_xi->Eval(xi);
368  const double v_y_F = i_F.s_v_y_vs_xi->Eval(xi);
369 
370  const double det = l_y_N * v_y_F - l_y_F * v_y_N;
371 
372  if (fitVtxY_) {
373  th_y = (y_eff_N * v_y_F - y_eff_F * v_y_N) / det;
374  vtx_y = (l_y_N * y_eff_F - l_y_F * y_eff_N) / det;
375  } else {
376  th_y = (y_eff_N / l_y_N + y_eff_F / l_y_F) / 2.;
377  vtx_y = 0.;
378  }
379 
380  // covariance matrix kept empty - not to compromise the speed of these special algorithms
381  }
382 
383  // print results
384  if (verbosity_)
385  os << " fit valid=" << valid << std::endl
386  << " xi=" << xi << ", th_x=" << th_x << ", th_y=" << th_y << ", vtx_y=" << vtx_y << ", chiSq = " << chi2
387  << std::endl;
388 
389  // save reco candidate
390  const double sign_z = (armId == 0) ? +1. : -1.; // CMS convention
391  const FP::Point vertex(0., vtx_y, 0.);
392  const double cos_th = sqrt(1. - th_x * th_x - th_y * th_y);
393  const double p = lhcInfo.energy() * (1. - xi);
394  const FP::Vector momentum(-p * th_x, // the signs reflect change LHC --> CMS convention
395  +p * th_y,
396  sign_z * p * cos_th);
397  signed int ndf = 2. * tracks.size() - ((fitVtxY_) ? 4. : 3.);
398 
399  return reco::ForwardProton(chi2, ndf, vertex, momentum, xi, cm, FP::ReconstructionMethod::multiRP, tracks, valid);
400 }
401 
402 //----------------------------------------------------------------------------------------------------
403 
405  const LHCInfo &lhcInfo,
406  std::ostream &os) const {
407  CTPPSDetId rpId(track->rpId());
408 
409  if (verbosity_)
410  os << "reconstructFromSingleRP(" << rpId.arm() * 100 + rpId.station() * 10 + rpId.rp() << ")" << std::endl;
411 
412  // make sure optics is available for the track
413  auto oit = m_rp_optics_.find(track->rpId());
414  if (oit == m_rp_optics_.end())
415  throw cms::Exception("ProtonReconstructionAlgorithm")
416  << "Optics data not available for RP " << track->rpId() << ", i.e. " << rpId << ".";
417 
418  // rough estimate of xi and th_y from each track
419  const double x_full = track->x() * 1E-1 + oit->second.x0; // conversion mm --> cm
420  const double xi = oit->second.s_xi_vs_x_d->Eval(x_full);
421  const double L_y = oit->second.s_L_y_vs_xi->Eval(xi);
422  const double th_y = track->y() * 1E-1 / L_y; // conversion mm --> cm
423 
424  const double ep_x = 1E-6;
425  const double dxi_dx = (oit->second.s_xi_vs_x_d->Eval(x_full + ep_x) - xi) / ep_x;
426  const double xi_unc = abs(dxi_dx) * track->xUnc() * 1E-1; // conversion mm --> cm
427 
428  const double ep_xi = 1E-4;
429  const double dL_y_dxi = (oit->second.s_L_y_vs_xi->Eval(xi + ep_xi) - L_y) / ep_xi;
430  const double th_y_unc = th_y * sqrt(pow(track->yUnc() / track->y(), 2.) + pow(dL_y_dxi * xi_unc / L_y, 2.));
431 
432  if (verbosity_)
433  os << " xi = " << xi << " +- " << xi_unc << ", th_y = " << th_y << " +- " << th_y_unc << "." << std::endl;
434 
435  using FP = reco::ForwardProton;
436 
437  // save proton candidate
438  const double sign_z = (CTPPSDetId(track->rpId()).arm() == 0) ? +1. : -1.; // CMS convention
439  const FP::Point vertex(0., 0., 0.);
440  const double cos_th = sqrt(1. - th_y * th_y);
441  const double p = lhcInfo.energy() * (1. - xi);
442  const FP::Vector momentum(0., p * th_y, sign_z * p * cos_th);
443 
445  cm((int)FP::Index::xi, (int)FP::Index::xi) = xi_unc * xi_unc;
446  cm((int)FP::Index::th_y, (int)FP::Index::th_y) = th_y_unc * th_y_unc;
447 
449  trk.push_back(track);
450 
451  return reco::ForwardProton(0., 0, vertex, momentum, xi, cm, FP::ReconstructionMethod::singleRP, trk, true);
452 }
ProtonReconstructionAlgorithm::RPOpticsData::x0
double x0
beam horizontal position, cm
Definition: ProtonReconstructionAlgorithm.h:60
ProtonReconstructionAlgorithm::ChiSquareCalculator::operator()
double operator()(const double *parameters) const
Definition: ProtonReconstructionAlgorithm.cc:126
LHCInfo::energy
const float energy() const
Definition: LHCInfo.cc:190
HIPAlignmentAlgorithm_cfi.verbosity
verbosity
Definition: HIPAlignmentAlgorithm_cfi.py:7
DDAxes::y
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
LHCOpticalFunctionsSet::exd
Definition: LHCOpticalFunctionsSet.h:16
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
funct::false
false
Definition: Factorize.h:34
LHCInterpolatedOpticalFunctionsSet::Kinematics
proton kinematics description
Definition: LHCInterpolatedOpticalFunctionsSet.h:28
ProtonReconstructionAlgorithm::RPOpticsData::optics
const LHCInterpolatedOpticalFunctionsSet * optics
Definition: ProtonReconstructionAlgorithm.h:58
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
HLTSiStripMonitoring_cff.Fitter
Fitter
Definition: HLTSiStripMonitoring_cff.py:207
hybridSuperClusters_cfi.xi
xi
Definition: hybridSuperClusters_cfi.py:10
reco::ForwardProton
Definition: ForwardProton.h:22
LHCInfo
Definition: LHCInfo.h:12
ProtonReconstructionAlgorithm::newtonGoalFcn
static double newtonGoalFcn(double xi, double x_N, double x_F, const RPOpticsData &i_N, const RPOpticsData &i_F)
Definition: ProtonReconstructionAlgorithm.cc:159
LHCOpticalFunctionsSet::getFcnValues
const std::vector< std::vector< double > > & getFcnValues() const
Definition: LHCOpticalFunctionsSet.h:29
ProtonReconstructionAlgorithm::init
void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions)
Definition: ProtonReconstructionAlgorithm.cc:53
CTPPSLocalTrackLite.h
DDAxes::x
edm::RefVector< CTPPSLocalTrackLiteCollection >
indexGen.s2
s2
Definition: indexGen.py:107
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
ProtonReconstructionAlgorithm::RPOpticsData::la1
double la1
slope for linear approximation of
Definition: ProtonReconstructionAlgorithm.h:65
year_2016_postTS2_cff.rpId
rpId
Definition: year_2016_postTS2_cff.py:23
ProtonReconstructionAlgorithm::RPOpticsData::la0
double la0
intercept for linear approximation of
Definition: ProtonReconstructionAlgorithm.h:64
ProtonReconstructionAlgorithm::reconstructFromMultiRP
reco::ForwardProton reconstructFromMultiRP(const CTPPSLocalTrackLiteRefVector &tracks, const LHCInfo &lhcInfo, std::ostream &os) const
run proton reconstruction using multiple-RP strategy
Definition: ProtonReconstructionAlgorithm.cc:170
validateGeometry_cfg.valid
valid
Definition: validateGeometry_cfg.py:21
ProtonReconstructionAlgorithm::RPOpticsData::s_x_d_vs_xi
std::shared_ptr< const TSpline3 > s_x_d_vs_xi
Definition: ProtonReconstructionAlgorithm.h:59
ProtonReconstructionAlgorithm::fitter_
std::unique_ptr< ROOT::Fit::Fitter > fitter_
fitter object
Definition: ProtonReconstructionAlgorithm.h:83
edm::Ref
Definition: AssociativeIterator.h:58
LHCOpticalFunctionsSet::eLx
Definition: LHCOpticalFunctionsSet.h:16
ProtonReconstructionAlgorithm::reconstructFromSingleRP
reco::ForwardProton reconstructFromSingleRP(const CTPPSLocalTrackLiteRef &track, const LHCInfo &lhcInfo, std::ostream &os) const
run proton reconstruction using single-RP strategy
Definition: ProtonReconstructionAlgorithm.cc:403
parameters
parameters
Definition: BeamSpot_PayloadInspector.cc:14
LHCOpticalFunctionsSet::getXiValues
const std::vector< double > & getXiValues() const
Definition: LHCOpticalFunctionsSet.h:28
Vector
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:29
LHCInterpolatedOpticalFunctionsSetCollection
Definition: LHCInterpolatedOpticalFunctionsSetCollection.h:10
ProtonReconstructionAlgorithm::chiSquareCalculator_
std::unique_ptr< ChiSquareCalculator > chiSquareCalculator_
object to calculate chi^2
Definition: ProtonReconstructionAlgorithm.h:86
CovarianceMatrix
math::Error< 5 >::type CovarianceMatrix
Definition: SeedToTrackProducer.h:48
Point
math::XYZPoint Point
Definition: TrackerDpgAnalysis.cc:107
ProtonReconstructionAlgorithm::doLinearFit
static void doLinearFit(const std::vector< double > &vx, const std::vector< double > &vy, double &b, double &a)
Definition: ProtonReconstructionAlgorithm.cc:98
LHCInterpolatedOpticalFunctionsSet::Kinematics::y
double y
Definition: LHCInterpolatedOpticalFunctionsSet.h:31
ProtonReconstructionAlgorithm::mrpaNewton
Definition: ProtonReconstructionAlgorithm.h:53
LHCOpticalFunctionsSet::eLy
Definition: LHCOpticalFunctionsSet.h:16
Fit
Definition: Fit.h:32
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
LHCOpticalFunctionsSet::evy
Definition: LHCOpticalFunctionsSet.h:16
ProtonReconstructionAlgorithm::RPOpticsData::s_xi_vs_x_d
std::shared_ptr< const TSpline3 > s_xi_vs_x_d
Definition: ProtonReconstructionAlgorithm.h:59
ProtonReconstructionAlgorithm::mrpaChi2
Definition: ProtonReconstructionAlgorithm.h:53
ProtonReconstructionAlgorithm::RPOpticsData::s_v_y_vs_xi
std::shared_ptr< const TSpline3 > s_v_y_vs_xi
Definition: ProtonReconstructionAlgorithm.h:59
ProtonReconstructionAlgorithm::mrpaAnalIter
Definition: ProtonReconstructionAlgorithm.h:53
LHCInterpolatedOpticalFunctionsSet::Kinematics::x
double x
Definition: LHCInterpolatedOpticalFunctionsSet.h:29
ProtonReconstructionAlgorithm::ProtonReconstructionAlgorithm
ProtonReconstructionAlgorithm(bool fit_vtx_y, bool improved_estimate, const std::string &multiRPAlgorithm, unsigned int verbosity)
Definition: ProtonReconstructionAlgorithm.cc:21
ProtonReconstructionAlgorithm::multi_rp_algorithm_
enum ProtonReconstructionAlgorithm::@957 multi_rp_algorithm_
b
double b
Definition: hdecay.h:118
ProtonReconstructionAlgorithm::RPOpticsData::s_L_x_vs_xi
std::shared_ptr< const TSpline3 > s_L_x_vs_xi
Definition: ProtonReconstructionAlgorithm.h:59
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
runTauDisplay.gp
gp
Definition: runTauDisplay.py:431
bphysicsOniaDQM_cfi.vertex
vertex
Definition: bphysicsOniaDQM_cfi.py:7
a
double a
Definition: hdecay.h:119
ProtonReconstructionAlgorithm::verbosity_
unsigned int verbosity_
Definition: ProtonReconstructionAlgorithm.h:50
ProtonReconstructionAlgorithm::RPOpticsData::ch0
double ch0
intercept for linear approximation of
Definition: ProtonReconstructionAlgorithm.h:62
ForwardProton.h
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
ProtonReconstructionAlgorithm::RPOpticsData::s_L_y_vs_xi
std::shared_ptr< const TSpline3 > s_L_y_vs_xi
Definition: ProtonReconstructionAlgorithm.h:59
createfilelist.int
int
Definition: createfilelist.py:10
ProtonReconstructionAlgorithm::ChiSquareCalculator
class for calculation of chi^2
Definition: ProtonReconstructionAlgorithm.h:72
ProtonReconstructionAlgorithm::RPOpticsData::s_y_d_vs_xi
std::shared_ptr< const TSpline3 > s_y_d_vs_xi
Definition: ProtonReconstructionAlgorithm.h:59
LHCOpticalFunctionsSet::eyd
Definition: LHCOpticalFunctionsSet.h:16
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
ProtonReconstructionAlgorithm::mrpaUndefined
Definition: ProtonReconstructionAlgorithm.h:53
edm::RefVector::push_back
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
ProtonReconstructionAlgorithm::RPOpticsData
optics data associated with 1 RP
Definition: ProtonReconstructionAlgorithm.h:57
ProtonReconstructionAlgorithm::fitVtxY_
bool fitVtxY_
Definition: ProtonReconstructionAlgorithm.h:51
ProtonReconstructionAlgorithm::RPOpticsData::y0
double y0
beam vertical position, cm
Definition: ProtonReconstructionAlgorithm.h:61
funct::D
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
ProtonReconstructionAlgorithm::initialized_
bool initialized_
Definition: ProtonReconstructionAlgorithm.h:54
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
ProtonReconstructionAlgorithm::RPOpticsData::ch1
double ch1
slope for linear approximation of
Definition: ProtonReconstructionAlgorithm.h:63
HLT_2018_cff.maxIterations
maxIterations
Definition: HLT_2018_cff.py:11849
ProtonReconstructionAlgorithm::ChiSquareCalculator::m_rp_optics
const std::map< unsigned int, RPOpticsData > * m_rp_optics
Definition: ProtonReconstructionAlgorithm.h:79
Exception
Definition: hltDiff.cc:246
CTPPSDetId.h
ProtonReconstructionAlgorithm::m_rp_optics_
std::map< unsigned int, RPOpticsData > m_rp_optics_
map: RP id --> optics data
Definition: ProtonReconstructionAlgorithm.h:69
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
ztail.d
d
Definition: ztail.py:151
mps_fire.result
result
Definition: mps_fire.py:303
cms::Exception
Definition: Exception.h:70
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LHCInterpolatedOpticalFunctionsSet::transport
void transport(const Kinematics &input, Kinematics &output, bool calculateAngles=false) const
transports proton according to the splines
Definition: LHCInterpolatedOpticalFunctionsSet.cc:17
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
ProtonReconstructionAlgorithm::release
void release()
Definition: ProtonReconstructionAlgorithm.cc:118
ProtonReconstructionAlgorithm::ChiSquareCalculator::tracks
const CTPPSLocalTrackLiteRefVector * tracks
Definition: ProtonReconstructionAlgorithm.h:78
ProtonReconstructionAlgorithm.h
ProtonReconstructionAlgorithm::useImprovedInitialEstimate_
bool useImprovedInitialEstimate_
Definition: ProtonReconstructionAlgorithm.h:52
ROOT
Definition: Transform3DPJ.h:35
OpticalFunctionsConfig_cfi.opticalFunctions
opticalFunctions
Definition: OpticalFunctionsConfig_cfi.py:16
g
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
LHCInterpolatedOpticalFunctionsSet::splines
const std::vector< std::shared_ptr< const TSpline3 > > & splines() const
Definition: LHCInterpolatedOpticalFunctionsSet.h:22