ATLAS Offline Software
Loading...
Searching...
No Matches
LArCalorimeter/LArSamplesMon/src/Residual.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include "TMath.h"
8#include <algorithm>
12
13#include "TRandom.h"
14#include "TH1D.h"
15
16#include <iostream>
17using std::cout;
18using std::endl;
19
20using namespace LArSamples;
21
22
23Residual::Residual(const TVectorD& deltas, int run, int event, double adcMax, double time)
26{
28}
29
30
35
36
37TVectorD Residual::scaledDeltas() const
38{
39 TVectorD vect = deltas();
40 vect *= 1/adcMax();
41 return vect;
42}
43
44
45double Residual::scaledDelta(short i) const
46{
47 if (!isInRange(i) || adcMax() <= 0) return 0;
48 return deltas()(i)/adcMax();
49}
50
51
53{
54 TVectorD vect(lwb(), upb() + 1);
55 vect.SetSub(lwb(), scaledDeltas());
56 vect(upb() + 1) = time();
57 return vect;
58}
59
60
61bool ResidualCompare::operator()(const Residual& r1, const Residual& r2) const
62{
63 if (!r1.hasSameRange(r2)) {
64 cout << "WARNING : residuals have different ranges, this may crash the sorting..." << endl;
65 return false;
66 }
67
68 // Very, very important! Since samples are integers, it *will* happen that some residuals will have *identical* deltas.
69 // It turns out that (at least in the way used here), operator<(double, double) is not a strict weak ordering: we will
70 // have cases where residual1[4] == residual2[4], and then we will have residual1[4] < residual2[4] and residual2[4] < residual1[4]
71 // both true, which will confuse std::sort and cause it to crash... so add this safety (lesson of 2 days of debugging)
72 // Note: if we write a<<b iff a < b - epsilon, then << is indeed a strict weak ordering so all is fine. It does mean there
73 // will be "sorting errors", i.e. residuals that differ by less than epsilon may be sorted the wrong way. But that's fine for our usage.
74 const double epsilon = 1E-5;
75// cout << "Comparing sample " << m_sampling << ", = "
76// << (m_sampling == r1.upb() + 1 ? r1.time() : r1.scaledDelta(m_sampling)) << " >/< "
77// << (m_sampling == r2.upb() + 1 ? r2.time() : r2.scaledDelta(m_sampling)) << " "
78// << ((m_sampling == r1.upb() + 1 ? r1.time() : r1.scaledDelta(m_sampling)) < (m_sampling == r2.upb() + 1 ? r2.time() : r2.scaledDelta(m_sampling))) << endl;
79
80 if (r1.isInRange(m_sampling)) return r1.scaledDelta(m_sampling) < r2.scaledDelta(m_sampling) - epsilon;
81 if (m_sampling == r1.upb() + 1) return r1.time() < r2.time() - epsilon;
82 cout << "WARNING : comparison on an unknown sample, this may crash the sorting..." << endl;
83 return false;
84}
85
86
87bool Residuals::medianVars(TVectorD& medians, TVectorD& widths) const
88{
89 if (size() == 0) return false;
90 medians.ResizeTo(lwb(), upb() + 1);
91 widths.ResizeTo (lwb(), upb() + 1);
92
93 double halfSigmaQuantile = TMath::Prob(1,1)/2;
94 double medianQuantile = 0.5;
95
96 std::vector<Residual> sortedResiduals = m_residuals;
97
98 for (short i = lwb(); i <= upb() + 1; i++) {
99// cout << "test-sorting sample " << i << " (n = " << m_residuals.size() << ")" << endl;
100// std::vector<const Residual*> testSort(m_residuals.size());
101// ResidualCompare comparer(i);
102// for (unsigned int i1 = 0; i1 < m_residuals.size(); i1++) {
103// const Residual& res1 = m_residuals[i1];
104// unsigned int nBefore = 0;
105// for (unsigned int i2 = 0; i2 < m_residuals.size(); i2++) {
106// if (i1 == i2) continue;
107// const Residual& res2 = m_residuals[i2];
108// //cout << "comparing " << i1 << " " << i2 << endl;
109// if (comparer(res2, res1)) nBefore++;
110// if ((i1 == 3617 || i1 == 3886) && comparer(res2, res1))
111// cout << i2 << " is bfore " << i1 << " " << res2.scaledDelta(i) << " < " << res1.scaledDelta(i)
112// << (res2.scaledDelta(i) < res1.scaledDelta(i)) << " " << (res2.scaledDelta(i) < res1.scaledDelta(i) - 1E-3) << endl;
113// }
114// cout << "index " << i1 << " : nBefore = " << nBefore << " " << res1.scaledDelta(i) << endl;
115// while (testSort[nBefore]) { nBefore++; if (nBefore == testSort.size()) cout << "ERROR!!! " << i1 << endl; }
116// testSort[nBefore] = &res1;
117// }
118// cout << "done sorting" << endl;
119// for (unsigned int i1 = 0; i1 < testSort.size(); i1++) {
120// const Residual* res1 = testSort[i1];
121// if (!res1) cout << "ERROR : null at index " << i1 << endl;
122// for (unsigned int i2 = 0; i2 < i1; i2++) {
123// const Residual* res2 = testSort[i2];
124// if (comparer(*res1, *res2)) cout << "Wrong order " << i1 << " " << i2 << endl;
125// }
126// }
127// cout << "done checking" << endl;
128 std::sort(sortedResiduals.begin(), sortedResiduals.end(), ResidualCompare(i));
129 Residual& medianRes = sortedResiduals[(unsigned int)(sortedResiduals.size()*medianQuantile)];
130 Residual& loHalfSigmaRes = sortedResiduals[(unsigned int)(sortedResiduals.size()*halfSigmaQuantile)];
131 Residual& hiHalfSigmaRes = sortedResiduals[(unsigned int)(sortedResiduals.size()*(1 - halfSigmaQuantile))];
132 double median = (i < upb() + 1 ? medianRes.scaledDelta(i) : medianRes.time()) ;
133 double loHalfSigma = (i < upb() + 1 ? loHalfSigmaRes.scaledDelta(i) : loHalfSigmaRes.time());
134 double hiHalfSigma = (i < upb() + 1 ? hiHalfSigmaRes.scaledDelta(i) : hiHalfSigmaRes.time());
135 //cout << loHalfSigma << " " << median << " " << hiHalfSigma << endl;
136 medians[i] = median;
137 widths[i] = (hiHalfSigma - median > median - loHalfSigma ? median - loHalfSigma : hiHalfSigma - median);
138 }
139 return true;
140}
141
142
143std::unique_ptr<Residuals> Residuals::truncate(double nWidthsRes, double nWidthsTime, unsigned int nMax) const
144{
145 if (size() == 0) return std::make_unique<Residuals>();
146 TVectorD medians, widths;
147
148 if (nMax > 0) {
149 auto original = std::make_unique<Residuals>();
150 for (unsigned int i = 0; i < nMax; i++) original->add(*residual(i));
151 if (!original->medianVars(medians, widths)) return nullptr;
152 }
153 else {
154 if (!medianVars(medians, widths)) return nullptr;
155 }
156 auto truncated = std::make_unique<Residuals>();
157
158 for (const Residual& residual : m_residuals) {
159 bool pass = true;
160 if (nMax > 0 && truncated->size() == nMax) break;
161 for (short i = lwb(); i <= upb(); i++) {
162 if (nWidthsRes > 0 && TMath::Abs(residual.scaledDelta(i) - medians[i]) > nWidthsRes*widths[i]) {
163 pass = false;
164 break;
165 }
166 }
167 if (!pass) continue;
168 if (nWidthsTime > 0 && TMath::Abs(residual.time() - medians[upb() + 1]) > nWidthsTime*widths[upb() + 1]) continue;
169 truncated->add(residual);
170 }
171 return truncated;
172}
173
174
175TH1D* Residuals::histogram(short sample, const TString& name, int nBins, double xMin, double xMax) const
176{
177 TH1D* h = new TH1D(name, "", nBins, xMin, xMax);
178 for (const Residual& residual : m_residuals)
179 h->Fill(sample <= upb() ? residual.scaledDelta(sample) : residual.time());
180 return h;
181}
182
183
185{
186 if (size() == 0) return nullptr;
187
188 ResidualCalculator* calc = new ResidualCalculator(lwb(), upb(), weigh);
189 for (const Residual& residual : m_residuals)
190 calc->fill(residual);
191
192 return calc;
193}
194
195
196int ResidualCalculator::find(int r, int e) const
197{
198 for (unsigned int i = 0; i < size(); i++)
199 if (event(i) == e && run(i) == r) return i;
200 return -1;
201}
202
203
205{
206 m_runs.push_back(r);
207 m_events.push_back(e);
208 return true;
209}
210
211
213{
214 if (find(residual.run(), residual.event()) >= 0) return true;
215 add(residual.run(), residual.event());
216 return fill_with_weight(residual, +1);
217}
218
219
221{
222 //std::pair<int, int> ev(residual.run(), residual.event());
223 //std::map< std::pair<int, int>, bool>::iterator findEv = m_events.find(ev);
224 //if (findEv == m_events.end()) return true;
225 //cout << "Will remove run = " << ev.first << ", event = " << ev.second << endl;
226 //m_events.erase(findEv);
227 int index = find(residual.run(), residual.event());
228 if (index == -1) return true;
229 m_events[index] = -1;
230 m_runs[index] = -1;
231 return fill_with_weight(residual, -1);
232}
233
234
235std::unique_ptr<ShapeErrorData> ResidualCalculator::shapeErrorData() const
236{
237 TVectorD xi = regresser()->means().GetSub(lwb(), upb(), "I");
238 CovMatrix xiErr = regresser()->meanErrorMatrix().GetSub(lwb(), upb(), lwb(), upb(), "I");
239 double tbar = regresser()->mean(upb() + 1);
240
241 double denom = regresser()->covarianceMatrix()(upb() + 1, upb() + 1);
242 if (denom < 1E-6) {
243 TVectorD xip(lwb(), upb());
244 CovMatrix xipErr(lwb(), upb());
245 // happens for size==2 if we are removing one of the residuals.
246 if (size() > 2) cout << "WARNING: variance of t < 1E-6, returning correction without derivative term. (V = " << denom << ", N = " << size() << ")" << endl;
247 return std::make_unique<ShapeErrorData>(xi, xip, xiErr, xipErr, tbar, regresser()->nEntries());
248 }
249
250 TVectorD xip = TVectorD(regresser()->covarianceMatrix()[upb() + 1]).GetSub(lwb(), upb(), "I");
251 xip *= -1/denom;
252
253 TVectorD xipErrVect = TVectorD(regresser()->covarianceMatrixErrors()[upb() + 1]).GetSub(lwb(), upb(), "I");
254 xipErrVect *= 1/denom;
255 CovMatrix xipErr(lwb(), upb());
256 for (int k1 = lwb(); k1 <= upb(); k1++) xipErr(k1, k1) = TMath::Power(xipErrVect(k1), 2);
257
258 return std::make_unique<ShapeErrorData>(xi, xip, xiErr, xipErr, tbar, regresser()->nEntries());
259}
260
261
262double ResidualCalculator::weight(const Residual& residual) const
263{
264 if (!weigh()) return 1;
265 return TMath::Power(residual.adcMax(), 2);
266}
267
268
269bool ResidualCalculator::fill_with_weight(const Residual& residual, double w)
270{
271 w *= weight(residual);
272 return m_regresser.fill(residual.scaledDeltasAndTime(), w);
273}
274
275
277{
278 if (!m_regresser.append(*other.regresser())) return false;
279 //for (std::map< std::pair<int, int >, bool >::const_iterator event = other.events().begin();
280 // event != other.events().end(); event++)
281 // m_events[event->first] = event->second;
282 for (unsigned int i = 0; i < other.size(); i++) add(other.run(i), other.event(i));
283 return true;
284}
285
286
288{
289 TString str = "\n";
290 /*for (std::map< std::pair<int, int>, bool>::const_iterator ev = m_events.begin();
291 ev != m_events.end(); ev++)
292 str += Form("run %6d event %10d\n", ev->first.first, ev->first.second);
293 */
294 for (unsigned int i = 0; i < size(); i++)
295 str += Form("run %6d event %10d\n", run(i), event(i));
296
297 return str;
298}
299
300
302{
303 Residuals testVect;
304 TRandom r;
305 for (unsigned int i = 0; i < 100000; i++) {
306 TVectorD d(5);
307 //for (unsigned short j = 0; j < 5; j++) d(j) = 0.01*i;
308 for (unsigned short j = 0; j < 5; j++) d(j) = r.Gaus(0, 1);
309 testVect.add(Residual(d, i, i, i, i));
310 }
311
312// TVectorD m,w;
313// if (!getMedianVars(testVect, m, w)) return false;
314// m.Print();
315// w.Print();
316
317 std::unique_ptr<Residuals> truncated = testVect.truncate(2);
318 if (!truncated) return false;
319
320 TH1D* h = truncated->histogram(0, "h", 100, -5, 5);
321 h->Draw();
322 return true;
323}
324
double mean(unsigned int i) const
Definition Averager.cxx:191
CovMatrix meanErrorMatrix() const
Definition Averager.cxx:138
TVectorD means() const
Definition Averager.cxx:112
CovMatrix covarianceMatrix() const
Definition Averager.cxx:153
bool isInRange(int i) const
Definition IndexRange.h:27
std::unique_ptr< ShapeErrorData > shapeErrorData() const
bool fill_with_weight(const Residual &residual, double w)
bool operator()(const Residual &r1, const Residual &r2) const
Residual(const TVectorD &deltas=TVectorD(), int run=0, int event=0, double adcMax=-1, double time=0)
bool medianVars(TVectorD &medians, TVectorD &widths) const
std::unique_ptr< Residuals > truncate(double nWidthsRes, double nWidthsTime=-1, unsigned int nMax=0) const
const Residual * residual(unsigned int i) const
TH1D * histogram(short sample, const TString &name, int nBins, double xMin, double xMax) const
ResidualCalculator * calculator(bool weigh=false) const
int r
Definition globals.cxx:22
Definition index.py:1
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.