15 #include "TMinuitMinimizer.h" 23 bool improved_estimate,
28 useImprovedInitialEstimate_(improved_estimate),
29 multi_rp_algorithm_(mrpaUndefined),
34 TMinuitMinimizer::UseStaticMinuit(
false);
48 double pStart[] = {0, 0, 0, 0};
71 vector<double> xiValues =
74 rpod.
s_xi_vs_x_d = make_shared<TSpline3>(
"", xDValues.data(), xiValues.data(), xiValues.size());
100 const std::vector<double> &
vy,
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) {
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;
141 oit->second.optics->transport(k_in, k_out);
144 const double x = k_out.
x - oit->second.x0;
145 const double y = k_out.
y - oit->second.y0;
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);
152 s2 += x_diff_norm * x_diff_norm + y_diff_norm * y_diff_norm;
173 std::ostream &
os)
const {
179 <<
"Optics data not available for RP " <<
it->rpId() <<
", i.e. " <<
CTPPSDetId(
it->rpId()) <<
".";
183 double xi_init = 0., th_x_init = 0.;
186 double x_N =
tracks[0]->x() * 1E-1,
187 x_F =
tracks[1]->x() * 1E-1;
192 const double a = i_F.ch1 * i_N.
la1 - i_N.
ch1 * i_F.la1;
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.;
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);
202 double s_xi0 = 0., s_1 = 0.;
205 double xi = oit->second.s_xi_vs_x_d->Eval(
track->x() * 1E-1 + oit->second.x0);
211 xi_init = s_xi0 / s_1;
214 if (!std::isfinite(xi_init))
216 if (!std::isfinite(th_x_init))
220 double y[2] = {0}, v_y[2] = {0}, L_y[2] = {0};
221 unsigned int y_idx = 0;
228 y[y_idx] =
track->y() * 1E-1 - oit->second.s_y_d_vs_xi->Eval(xi_init);
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);
235 double vtx_y_init = 0.;
236 double th_y_init = 0.;
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.;
244 th_y_init = (
y[1] / L_y[1] +
y[0] / L_y[0]) / 2.;
247 if (!std::isfinite(vtx_y_init))
249 if (!std::isfinite(th_y_init))
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;
262 double xi = 0., th_x = 0., th_y = 0., vtx_y = 0.;
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);
275 fitter_->Config().ParSettings(3).Fix();
294 map<unsigned int, signed int> index_map = {
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},
302 for (
unsigned int i = 0;
i < (
unsigned int)FP::Index::num_indices; ++
i) {
303 signed int fit_i = index_map[
i];
305 for (
unsigned int j = 0;
j < (
unsigned int)FP::Index::num_indices; ++
j) {
306 signed int fit_j = index_map[
j];
308 cm(
i,
j) = (fit_i >= 0 && fit_j >= 0) ?
result.CovMatrix(fit_i, fit_j) : 0.;
316 const double maxXiDiff = 1E-6;
317 const double xi_ep = 1E-5;
320 double x_N =
tracks[0]->x() * 1E-1,
321 x_F =
tracks[1]->x() * 1E-1;
323 double y_N =
tracks[0]->y() * 1E-1,
324 y_F =
tracks[1]->y() * 1E-1;
331 double xi_prev = xi_init;
335 const double gp = (
newtonGoalFcn(xi_prev + xi_ep, x_N, x_F, i_N, i_F) -
g) / xi_ep;
337 xi = xi_prev -
g /
gp;
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;
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);
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);
350 if (
abs(
xi - xi_prev) < maxXiDiff) {
358 th_x = (x_F - i_F.s_x_d_vs_xi->Eval(
xi)) / i_F.s_L_x_vs_xi->Eval(
xi);
362 const double y_eff_F = y_F - i_F.s_y_d_vs_xi->Eval(
xi);
365 const double l_y_F = i_F.s_L_y_vs_xi->Eval(
xi);
368 const double v_y_F = i_F.s_v_y_vs_xi->Eval(
xi);
370 const double det = l_y_N * v_y_F - l_y_F * v_y_N;
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;
376 th_y = (y_eff_N / l_y_N + y_eff_F / l_y_F) / 2.;
385 os <<
" fit valid=" <<
valid << std::endl
386 <<
" xi=" <<
xi <<
", th_x=" << th_x <<
", th_y=" << th_y <<
", vtx_y=" << vtx_y <<
", chiSq = " <<
chi2 390 const double sign_z = (armId == 0) ? +1. : -1.;
392 const double cos_th_sq = 1. - th_x * th_x - th_y * th_y;
393 const double cos_th = (cos_th_sq > 0.) ?
sqrt(cos_th_sq) : 1.;
397 sign_z *
p * cos_th);
407 std::ostream &
os)
const {
411 os <<
"reconstructFromSingleRP(" << rpId.arm() * 100 + rpId.station() * 10 + rpId.rp() <<
")" << std::endl;
417 <<
"Optics data not available for RP " <<
track->rpId() <<
", i.e. " << rpId <<
".";
420 const double x_full =
track->x() * 1E-1 + oit->second.x0;
421 const double xi = oit->second.s_xi_vs_x_d->Eval(x_full);
422 const double L_y = oit->second.s_L_y_vs_xi->Eval(
xi);
423 const double th_y =
track->y() * 1E-1 / L_y;
425 const double ep_x = 1E-6;
426 const double dxi_dx = (oit->second.s_xi_vs_x_d->Eval(x_full + ep_x) -
xi) / ep_x;
427 const double xi_unc =
abs(dxi_dx) *
track->xUnc() * 1E-1;
429 const double ep_xi = 1E-4;
430 const double dL_y_dxi = (oit->second.s_L_y_vs_xi->Eval(
xi + ep_xi) - L_y) / ep_xi;
431 const double th_y_unc_sq = th_y * th_y * (
pow(
track->yUnc() /
track->y(), 2.) +
pow(dL_y_dxi * xi_unc / L_y, 2.));
434 os <<
" xi = " <<
xi <<
" +- " << xi_unc <<
", th_y = " << th_y <<
" +- " <<
sqrt(th_y_unc_sq) <<
"." 442 const double cos_th_sq = 1. - th_y * th_y;
443 const double cos_th = (cos_th_sq > 0.) ?
sqrt(cos_th_sq) : 1.;
445 const FP::Vector momentum(0.,
p * th_y, sign_z *
p * cos_th);
449 cm((
int)FP::Index::th_y, (
int)FP::Index::th_y) = th_y_unc_sq;
std::unique_ptr< ROOT::Fit::Fitter > fitter_
fitter object
std::shared_ptr< const TSpline3 > s_xi_vs_x_d
optics data associated with 1 RP
class for calculation of chi^2
reco::ForwardProton reconstructFromSingleRP(const CTPPSLocalTrackLiteRef &track, const float energy, std::ostream &os) const
run proton reconstruction using single-RP strategy
double operator()(const double *parameters) const
std::map< unsigned int, RPOpticsData > m_rp_optics_
map: RP id –> optics data
double y0
beam vertical position, cm
std::shared_ptr< const TSpline3 > s_x_d_vs_xi
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
const std::map< unsigned int, RPOpticsData > * m_rp_optics
void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions)
std::unique_ptr< ChiSquareCalculator > chiSquareCalculator_
object to calculate chi^2
bool useImprovedInitialEstimate_
double ch0
intercept for linear approximation of
ProtonReconstructionAlgorithm(bool fit_vtx_y, bool improved_estimate, const std::string &multiRPAlgorithm, unsigned int verbosity)
double x0
beam horizontal position, cm
static double newtonGoalFcn(double xi, double x_N, double x_F, const RPOpticsData &i_N, const RPOpticsData &i_F)
double la1
slope for linear approximation of
static void doLinearFit(const std::vector< double > &vx, const std::vector< double > &vy, double &b, double &a)
Abs< T >::type abs(const T &t)
void transport(const Kinematics &input, Kinematics &output, bool calculateAngles=false) const
transports proton according to the splines
const std::vector< std::shared_ptr< const TSpline3 > > & splines() const
enum ProtonReconstructionAlgorithm::@946 multi_rp_algorithm_
reco::ForwardProton reconstructFromMultiRP(const CTPPSLocalTrackLiteRefVector &tracks, const float energy, std::ostream &os) const
run proton reconstruction using multiple-RP strategy
proton kinematics description
const CTPPSLocalTrackLiteRefVector * tracks
double la0
intercept for linear approximation of
const std::vector< std::vector< double > > & getFcnValues() const
std::shared_ptr< const TSpline3 > s_v_y_vs_xi
DecomposeProduct< arg, typename Div::arg > D
Set of optical functions corresponding to one scoring plane along LHC, including splines for interpol...
const std::vector< double > & getXiValues() const
Base class for CTPPS detector IDs.
std::shared_ptr< const TSpline3 > s_L_y_vs_xi
Eigen::Matrix< float, EcalPulseShape::TEMPLATESAMPLES, EcalPulseShape::TEMPLATESAMPLES > CovarianceMatrix
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
const LHCInterpolatedOpticalFunctionsSet * optics
double ch1
slope for linear approximation of
std::shared_ptr< const TSpline3 > s_L_x_vs_xi
std::shared_ptr< const TSpline3 > s_y_d_vs_xi
Power< A, B >::type pow(const A &a, const B &b)