ATLAS Offline Software
Loading...
Searching...
No Matches
HLT::MET::PufitUtils Namespace Reference

Classes

struct  CovarianceSum
 Helper struct to hold the sum over pileup objects and its covariance. More...

Functions

void trimmedMeanAndVariance (const std::vector< double > &sorted, double trimFraction, double &mean, double &variance)
 Calculate the trimmed mean and variance for a vector of tower sumEts.
void trimmedMeanAndVariance (const PufitGrid &grid, double trimFraction, double &mean, double &variance)
 Calculate the trimmed mean and variance for a vector of tower sumEts.
void unmaskedMeanAndVariance (const PufitGrid &grid, double &mean, double &variance)
 Calculate the mean and variance of unmasked towers.
GridDisplacement selectGrid (const PufitGridSet &grids)
 Select the grid with the highest masked sumEt.
Eigen::VectorXd pufit (const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, const Eigen::VectorXd &towerExpectations, const Eigen::VectorXd &towerVariances, const Eigen::VectorXd &correctionDirections, double constraintImportance=1)
 Perform the pile-up fit.
Eigen::VectorXd pufit (const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, double towerMean, double towerVariance, const Eigen::VectorXd &correctionDirections, double constraintImportance=1)
 Perform the pile-up fit.
Eigen::VectorXd pufit (const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, const Eigen::VectorXd &towerExpectations, const Eigen::VectorXd &towerVariances, const Eigen::Matrix< double, 2, Eigen::Dynamic > &cosSin, double constraintImportance=1)
 Perform the pile-up fit.
Eigen::VectorXd pufit (const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, double towerMean, double towerVariance, const Eigen::Matrix< double, 2, Eigen::Dynamic > &cosSin, double constraintImportance=1)
 Perform the pile-up fit.
std::vector< SignedKinematicspufit (const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, const std::vector< double > &towerExpectations, const std::vector< double > &towerVariances, const std::vector< SignedKinematics > &toCorrect, double constraintImportance=1)
 Perform the pile-up fit.
std::vector< SignedKinematicspufit (const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, double towerMean, double towerVariance, const std::vector< SignedKinematics > &toCorrect, double constraintImportance=1)
 Perform the pile-up fit.
template<std::size_t N>
void trimmedMeanAndVariance (const PufitMultiGrid< N > &grid, std::size_t type, double trimFraction, double &mean, double &variance)
 Calculate the trimmed mean and variance for a vector of tower sumEts.
template<std::size_t N>
void unmaskedMeanAndVariance (const PufitMultiGrid< N > &grid, int type, double &mean, double &variance)
 Calculate the mean and variance of unmasked towers.
template<typename Grid>
GridDisplacement selectGrid (const PufitMultiGridSet< Grid > &grids, std::size_t type)
 Select the grid with the highest masked sumEt.

Function Documentation

◆ pufit() [1/6]

Eigen::VectorXd HLT::MET::PufitUtils::pufit ( const Eigen::Vector2d & pileupSum,
const Eigen::Matrix2d & pileupCovariance,
const Eigen::VectorXd & towerExpectations,
const Eigen::VectorXd & towerVariances,
const Eigen::Matrix< double, 2, Eigen::Dynamic > & cosSin,
double constraintImportance = 1 )

Perform the pile-up fit.

Parameters
pileupSum2-vector sum over the pileup contributions
pileupCovarianceCovariance matrix of the pileup sum
towerExpectationsExpected values for each tower based on the pileup energy density
towerVariancesExpected values for the variances of each tower based on variances of background towers
cosSinThe cosines and sines of the objects to correct
Therelative importance of the constraints.
Returns
A vector containing the magnitudes of the corrections

This variant directly provides the cos-sin matrix for the correction, potentially allowing many cos/sin calls to be avoided. This variant provides individual estimates of the expected energy in each tower and the corresponding variances - these are for cases where the areas of each object are different.

Definition at line 163 of file PufitUtils.cxx.

