35 #define MAX_ITERATIONS_NEWTON 10
37 #define MAX_ITERATIONS_RUNGEKUTTA 5000
96 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.order",
97 po::value<int>(0)->default_value(5),
98 "Runge-Kutta method order (higher than 4 implies adaptive step)");
100 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.step",
101 po::value<float>(0)->default_value(1.0),
102 "Runge-Kutta fixed step size (for order 4 or less) [cm]");
104 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.min_step",
105 po::value<float>(0)->default_value(0.1),
106 "Runge-Kutta minimum step size (for adaptive step) [cm]");
108 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.max_step",
109 po::value<float>(0)->default_value(10.0),
110 "Runge-Kutta maximum step size (for adaptive step) [cm]");
112 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.tolerance",
113 po::value<float>(0)->default_value(1
e-10),
114 "Runge-Kutta adaptive step maximum allowable truncation error");
116 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.momentum_step",
117 po::value<float>(0)->default_value(30.0),
118 "Newton's method momentum step [MeV]");
123 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.initial_momemtum",
124 po::value<float>(0)->default_value(1.25),
125 "Newton's method position step [GeV]");
128 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.start_position",
129 po::value<float>(0)->default_value(-250.0),
130 "Magnetic field swimming start position [cm]");
132 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.end_position",
133 po::value<float>(0)->default_value(+577.0),
134 "Magnetic field swimming end position [cm]");
137 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.optim_variable",
138 po::value<int>(0)->default_value(0),
139 "Newton's method optimization variable: 0 for position.Perp(), 1 for direction.Theta()");
141 options.
AddOptions(
"Momentum reconstruction")(
"QwRayTracer.optim_resolution",
142 po::value<float>(0)->default_value(0.3),
143 "Newton's method optimization cut-off resolution [cm for position, rad for angles]");
165 switch (fOptimizationVariable) {
166 case 0: fOptimizationResolution *=
Qw::cm;
167 case 1: fOptimizationResolution *=
Qw::rad;
168 default: fOptimizationResolution *=
Qw::cm;
206 TVector3 position = start_position;
207 TVector3 direction = start_direction;
214 int iterations_newton = 0;
215 int iterations_rungekutta = 0;
225 Double_t r[2] = {0.0, 0.0};
228 position = start_position;
229 direction = start_direction;
234 position = start_position;
235 direction = start_direction;
245 position = start_position;
246 direction = start_direction;
281 track =
new QwTrack(front,back);
289 position = start_position;
290 direction = start_direction;
296 position = start_position;
297 direction = start_direction;
310 TVector3 position_diff = position - end_position;
314 track->
fPositionRoff = position.Perp() - end_position.Perp();
319 TVector3 direction_diff = direction - end_direction;
355 const double alpha = 0.5;
357 TMatrixD A(s,s); A.Zero();
360 TMatrixD b(q,s); b.Zero();
361 b[0][0] = 1.0 - 1.0 / (2.0 * alpha);
362 b[0][1] = 1.0 / (2.0 * alpha);
364 return IntegrateRK(A, b, r, v, p, z, h, tolerance);
373 TMatrixD A(s,s); A.Zero();
374 A[1][0] = A[2][1] = 1.0 / 2.0;
377 TMatrixD b(q,s); b.Zero();
378 b[0][0] = b[0][3] = 1.0 / 6.0;
379 b[0][1] = b[0][2] = 1.0 / 3.0;
381 return IntegrateRK(A, b, r, v, p, z, h, tolerance);
390 TMatrixD A(s,s); A.Zero();
392 A[2][0] = 3.0 / 32.0; A[2][1] = 9.0 / 32.0;
393 A[3][0] = 1932.0 / 2197.0; A[3][1] = -7200.0 / 2197.0; A[3][2] = 7296.0 / 2197.0;
394 A[4][0] = 439.0 / 216.0; A[4][1] = -8.0; A[4][2] = 3680.0 / 513.0; A[4][3] = -845.0 / 4104.0;
395 A[5][0] = -8.0 / 27.0; A[5][1] = 2.0; A[5][2] = -3544.0 / 2565.0; A[5][3] = 1859.0 / 4104.0; A[5][4] = -11.0 / 40.0;
397 TMatrixD b(q,s); b.Zero();
398 b[0][0] = 16.0 / 135.0;
399 b[0][2] = 6656.0 / 12825.0;
400 b[0][3] = 28561.0 / 56430.0;
401 b[0][4] = -9.0 / 50.0;
402 b[0][5] = 2.0 / 55.0;
403 b[1][0] = 25.0 / 216.0;
404 b[1][2] = 1408.0 / 2565.0;
405 b[1][3] = 2197.0 / 4104.0;
406 b[1][4] = -1.0 / 5.0;
408 return IntegrateRK(A, b, r, v, p, z, h, tolerance);
416 TMatrixD A(s,s); A.Zero();
418 A[2][0] = 4.0 / 75.0; A[2][1] = 16.0 / 75.0;
419 A[3][0] = 5.0 / 6.0; A[3][1] = -8.0 / 3.0; A[3][2] = 5.0 / 2.0;
420 A[4][0] = -8.0 / 5.0; A[4][1] = 144.0 / 25.0; A[4][2] = -4.0; A[4][3] = 16.0 / 25.0;
421 A[5][0] = 361.0 / 320.0; A[5][1] = -18.0 / 5.0; A[5][2] = 407.0 / 128.0; A[5][3] = -11.0 / 80.0; A[5][4] = 55.0 / 128.0;
422 A[6][0] = -11.0 / 640.0; A[6][1] = 0.0; A[6][2] = 11.0 / 256.0; A[6][3] = -11.0 / 160.0; A[6][4] = 11.0 / 256.0; A[6][5] = 0.0;
423 A[7][0] = 93.0 / 640.0; A[7][1] = -18.0 / 5.0; A[7][2] = 803.0 / 256.0; A[7][3] = -11.0 / 160.0; A[7][4] = 99.0 / 256.0; A[7][5] = 0.0; A[7][6] = 1.0;
425 TMatrixD b(q,s); b.Zero();
426 b[0][0] = 31.0 / 384.0;
427 b[0][2] = 1125.0 / 2816.0;
428 b[0][3] = 9.0 / 32.0;
429 b[0][4] = 125.0 / 768.0;
430 b[0][5] = 5.0 / 66.0;
431 b[1][0] = 7.0 / 1408.0;
432 b[1][2] = 1125.0 / 2816.0;
433 b[1][3] = 9.0 / 32.0;
434 b[1][4] = 125.0 / 768.0;
435 b[1][6] = 5.0 / 66.0;
436 b[1][7] = 5.0 / 66.0;
444 TMatrixD A(s,s); A.Zero();
445 A[1][0] = 2.0 / 27.0;
446 A[2][0] = 1.0 / 36.0; A[2][1] = 1.0 / 12.0;
447 A[3][0] = 1.0 / 24.0; A[3][1] = 0.0; A[3][2] = 1.0 / 8.0;
448 A[4][0] = 5.0 / 12.0; A[4][1] = 0.0; A[4][2] = -25.0 / 16.0; A[4][3] = 25.0 / 16.0;
449 A[5][0] = 1.0 / 20.0; A[5][1] = 0.0; A[5][2] = 0.0; A[5][3] = 1.0 / 4.0; A[5][4] = 1.0 / 5.0;
450 A[6][0] = -25.0 / 108.0; A[6][1] = 0.0; A[6][2] = 0.0; A[6][3] = 125.0 / 108.0; A[6][4] = -65.0 / 27.0; A[6][5] = 125.0 / 54.0;
451 A[7][0] = 31.0 / 300.0; A[7][1] = 0.0; A[7][2] = 0.0; A[7][3] = 0.0; A[7][4] = 61.0 / 225.0; A[7][5] = -2.0 / 9.0; A[7][6] = 13.0 / 900.0;
452 A[8][0] = 2.0; A[8][1] = 0.0; A[8][2] = 0.0; A[8][3] = -53.0 / 6.0; A[8][4] = 704.0 / 45.0; A[8][5] = -107.0 / 9.0; A[8][6] = 67.0 / 90.0; A[8][7] = 3.0;
453 A[9][0] = -91.0 / 108.0; A[9][1] = 0.0; A[9][2] = 0.0; A[9][3] = 23.0 / 108.0; A[9][4] = -976.0 / 135.0; A[9][5] = 311.0 / 54.0; A[9][6] = -19.0 / 60.0; A[9][7] = 17.0 / 6.0; A[9][9] = -1.0 / 12.0;
454 A[10][0] = 2383.0 / 4100.0; A[10][1] = 0.0; A[10][2] = 0.0; A[10][3] = -341.0 / 164.0; A[10][4] = 4496.0 / 1025.0; A[10][5] = -301.0 / 82.0; A[10][6] = 2133.0 / 4100.0; A[10][7] = 45.0 / 82.0; A[10][8] = 45.0 / 164.0; A[10][9] = 18.0 / 41.0;
455 A[11][0] = 3.0 / 205.0; A[11][1] = 0.0; A[11][2] = 0.0; A[11][3] = 0.0; A[11][4] = 0.0; A[11][5] = -6.0 / 41.0; A[11][6] = -3.0 / 205.0; A[11][7] = -3.0 / 41.0; A[11][8] = 3.0 / 41.0; A[11][9] = 6.0 / 41.0; A[11][10] = 0.0;
456 A[12][0] = -1777.0 / 4100.0; A[12][1] = 0.0; A[12][2] = 0.0; A[12][3] = -341.0 / 164.0; A[12][4] = 4496.0 / 1025.0; A[12][5] = -289.0 / 82.0; A[12][6] = 2193.0 / 4100.0; A[12][7] = 51.0 / 82.0; A[12][8] = 33.0 / 164.0; A[12][9] = 12.0 / 41.0; A[12][10] = 0.0; A[12][11] = 1.0;
458 TMatrixD b(q,s); b.Zero();
459 b[0][0] = 41.0 / 840.0;
460 b[0][5] = 34.0 / 105.0;
461 b[0][6] = 9.0 / 35.0;
462 b[0][7] = 9.0 / 35.0;
463 b[0][8] = 9.0 / 280.0;
464 b[0][9] = 9.0 / 280.0;
465 b[0][10] = 41.0 / 840.0;
466 b[1][5] = 34.0 / 105.0;
467 b[1][6] = 9.0 / 35.0;
468 b[1][7] = 9.0 / 35.0;
469 b[1][8] = 9.0 / 280.0;
470 b[1][9] = 9.0 / 280.0;
471 b[1][11] = 41.0 / 840.0;
472 b[1][12] = 41.0 / 840.0;
518 const double tolerance)
521 const int s = A.GetNrows();
524 const int q = b.GetNrows();
531 for (
int i = 0; i < s; i++)
532 for (
int j = 0; j < s; j++)
539 TVector3 r_new[2], r_old[2];
540 TVector3 v_new[2], v_old[2];
541 r_old[0] = r_old[1] = r;
542 v_old[0] = v_old[1] = v;
546 bool last_iteration =
false;
551 if (fabs(r_old[0].Z() - z) < h) {
552 TVector3 r_diff = (r_old[0].Z() - z) / v_old[0].Z() * v_old[0];
554 last_iteration =
true;
561 for (
int i = 0; i < s; i++) {
564 for (
int j = 0; j < i; j++) {
565 r_new[0] += A[i][j] * k_r[j];
566 v_new[0] += A[i][j] * k_v[j];
570 k_r[i] = h * v_new[0];
571 k_v[i] = beta * k_r[i].Cross(B);
575 for (
int p = 0; p < q; p++) {
578 for (
int i = 0; i < s; i++) {
579 r_new[p] += b[p][i] * k_r[i];
580 v_new[p] += b[p][i] * k_v[i];
589 TVector3 r_diff = r_new[0] - r_new[1];
590 h *= pow(tolerance * h / (2 * r_diff.Mag()), 0.25);
591 if (h < h_min) h = h_min;
592 if (h > h_max) h = h_max;
596 r_old[0] = r_old[1] = r_new[0];
597 v_old[0] = v_old[1] = v_new[0];
606 else return iterations;
double fDirectionYoff
Difference in Y momentum at matching plane.
bool LoadMagneticFieldMap(QwOptions &options)
Load the magnetic field based on config file options.
#define QwMessage
Predefined log drain for regular messages.
double fChi
Combined chi, i.e. the sum of the chi for front and back partial track.
static QwMagneticField * fBfield
Magnetic field (static)
static void DefineOptions(QwOptions &options)
Define command line and config file options.
int fIterationsRungeKutta
Number of iterations in Runge-Kutta method.
Magnetic field map object.
double fDirectionThetaoff
Difference in momentum polar angle at matching plane.
Definition of the ray-tracing bridging method for R2/R3 partial tracks.
double fOptimizationResolution
Double_t fBdly
y component of the field integral
double fPositionPhioff
Difference in azimuthal angle phi at matching plane.
int fIterationsRK4
Number of iterations using Runge-Kutta 4th order.
int fIterationsNewton
Number of iterations in Newton's method.
double fStartPosition
Starting position for magnetic field swimming.
double fIntegrationStep
Runge-Kutta fixed step size for order 4 or less.
int fOptimizationVariable
Newton's method optimization variable.
TVector3 fStartPosition
Position and direction before and after swimming.
int GetOctant() const
Get the octant number.
TVector3 fEndPositionActualRK4
Actual position of track at back plane using Runge-Kutta 4th order.
#define QwVerbose
Predefined log drain for verbose messages.
TVector3 fDirectionDiff
Difference in momentum at matching plane.
Double_t fChi
combined chi square
#define MAX_ITERATIONS_NEWTON
void SetPackage(EQwDetectorPackage package)
Set the package.
const TVector3 GetPosition(const double z) const
Return the vertex at position z.
Double_t fBdl
scalar field integrals
double fDirectionPhioff
Difference in momentum azimuthal angle at matching plane.
Definition of the track class.
double fPositionRoff
Difference in radial position at matching plane.
po::options_description_easy_init AddOptions(const std::string &blockname="Specialized options")
Add an option to a named block or create new block.
void SetMagneticFieldIntegral(double bdl, double bdlx, double bdly, double bdlz)
Set magnetic field integral.
TVector3 fEndDirectionActual
Actual direction of track at back plane.
TVector3 fStartDirection
Start direction of front track.
T GetValue(const std::string &key)
Get a templated value.
void GetCartesianFieldValue(const double point_xyz[3], double field_xyz[3]) const
Get the cartesian components of the field value.
A logfile class, based on an identical class in the Hermes analyzer.
const TVector3 GetMomentumDirection() const
Return the direction.
QwRayTracer(QwOptions &options)
Default constructor.
double fIntegrationMaximumStep
Maximum step size.
Interface to the various bridging methods.
TVector3 fEndPositionActual
Actual position of track at back plane.
Contains the complete track as a concatenation of partial tracks.
void ProcessOptions(QwOptions &options)
Process command line and config file options.
double fDirectionXoff
Difference in X momentum at matching plane.
double fInitialMomentum
Newton's method initial momentum.
void SetOctant(int octant)
Set the octant number.
#define MAX_ITERATIONS_RUNGEKUTTA
double fDirectionZoff
Difference in Z momentum at matching plane.
TVector3 fEndPositionActualRKF45
Actual position of track at back plane using Runge-Kutta-Fehlberg.
const QwTrack * Bridge(const QwPartialTrack *front, const QwPartialTrack *back)
Bridge from the front to back partial track.
TVector3 fPositionDiff
Matching of front and back track position and direction at matching plane.
Double_t fBdlx
x component of the field integral
TVector3 fEndDirectionGoal
Goal direction of back track.
int fIterationsRKF45
Number of iterations using Runge-Kutta-Fehlberg.
static std::ostream & endl(std::ostream &)
End of the line.
TVector3 fEndPositionGoal
Goal position of back track.
double fPositionXoff
Difference in X position at matching plane.
int fIntegrationOrder
Runge-Kutta method order.
double fEndPosition
Ending position for magnetic field swimming.
TVector3 fEndDirectionActualRK4
Actual direction of track at back plane using Runge-Kutta 4th order.
double fPositionYoff
Difference in Y position at matching plane.
double fIntegrationMinimumStep
Runge-Kutta minimum and maximum step size for adaptive step.
TVector3 fEndDirectionActualRKF45
Actual direction of track at back plane using Runge-Kutta-Fehlberg.
double fPositionThetaoff
Difference in polar angle theta at matching plane.
double GetOptimizationVariable(const TVector3 &position, const TVector3 &direction) const
Optimization variable from position and direction.
double fIntegrationTolerance
Allowable truncation error in adaptive step.
double fMomentum
Spectrometer momentum.
bool TestFieldMap()
Test the field map.
double fMomentumStep
Newton's method step size in momentum.
Contains the straight part of a track in one region only.
int IntegrateRK(const TMatrixD &A, const TMatrixD &b, TVector3 &r, TVector3 &v, const double p, const double z, const double h, const double epsilon)
Runge-Kutta numerical integration by Butcher tableau.
Double_t fBdlz
z component of the field integral
static const double c
Physical constants.
EQwDetectorPackage GetPackage() const
Get the package.
static const double cm
Length units: base unit is mm.
virtual ~QwRayTracer()
Destructor.