ATLAS Offline Software
PufitUtils.cxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3  */
4 
6 #include <algorithm>
7 #include <numeric>
8 
9 namespace HLT
10 {
11  namespace MET
12  {
13  namespace PufitUtils
14  {
15 
16  CovarianceSum::CovarianceSum() : sum(Eigen::Vector2d::Zero()),
17  covariance(Eigen::Matrix2d::Zero())
18  {
19  }
20 
22  const Eigen::Vector2d &sum,
23  const Eigen::Matrix2d &covariance) : sum(sum), covariance(covariance) {}
24 
26  const SignedKinematics &kin,
27  double sigma)
28  {
29  Eigen::Vector2d cosSin(kin.cosPhi(), kin.sinPhi());
30  sum += Eigen::Vector2d(kin.px(), kin.py());
31  covariance += sigma * cosSin * cosSin.transpose();
32  return *this;
33  }
34 
36  const std::vector<double> &sorted,
37  double trimFraction,
38  double &mean,
39  double &variance)
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  }
56 
58  const PufitGrid &grid,
59  double trimFraction,
60  double &mean,
61  double &variance)
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  }
71 
73  const PufitGrid &grid,
74  double &mean,
75  double &variance)
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  }
94 
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  }
113 
114  Eigen::VectorXd pufit(
115  const Eigen::Vector2d &pileupSum,
116  const Eigen::Matrix2d &pileupCovariance,
117  const Eigen::VectorXd &towerExpectations,
118  const Eigen::VectorXd &towerVariances,
119  const Eigen::VectorXd &correctionDirections,
120  double constraintImportance)
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  }
139 
140  Eigen::VectorXd pufit(
141  const Eigen::Vector2d &pileupSum,
142  const Eigen::Matrix2d &pileupCovariance,
143  double towerMean,
144  double towerVariance,
145  const Eigen::VectorXd &correctionDirections,
146  double constraintImportance)
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  }
157 
158  Eigen::VectorXd pufit(
159  const Eigen::Vector2d &pileupSum,
160  const Eigen::Matrix2d &pileupCovariance,
161  const Eigen::VectorXd &towerExpectations,
162  const Eigen::VectorXd &towerVariances,
163  const Eigen::Matrix<double, 2, Eigen::Dynamic> &cosSin,
164  double constraintImportance)
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  }
198 
199  Eigen::VectorXd pufit(
200  const Eigen::Vector2d &pileupSum,
201  const Eigen::Matrix2d &pileupCovariance,
202  double towerMean,
203  double towerVariance,
204  const Eigen::Matrix<double, 2, Eigen::Dynamic> &cosSin,
205  double constraintImportance)
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  }
216 
217  std::vector<SignedKinematics> pufit(
218  const Eigen::Vector2d &pileupSum,
219  const Eigen::Matrix2d &pileupCovariance,
220  const std::vector<double> &towerExpectations,
221  const std::vector<double> &towerVariances,
222  const std::vector<SignedKinematics> &toCorrect,
223  double constraintImportance)
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  }
259  std::vector<SignedKinematics> pufit(
260  const Eigen::Vector2d &pileupSum,
261  const Eigen::Matrix2d &pileupCovariance,
262  double towerMean,
263  double towerVariance,
264  const std::vector<SignedKinematics> &toCorrect,
265  double constraintImportance)
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  }
281  } // namespace PufitUtils
282  } // namespace MET
283 } // namespace HLT
HLT::MET::SignedKinematics::cosPhi
double cosPhi() const
Definition: SignedKinematics.cxx:79
HLT::MET::SignedKinematics::px
double px() const
Definition: SignedKinematics.cxx:121
HLT::MET::PufitGrid
Bins energy deposits into a grid.
Definition: PufitGrid.h:50
HLT::MET::PufitUtils::pufit
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.
Definition: PufitUtils.cxx:114
HLT::MET::PufitUtils::CovarianceSum::sum
Eigen::Vector2d sum
The sum.
Definition: PufitUtils.h:61
pdg_comparison.sigma
sigma
Definition: pdg_comparison.py:324
mean
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="")
Definition: dependence.cxx:254
hist_file_dump.d
d
Definition: hist_file_dump.py:137
accumulate
bool accumulate(AccumulateMap &map, std::vector< module_t > const &modules, FPGATrackSimMatrixAccumulator const &acc)
Accumulates an accumulator (e.g.
Definition: FPGATrackSimMatrixAccumulator.cxx:22
HLT::MET::PufitGridSet
Helper struct to contain a full set of grids.
Definition: PufitGrid.h:252
HLT::MET::PufitUtils::CovarianceSum::covariance
Eigen::Matrix2d covariance
The covariance matrix.
Definition: PufitUtils.h:63
HLT::MET::PufitUtils::CovarianceSum
Helper struct to hold the sum over pileup objects and its covariance.
Definition: PufitUtils.h:44
HLT::MET::SignedKinematics::sinPhi
double sinPhi() const
Provide accessors for sin and cos phi.
Definition: SignedKinematics.cxx:73
HLT::MET::NoDisplacement
@ NoDisplacement
The grid is not shifted.
Definition: PeriodicGridBase.h:25
HLT::MET::GridDisplacement
GridDisplacement
Enum to describe the positioning of the grid.
Definition: PeriodicGridBase.h:23
convertTimingResiduals.sum
sum
Definition: convertTimingResiduals.py:55
HLT::MET::SignedKinematics
Class to describe the kinematics of an object that can have negative energies.
Definition: SignedKinematics.h:42
HLT
It used to be useful piece of code for replacing actual SG with other store of similar functionality ...
Definition: HLTResultReader.h:26
beamspotman.n
n
Definition: beamspotman.py:731
PufitUtils.h
HLT::MET::PufitGrid::Tower
Describes a single element of the grid.
Definition: PufitGrid.h:60
TCS::MET
@ MET
Definition: Trigger/TrigT1/L1Topo/L1TopoCommon/L1TopoCommon/Types.h:16
HLT::MET::PufitUtils::selectGrid
GridDisplacement selectGrid(const PufitGridSet &grids)
Select the grid with the highest masked sumEt.
Definition: PufitUtils.cxx:95
DerivationFramework::TriggerMatchingUtils::sorted
std::vector< typename T::value_type > sorted(T begin, T end)
Helper function to create a sorted vector from an unsorted one.
HLT::MET::PufitUtils::trimmedMeanAndVariance
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.
Definition: PufitUtils.cxx:35
HLT::MET::PufitUtils::unmaskedMeanAndVariance
void unmaskedMeanAndVariance(const PufitGrid &grid, double &mean, double &variance)
Calculate the mean and variance of unmasked towers.
Definition: PufitUtils.cxx:72
HLT::MET::PufitUtils::CovarianceSum::add
CovarianceSum & add(const SignedKinematics &kin, double sigma)
Add a new contribution to the sum.
Definition: PufitUtils.cxx:25
constants
Definition: Calorimeter/CaloClusterCorrection/python/constants.py:1
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
HLT::MET::PufitUtils::CovarianceSum::CovarianceSum
CovarianceSum()
Default constructor - zero initialize everything.
Definition: PufitUtils.cxx:16
HLT::MET::PeriodicGridBase::nTowers
std::size_t nTowers() const
The number of bins.
Definition: PeriodicGridBase.cxx:142
HLT::MET::SignedKinematics::fromEtEtaPhi
static SignedKinematics fromEtEtaPhi(double et, double eta, double phi)
Factory function to construct from et, eta, phi (massless)
Definition: SignedKinematics.cxx:46
HLT::MET::SignedKinematics::py
double py() const
Definition: SignedKinematics.cxx:124