Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

refactor(gx2f): make function to update parameter #3840

Merged
merged 5 commits into from
Nov 14, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
42 changes: 17 additions & 25 deletions Core/include/Acts/TrackFitting/GlobalChiSquareFitter.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -694,6 +694,20 @@ std::size_t countMaterialStates(
return nMaterialSurfaces;
}

/// @brief Update parameters (and scattering angles if applicable)
///
/// @param params Parameters to be updated
/// @param deltaParamsExtended Delta parameters for bound parameter and scattering angles
/// @param nMaterialSurfaces Number of material surfaces in the track
/// @param scatteringMap Map of geometry identifiers to scattering properties,
/// containing all scattering angles and covariances
/// @param geoIdVector Vector of geometry identifiers corresponding to material surfaces
void updateGx2fParams(
BoundTrackParameters& params, const Eigen::VectorXd& deltaParamsExtended,
const std::size_t nMaterialSurfaces,
std::unordered_map<GeometryIdentifier, ScatteringProperties>& scatteringMap,
const std::vector<GeometryIdentifier>& geoIdVector);

/// @brief Calculate and update the covariance of the fitted parameters
///
/// This function calculates the covariance of the fitted parameters using
Expand All @@ -704,8 +718,6 @@ std::size_t countMaterialStates(
///
/// @param fullCovariancePredicted The covariance matrix to update
/// @param extendedSystem All parameters of the current equation system
///
/// @return deltaParams The calculated delta parameters.
void updateGx2fCovarianceParams(BoundMatrix& fullCovariancePredicted,
Gx2fSystem& extendedSystem);

Expand Down Expand Up @@ -1234,7 +1246,6 @@ class Gx2Fitter {
using PropagatorOptions = typename propagator_t::template Options<Actors>;

start_parameters_t params = sParameters;
BoundVector deltaParams = BoundVector::Zero();
double chi2sum = 0;
double oldChi2sum = std::numeric_limits<double>::max();

Expand Down Expand Up @@ -1272,10 +1283,6 @@ class Gx2Fitter {
for (nUpdate = 0; nUpdate < gx2fOptions.nUpdateMax; nUpdate++) {
ACTS_DEBUG("nUpdate = " << nUpdate + 1 << "/" << gx2fOptions.nUpdateMax);

// update params
params.parameters() += deltaParams;
ACTS_VERBOSE("Updated parameters: " << params.parameters().transpose());

// set up propagator and co
Acts::GeometryContext geoCtx = gx2fOptions.geoContext;
Acts::MagneticFieldContext magCtx = gx2fOptions.magFieldContext;
Expand Down Expand Up @@ -1392,14 +1399,10 @@ class Gx2Fitter {
extendedSystem.aMatrix().colPivHouseholderQr().solve(
extendedSystem.bVector());

deltaParams = deltaParamsExtended.topLeftCorner<eBoundSize, 1>().eval();

ACTS_VERBOSE("aMatrix:\n"
<< extendedSystem.aMatrix() << "\n"
<< "bVector:\n"
<< extendedSystem.bVector() << "\n"
<< "deltaParams:\n"
<< deltaParams << "\n"
<< "deltaParamsExtended:\n"
<< deltaParamsExtended << "\n"
<< "oldChi2sum = " << oldChi2sum << "\n"
Expand Down Expand Up @@ -1440,19 +1443,9 @@ class Gx2Fitter {
break;
}

if (multipleScattering) {
// update the scattering angles
for (std::size_t matSurface = 0; matSurface < nMaterialSurfaces;
matSurface++) {
const std::size_t deltaPosition = eBoundSize + 2 * matSurface;
const GeometryIdentifier geoId = geoIdVector[matSurface];
const auto scatteringMapId = scatteringMap.find(geoId);
assert(scatteringMapId != scatteringMap.end() &&
"No scattering angles found for material surface.");
scatteringMapId->second.scatteringAngles().block<2, 1>(2, 0) +=
deltaParamsExtended.block<2, 1>(deltaPosition, 0).eval();
}
}
updateGx2fParams(params, deltaParamsExtended, nMaterialSurfaces,
scatteringMap, geoIdVector);
ACTS_VERBOSE("Updated parameters: " << params.parameters().transpose());
AJPfleger marked this conversation as resolved.
Show resolved Hide resolved

oldChi2sum = extendedSystem.chi2();
}
Expand All @@ -1478,7 +1471,6 @@ class Gx2Fitter {
// step, we will not ignore the boundary checks for measurement surfaces. We
// want to create trackstates only on surfaces, that we actually hit.
if (gx2fOptions.nUpdateMax > 0) {
ACTS_VERBOSE("Final delta parameters: " << deltaParams.transpose());
ACTS_VERBOSE("Propagate with the final covariance.");
// update covariance
params.covariance() = fullCovariancePredicted;
Expand Down
24 changes: 24 additions & 0 deletions Core/src/TrackFitting/GlobalChiSquareFitter.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,30 @@

#include "Acts/Definitions/TrackParametrization.hpp"

void Acts::Experimental::updateGx2fParams(
BoundTrackParameters& params, const Eigen::VectorXd& deltaParamsExtended,
const std::size_t nMaterialSurfaces,
std::unordered_map<GeometryIdentifier, ScatteringProperties>& scatteringMap,
const std::vector<GeometryIdentifier>& geoIdVector) {
// update params
params.parameters() +=
deltaParamsExtended.topLeftCorner<eBoundSize, 1>().eval();

// update the scattering angles.
for (std::size_t matSurface = 0; matSurface < nMaterialSurfaces;
matSurface++) {
const std::size_t deltaPosition = eBoundSize + 2 * matSurface;
const GeometryIdentifier geoId = geoIdVector[matSurface];
const auto scatteringMapId = scatteringMap.find(geoId);
assert(scatteringMapId != scatteringMap.end() &&
"No scattering angles found for material surface.");
scatteringMapId->second.scatteringAngles().block<2, 1>(2, 0) +=
deltaParamsExtended.block<2, 1>(deltaPosition, 0).eval();
}

return;
AJPfleger marked this conversation as resolved.
Show resolved Hide resolved
}

void Acts::Experimental::updateGx2fCovarianceParams(
BoundMatrix& fullCovariancePredicted, Gx2fSystem& extendedSystem) {
// make invertible
Expand Down
Loading