170 {
171 // Keep track of the number of corrections to derive
172 std::size_t nCorr = cosSin.cols();
173 if (nCorr == 0)
174 // No need to do anything to derive corrections to nothing
175 return {};
176 // Eigen internally uses somewhat complicated 'expression templates' to hold
177 // the mathematics of a matrix calculation without actually doing the loop.
178 // This allows it to avoid doing more loops than it actually needs to. It
179 // will do this by default unless you force it to evaluate an expression by
180 // binding it to a vector matrix type. Therefore we use 'auto' to store the
181 // individual expression pieces, rather than forcing an early evaluation.
182 if (pileupCovariance.determinant() == 0)
183 {
184 // In very rare circumstances the pileup covariance can be 0. This essentially means that
185 // there is nothing interpreted as pileup. The calculation will fail in this case so we
186 // cannot run the fit. Just return the expectations
187 return towerExpectations;
188 }
189 auto constants =
190 // Constant part derived from the uniformity constraint
191 towerExpectations.cwiseQuotient(towerVariances) -
192 // Constant part derived from the pileup balance constraint
193 constraintImportance * cosSin.transpose() * pileupCovariance.inverse() * pileupSum;
194 Eigen::MatrixXd diagonal = towerVariances.cwiseInverse().asDiagonal();
195 auto coeffMatrix =
196 // Matrix part derived from the uniformity constraint
197 diagonal +
198 // Matrix part derived from the pileup balance constraint
199 constraintImportance * (cosSin.transpose() * pileupCovariance.inverse() * cosSin);
200 // Now return the actual corrections
201 return coeffMatrix.inverse() * constants;
202 }

◆ pufit() [2/6]

Eigen::VectorXd HLT::MET::PufitUtils::pufit ( const Eigen::Vector2d & pileupSum,
const Eigen::Matrix2d & pileupCovariance,
const Eigen::VectorXd & towerExpectations,
const Eigen::VectorXd & towerVariances,
const Eigen::VectorXd & correctionDirections,
double constraintImportance = 1 )

Perform the pile-up fit.

Parameters
pileupSum2-vector sum over the pileup contributions
pileupCovarianceCovariance matrix of the pileup sum
towerExpectationsExpected values for each tower based on the pileup energy density
towerVariancesExpected values for the variances of each tower based on variances of background towers
correctionDirectionsThe directions of the objects to correct
Therelative importance of the constraints.
Returns
A vector containing the magnitudes of the corrections

This variant provides individual estimates of the expected energy in each tower and the corresponding variances - these are for cases where the areas of each object are different.

Definition at line 119 of file PufitUtils.cxx.

126 {
127 // Keep track of the number of corrections to derive
128 std::size_t nCorr = correctionDirections.size();
129 if (nCorr == 0)
130 // No need to do anything to derive corrections to nothing
131 return {};
132 // Turn the directions into cos/sin pairs
133 Eigen::Matrix<double, 2, Eigen::Dynamic> cosSin(2, nCorr);
134 cosSin.row(0) = correctionDirections.array().cos();
135 cosSin.row(1) = correctionDirections.array().sin();
136 return pufit(
137 pileupSum,
138 pileupCovariance,
139 towerExpectations,
140 towerVariances,
141 cosSin,
142 constraintImportance);
143 }
Eigen::VectorXd pufit(const Eigen::Vector2d &pileupSum, const Eigen::Matrix2d &pileupCovariance, const Eigen::VectorXd &towerExpectations, const Eigen::VectorXd &towerVariances, const Eigen::VectorXd &correctionDirections, double constraintImportance)
Perform the pile-up fit.

◆ pufit() [3/6]

std::vector< SignedKinematics > HLT::MET::PufitUtils::pufit ( const Eigen::Vector2d & pileupSum,
const Eigen::Matrix2d & pileupCovariance,
const std::vector< double > & towerExpectations,
const std::vector< double > & towerVariances,
const std::vector< SignedKinematics > & toCorrect,
double constraintImportance = 1 )

