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 158 of file PufitUtils.cxx.

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

◆ 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 114 of file PufitUtils.cxx.

121 {
122 // Keep track of the number of corrections to derive
123 std::size_t nCorr = correctionDirections.size();
124 if (nCorr == 0)
125 // No need to do anything to derive corrections to nothing
126 return {};
127 // Turn the directions into cos/sin pairs
128 Eigen::Matrix<double, 2, Eigen::Dynamic> cosSin(2, nCorr);
129 cosSin.row(0) = correctionDirections.array().cos();
130 cosSin.row(1) = correctionDirections.array().sin();
131 return pufit(
132 pileupSum,
133 pileupCovariance,
134 towerExpectations,
135 towerVariances,
136 cosSin,
137 constraintImportance);
138 }
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 217 of file PufitUtils.cxx.

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

206 {
207 std::size_t nCorr = cosSin.cols();
208 return pufit(
209 pileupSum,
210 pileupCovariance,
211 Eigen::VectorXd::Constant(nCorr, towerMean),
212 Eigen::VectorXd::Constant(nCorr, towerVariance),
213 cosSin,
214 constraintImportance);
215 }

◆ 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 140 of file PufitUtils.cxx.

147 {
148 std::size_t nCorr = correctionDirections.size();
149 return pufit(
150 pileupSum,
151 pileupCovariance,
152 Eigen::VectorXd::Constant(nCorr, towerMean),
153 Eigen::VectorXd::Constant(nCorr, towerVariance),
154 correctionDirections,
155 constraintImportance);
156 }

◆ 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 259 of file PufitUtils.cxx.

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

◆ selectGrid() [1/2]

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

Select the grid with the highest masked sumEt.

Definition at line 95 of file PufitUtils.cxx.

96 {
98 double maxSum = 0;
99 for (std::size_t d = 0; d < 4; ++d)
100 {
101 double sum = 0;
102 for (const PufitGrid::Tower &tower : grids[GridDisplacement(d)])
103 if (tower.masked())
104 sum += tower.sumEt();
105 if (sum > maxSum)
106 {
107 maximum = GridDisplacement(d);
108 maxSum = sum;
109 }
110 }
111 return maximum;
112 }
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 57 of file PufitUtils.cxx.

62 {
63 std::vector<double> sorted;
64 sorted.reserve(grid.nTowers());
65 for (const PufitGrid::Tower &tower : grid)
66 sorted.insert(
67 std::lower_bound(sorted.begin(), sorted.end(), tower.sumEt()),
68 tower.sumEt());
69 trimmedMeanAndVariance(sorted, trimFraction, mean, variance);
70 }
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 35 of file PufitUtils.cxx.

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

◆ 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 72 of file PufitUtils.cxx.

76 {
77 double sum = 0;
78 double squaredSum = 0;
79 std::size_t n = 0;
80 for (const PufitGrid::Tower &tower : grid)
81 {
82 if (tower.masked())
83 continue;
84 ++n;
85 sum += tower.sumEt();
86 squaredSum += tower.sumEt() * tower.sumEt();
87 }
88 mean = sum / n;
89 // Note that this could result in catastrophic cancellation in some cases.
90 // However the amount of precision present in a double is around 10^15 so we
91 // can afford to lose some significant figures
92 variance = squaredSum / n - mean * mean;
93 }

◆ 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