Perform the pile-up fit.

Parameters
pileupSum2-vector sum over the pileup contributions
pileupCovarianceCovariance matrix of the pileup sum
towerExpectationsExpected values for each tower based on the pileup energy density
towerVariancesExpected values for the variances of each tower based on variances of background towers
toCorrectKinematics of the objects to correct
Therelative importance of the constraints.
Returns
The corrections to apply

Note that this version still returns the corrections that should be applied, not the corrected objects. The corrections are returned using the negatives of the value calculated by the fit - i.e. you should add the returned objects and let them act as negative energy objects (unless the sign of the correction was negative in which case it will act as a positive energy object). This variant provides individual estimates of the expected energy in each tower and the corresponding variances - these are for cases where the areas of each object are different.

Definition at line 222 of file PufitUtils.cxx.

229 {
230 // Keep track of the number of corrections to derive
231 std::size_t nCorr = toCorrect.size();
232 if (nCorr == 0)
233 // No need to do anything to derive corrections to nothing
234 return {};
235
236 // Get the cos-sin matrix
237 Eigen::Matrix<double, 2, Eigen::Dynamic> cosSin(2, nCorr);
238 Eigen::VectorXd vecTowerExpectations(nCorr);
239 Eigen::VectorXd vecTowerVariances(nCorr);
240 for (std::size_t ii = 0; ii < nCorr; ++ii)
241 {
242 cosSin(0, ii) = toCorrect.at(ii).cosPhi();
243 cosSin(1, ii) = toCorrect.at(ii).sinPhi();
244 vecTowerExpectations(ii) = towerExpectations.at(ii);
245 vecTowerVariances(ii) = towerVariances.at(ii);
246 }
247 Eigen::VectorXd corrections = pufit(
248 pileupSum,
249 pileupCovariance,
250 vecTowerExpectations,
251 vecTowerVariances,
252 cosSin,
253 constraintImportance);
254 // Now convert these into kinematics
255 std::vector<SignedKinematics> kinCorrections;
256 kinCorrections.reserve(nCorr);
257 for (std::size_t ii = 0; ii < nCorr; ++ii)
258 kinCorrections.push_back(SignedKinematics::fromEtEtaPhi(
259 -corrections.coeff(ii),
260 toCorrect.at(ii).eta(),
261 toCorrect.at(ii).phi()));
262 return kinCorrections;
263 }
static SignedKinematics fromEtEtaPhi(double et, double eta, double phi)
Factory function to construct from et, eta, phi (massless)

◆ pufit() [4/6]

Eigen::VectorXd HLT::MET::PufitUtils::pufit ( const Eigen::Vector2d & pileupSum,
const Eigen::Matrix2d & pileupCovariance,
double towerMean,
double towerVariance,
const Eigen::Matrix< double, 2, Eigen::Dynamic > & cosSin,
double constraintImportance = 1 )

Perform the pile-up fit.

Parameters
pileupSum2-vector sum over the pileup contributions
pileupCovarianceCovariance matrix of the pileup sum
towerMeanEstimate of the mean pileup tower transverse energy
towerVarianceEstimate of the variance of pileup tower transverse energies
cosSinThe cosines and sines of the objects to correct
Therelative importance of the constraints.
Returns
A vector containing the magnitudes of the corrections

This variant directly provides the cos-sin matrix for the correction, potentially allowing many cos/sin calls to be avoided.

Definition at line 204 of file PufitUtils.cxx.

211 {
212 std::size_t nCorr = cosSin.cols();
213 return pufit(
214 pileupSum,
215 pileupCovariance,
216 Eigen::VectorXd::Constant(nCorr, towerMean),
217 Eigen::VectorXd::Constant(nCorr, towerVariance),
218 cosSin,
219 constraintImportance);
220 }

◆ pufit() [5/6]

Eigen::VectorXd HLT::MET::PufitUtils::pufit ( const Eigen::Vector2d & pileupSum,
const Eigen::Matrix2d & pileupCovariance,
double towerMean,
double towerVariance,
const Eigen::VectorXd & correctionDirections,
double constraintImportance = 1 )

Perform the pile-up fit.

Parameters
pileupSum2-vector sum over the pileup contributions
pileupCovarianceCovariance matrix of the pileup sum
towerMeanEstimate of the mean pileup tower transverse energy
towerVarianceEstimate of the variance of pileup tower transverse energies
correctionDirectionsThe directions of the objects to correct
Therelative importance of the constraints.
Returns
A vector containing the magnitudes of the corrections

Definition at line 145 of file PufitUtils.cxx.

152 {
153 std::size_t nCorr = correctionDirections.size();
154 return pufit(
155 pileupSum,
156 pileupCovariance,
157 Eigen::VectorXd::Constant(nCorr, towerMean),
158 Eigen::VectorXd::Constant(nCorr, towerVariance),
159 correctionDirections,
160 constraintImportance);
161 }

◆ pufit() [6/6]

std::vector< SignedKinematics > HLT::MET::PufitUtils::pufit ( const Eigen::Vector2d & pileupSum,
const Eigen::Matrix2d & pileupCovariance,
double towerMean,
double towerVariance,
const std::vector< SignedKinematics > & toCorrect,
double constraintImportance = 1 )

Perform the pile-up fit.

Parameters
pileupSum2-vector sum over the pileup contributions
pileupCovarianceCovariance matrix of the pileup sum
towerMeanEstimate of the mean pileup tower transverse energy
towerVarianceEstimate of the variance of pileup tower transverse energies
toCorrectKinematics of the objects to correct
Therelative importance of the constraints.
Returns
The corrections to apply

Note that this version still returns the corrections that should be applied, not the corrected objects. The corrections are returned using the negatives of the value calculated by the fit - i.e. you should add the returned objects and let them act as negative energy objects (unless the sign of the correction was negative in which case it will act as a positive energy object).

Definition at line 264 of file PufitUtils.cxx.

271 {
272 // Keep track of the number of corrections to derive
273 std::size_t nCorr = toCorrect.size();
274 if (nCorr == 0)
275 // No need to do anything to derive corrections to nothing
276 return {};
277
278 return pufit(
279 pileupSum,
280 pileupCovariance,
281 std::vector<double>(nCorr, towerMean),
282 std::vector<double>(nCorr, towerVariance),
283 toCorrect,
284 constraintImportance);
285 }

◆ selectGrid() [1/2]

GridDisplacement HLT::MET::PufitUtils::selectGrid ( const PufitGridSet & grids)

Select the grid with the highest masked sumEt.

Definition at line 100 of file PufitUtils.cxx.

101 {
103 double maxSum = 0;
104 for (std::size_t d = 0; d < 4; ++d)
105 {
106 double sum = 0;
107 for (const PufitGrid::Tower &tower : grids[GridDisplacement(d)])
108 if (tower.masked())
109 sum += tower.sumEt();
110 if (sum > maxSum)
111 {
112 maximum = GridDisplacement(d);
113 maxSum = sum;
114 }
115 }
116 return maximum;
117 }
Describes a single element of the grid.
Definition PufitGrid.h:48
GridDisplacement
Enum to describe the positioning of the grid.
@ NoDisplacement
The grid is not shifted.

◆ selectGrid() [2/2]

template<typename Grid>
GridDisplacement HLT::MET::PufitUtils::selectGrid ( const PufitMultiGridSet< Grid > & grids,
std::size_t type )

Select the grid with the highest masked sumEt.

◆ trimmedMeanAndVariance() [1/3]

void HLT::MET::PufitUtils::trimmedMeanAndVariance ( const PufitGrid & grid,
double trimFraction,
double & mean,
double & variance )

Calculate the trimmed mean and variance for a vector of tower sumEts.

Parameters
gridThe grid to calculate
trimFractionThe fraction of towers to use in the calculation
[out]meanThe calculated mean
[out]varianceThe calculated variance

Definition at line 58 of file PufitUtils.cxx.

63 {
64 std::vector<double> sorted;
65 sorted.reserve(grid.nTowers());
66 for (const PufitGrid::Tower &tower : grid)
67 sorted.insert(
68 std::lower_bound(sorted.begin(), sorted.end(), tower.sumEt()),
69 tower.sumEt());
70 trimmedMeanAndVariance(sorted, trimFraction, mean, variance);
71 }
std::size_t nTowers() const
The number of bins.
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
void trimmedMeanAndVariance(const std::vector< double > &sorted, double trimFraction, double &mean, double &variance)
Calculate the trimmed mean and variance for a vector of tower sumEts.

◆ trimmedMeanAndVariance() [2/3]

template<std::size_t N>
void HLT::MET::PufitUtils::trimmedMeanAndVariance ( const PufitMultiGrid< N > & grid,
std::size_t type,
double trimFraction,
double & mean,
double & variance )

Calculate the trimmed mean and variance for a vector of tower sumEts.

Parameters
gridThe grid to calculate
typeThe int mask of subgrid types
trimFractionThe fraction of twoers to use in the calculation
[out]meanThe calculated mean
[out]varianceThe calculated variance

◆ trimmedMeanAndVariance() [3/3]

void HLT::MET::PufitUtils::trimmedMeanAndVariance ( const std::vector< double > & sorted,
double trimFraction,
double & mean,
double & variance )

Calculate the trimmed mean and variance for a vector of tower sumEts.

Parameters
sortedThe sumEts for all towers, sorted in ascending order
trimFractionThe fraction of towers to use in the calculation
[out]meanThe calculated mean
[out]varianceThe calculated variance

Definition at line 36 of file PufitUtils.cxx.

41 {
42 // The number to trim from each side
43 std::size_t nTrim = sorted.size() * (1 - trimFraction) / 2;
44
45 mean = std::accumulate(sorted.begin() + nTrim, sorted.end() - nTrim, 0.);
46 mean /= (sorted.size() - 2 * nTrim);
47
48 // Now get the variance
49 variance = 0;
50 for (std::size_t ii = 0; ii < sorted.size() - nTrim; ++ii)
51 {
52 double val = sorted.at(ii) - mean;
53 variance += val * val * (ii < nTrim ? 2 : 1);
54 }
55 variance /= sorted.size();
56 }

◆ unmaskedMeanAndVariance() [1/2]

void HLT::MET::PufitUtils::unmaskedMeanAndVariance ( const PufitGrid & grid,
double & mean,
double & variance )

Calculate the mean and variance of unmasked towers.

Parameters
gridThe grid to calculate
[out]meanThe calculated mean
[out]varianceThe calculated variance

Definition at line 73 of file PufitUtils.cxx.

77 {
78 double sum = 0;
79 double squaredSum = 0;
80 std::size_t n = 0;
81 for (const PufitGrid::Tower &tower : grid)
82 {
83 if (tower.masked())
84 continue;
85 ++n;
86 sum += tower.sumEt();
87 squaredSum += tower.sumEt() * tower.sumEt();
88 }
89 //if n is still zero here, something went very wrong
90 if (n == 0)[[unlikely]]{
91 throw std::runtime_error("unmaskedMeanAndVariance: n is zero in denominator.");
92 }
93 mean = sum / n;
94 // Note that this could result in catastrophic cancellation in some cases.
95 // However the amount of precision present in a double is around 10^15 so we
96 // can afford to lose some significant figures
97 variance = squaredSum / n - mean * mean;
98 }
#define unlikely(x)

◆ unmaskedMeanAndVariance() [2/2]

template<std::size_t N>
void HLT::MET::PufitUtils::unmaskedMeanAndVariance ( const PufitMultiGrid< N > & grid,
int type,
double & mean,
double & variance )

Calculate the mean and variance of unmasked towers.

Parameters
gridThe grid to calculate
typeThe int mask of subgrid types
[out]meanThe calculated mean
[out]varianceThe calculated variance