ATLAS Offline Software
Loading...
Searching...
No Matches
DigitMonitor.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
10#include "LArSamplesMon/OFC.h"
17#include "LArCafJobs/Geometry.h"
18
19#include "TFile.h"
20#include "TCanvas.h"
21#include "TF1.h"
22#include "TH1I.h"
23#include "TH2D.h"
24#include "TLatex.h"
25#include "TMath.h"
26#include "TPaveText.h"
27#include "TROOT.h"
28#include "TVectorD.h"
29#include <map>
30#include <utility>
31
32#include <vector>
33
34#include <iostream>
35using std::cout;
36using std::endl;
37
38using namespace LArSamples;
39
40
41TH1D* DigitMonitor::energyDist(const TString& name, int nBins, double eMax) const
42{return dist(&Data::_energy, DataFuncArgs(), name, nBins, 0, eMax, "Energy distribution", "Energy (MeV)"); }
43
44TH1D* DigitMonitor::timeDist(const TString& name, int nBins, double tMin, double tMax) const
45{ return dist(&Data::_ofcTime, DataFuncArgs(), name, nBins, tMin, tMax, "Time distribution", "Time (ns)"); }
46
47TH1D* DigitMonitor::adcMaxDist(const TString& name, int nBins, double aMin, double aMax) const
48{ return dist(&Data::_adcMax, DataFuncArgs(), name, nBins, aMin, aMax, "ADCMax distribution", "ADCMax (counts)"); }
49
50TH1D* DigitMonitor::noiseDist(const TString& name, int nBins, double max) const
51{ return dist(&Data::_noise, DataFuncArgs(), name, nBins, 0, max, "Noise distribution", "Noise (ADC counts)"); }
52
53TH2D* DigitMonitor::maxValueMap(TString name, PartitionId partition) const
54{ return partitionMap(&Data::_maxValue, DataFuncArgs(), std::move(name), partition, "Max sample value (ADC counts)", MaxValue); }
55
56TH2D* DigitMonitor::minValueMap(TString name, PartitionId partition) const
57{ return partitionMap(&Data::_minValue, DataFuncArgs(), std::move(name), partition, "Min sample value (ADC counts)", MinValue); }
58
59
60TH1D* DigitMonitor::chi2Dist(const TString& name, int nBins, double max, double kFactor,
61 double fitMax, int lwb, int upb, unsigned int chi2Pars,
62 ShapeErrorType shapeErrorType, unsigned int nDof) const
63{
64 TH1D* h = new TH1D(name, "#chi^{2} distribution", nBins, 0, max);
65 unsigned int nDofEff = 0, nDofLast = 0;
66 std::unique_ptr<UniformShapeErrorGetter> kFactorGetter;
67 std::unique_ptr<CombinedShapeErrorGetter> shapeErrorGetter;
68 if (kFactor > 0) {
69 kFactorGetter = std::make_unique<UniformShapeErrorGetter>(kFactor);
70 shapeErrorGetter = std::make_unique<CombinedShapeErrorGetter>();
71 if (interface().shapeErrorGetter()) shapeErrorGetter->add(*interface().shapeErrorGetter());
72 shapeErrorGetter->add(*kFactorGetter);
73 }
74 for (unsigned int i = 0; i < nChannels(); i++) {
75 if ((i+1) % 10000 == 0) cout << "Processing entry # " << i+1 << endl;
76 const History* history = cellHistory(i);
77 if (!history) continue;
78 if (shapeErrorGetter) history->setShapeErrorGetter(shapeErrorGetter.get()); // Do it at history level so it doesn't "stick" afterwards,,,
79 for (unsigned int j = 0; j < history->nData(); j++) {
80 h->Fill(history->chi2(j, lwb, upb, chi2Pars, shapeErrorType, &nDofEff));
81 if (nDofLast > 0 && nDofLast != nDofEff) cout << "WARNING @ hash = " << i << ", index = " << j << " : nDof varied from " << nDofLast << " to " << nDofEff << endl;
82 nDofLast = nDofEff;
83 }
84 }
85
86 if (nDof > 0 && nDof != nDofEff) cout << "WARNING : fixing nDof = " << nDof << ", but computer nDofEff = " << nDofEff << endl;
87 h->Draw();
88 if (fitMax < 0) return h;
89 TF1* fChi2 = fitChi2(*h, "chi2", 0, max, nDof, 0, fitMax);
90 TPaveText* p = new TPaveText(0.6, 0.85, 0.85, 0.65, "NDC");
91 p->AddText(Form("k = %.1f%%", kFactor*100));
92 p->AddText(Form("Fit #chi^{2} = %.1f", fChi2->GetChisquare()));
93 p->AddText(Form("Fit nDOF = %.1f", fChi2->GetParameter(1)));
94 p->Draw();
95 return h;
96}
97
98
99TF1* DigitMonitor::chi2Func(const char* name, double xMin, double xMax)
100{
101 TF1* fChi2 = new TF1(name, "[0]*TMath::Power(x/2, [1]/2)/x*TMath::Exp(-x/2)/TMath::Gamma([1]/2)", xMin, xMax);
102 fChi2->SetParameter(0, 100);
103 fChi2->SetParameter(1, 3);
104 fChi2->SetParLimits(1, 1, 7);
105 return fChi2;
106}
107
108
109TF1* DigitMonitor::fitChi2(TH1D& h, const char* name, double xMin, double xMax,
110 double nDof, double fitMin, double fitMax)
111{
112 if (fitMin == Definitions::none) fitMin = xMin;
113 if (fitMax == Definitions::none) fitMin = xMax;
114
115 TF1* fChi2 = chi2Func(name, xMin, xMax);
116 fChi2->SetParameter(0, 4*h.GetMaximum());
117 if (nDof > 0) fChi2->FixParameter(1, nDof);
118
119 h.Fit(fChi2, "", "", fitMin, fitMax);
120 fChi2->SetRange(xMin, xMax);
121 return fChi2;
122}
123
124
125TH1D* DigitMonitor::bestChi2Dist(const TString& name, int nBins, double max,
126 double refErrMin, double refErrMax, unsigned int refErrNBins,
127 int lwb, int upb, unsigned int chi2Pars, unsigned int nDof) const
128{
129 TF1* fChi2 = chi2Func("fChi2", 0, max);
130 if (nDof > 0) fChi2->FixParameter(1, nDof);
131
132 double delta = (refErrMax - refErrMin)/refErrNBins;
133 TH1D* fitResults = new TH1D("fitResults", "", refErrNBins, refErrMin - delta/2, refErrMax - delta/2);
134
135 for (unsigned int k = 0; k < refErrNBins; k++) {
136 double refErr = refErrMin + delta*k;
137 TH1D* h = new TH1D(name, "#chi^{2} distribution", nBins, 0, max);
138 UniformShapeErrorGetter kFactorGetter(refErr);
139 CombinedShapeErrorGetter shapeErrorGetter;
140 if (interface().shapeErrorGetter()) shapeErrorGetter.add(*interface().shapeErrorGetter());
141 shapeErrorGetter.add(kFactorGetter);
142 for (unsigned int i = 0; i < nChannels(); i++) {
143 const History* history = cellHistory(i);
144 if (!history) continue;
145 history->setShapeErrorGetter(&shapeErrorGetter);
146 for (unsigned int j = 0; j < history->nData(); j++) {
147 h->Fill(history->chi2(j, lwb, upb, chi2Pars));
148 }
149 }
150 h->Fit(fChi2);
151 fitResults->SetBinContent(k+1, fChi2->GetChisquare());
152 delete h;
153 }
154
155 TF1* fPol2 = new TF1("fPol2", "[0]*(x - [1])*(x - [1]) + [2]", refErrMin, refErrMax);
156
157 double center = fitResults->GetBinCenter(fitResults->GetMinimumBin());
158 fPol2->SetParameter(0, (fitResults->GetMaximum() - fitResults->GetMinimum())/center/center);
159 fPol2->SetParameter(1, center);
160 fPol2->SetParameter(2, fitResults->GetMinimum());
161
162 fitResults->Draw();
163 fitResults->Fit(fPol2);
164
165 TPaveText* p = new TPaveText(0.6, 0.85, 0.85, 0.65, "NDC");
166 p->AddText(Form("Fit k_{best} = %.1f%%", fPol2->GetParameter(1)*100));
167 p->AddText(Form("Fit #delta k_{best} = %.2f%%", 1/TMath::Sqrt(fPol2->GetParameter(0))*100));
168 p->AddText(Form("Fit #chi^{2}_{min} = %.2f", fPol2->GetParameter(2)));
169 p->Draw();
170
171 return fitResults;
172}
173
174
175
176TH1D* DigitMonitor::gainDist(const TString& name) const
177{
178 TH1D* h = new TH1D(name, "gain distribution", 3, -0.5, 2.5);
179 h->GetXaxis()->SetBinLabel(1, "Low");
180 h->GetXaxis()->SetBinLabel(2, "Medium");
181 h->GetXaxis()->SetBinLabel(3, "High");
182 for (unsigned int i = 0; i < nChannels(); i++) {
183 const History* history = cellHistory(i);
184 if (!history) continue;
185 for (unsigned int j = 0; j < history->nData(); j++) h->Fill(history->data(j)->gain());
186 }
187 return h;
188}
189
190
191TH1D* DigitMonitor::layerDist(const TString& name) const
192{
193 TH1D* h = new TH1D(name, "layer distribution", 4, -0.5, 3.5);
194 h->GetXaxis()->SetBinLabel(1, "PS");
195 h->GetXaxis()->SetBinLabel(2, "Strips");
196 h->GetXaxis()->SetBinLabel(3, "Middle");
197 h->GetXaxis()->SetBinLabel(4, "Back");
198 for (unsigned int i = 0; i < nChannels(); i++) {
199 const History* history = cellHistory(i);
200 if (!history) continue;
201 for (unsigned int j = 0; j < history->nData(); j++) h->Fill(history->cellInfo()->layer());
202 }
203 return h;
204}
205
206
207TH1D* DigitMonitor::residualDist(unsigned int k, const TString& name, int nBins, double rMin, double rMax, bool norm) const
208{
209 TH1D* h = new TH1D(name, "Residual distribution", nBins, rMin, rMax);
210 h->GetXaxis()->SetTitle(norm ? "Residual/ADCMax (%)" : "Residual (ADC counts)");
211
212 for (unsigned int i = 0; i < nChannels(); i++) {
213 const History* history = cellHistory(i);
214 if (!history) continue;
215 for (unsigned int j = 0; j < history->nData(); j++)
216 h->Fill(history->data(j)->delta(k)*(norm ? 100/history->data(j)->adcMax() : 1));
217 }
218 return h;
219}
220
221
222bool DigitMonitor::residualPlots(CaloId calo, unsigned int layer, CaloGain::CaloGain gain, bool xip, bool ring,
223 double lo, double hi, const TString& fileName)
224{
225 TH1D* h1[5];
226 TH2D* hEtaPhi[5];
227 TString globalTitle = TString(xip ? "Derivative" : "Offset") + " correction";
228 for (unsigned short k = 0; k < 5; k++) {
229 TString title = globalTitle + " for sampling " + Form("%d", k);
230 h1[k] = new TH1D(Form("h%d", k), title, 100, lo, hi);
231 h1[k]->GetXaxis()->SetTitle("correction");
232 hEtaPhi[k] = Geo::etaPhiHist(calo, layer, Form("hEtaPhi%d", k), title);
233 hEtaPhi[k]->SetMinimum(lo);
234 hEtaPhi[k]->SetMaximum(hi);
235 }
236 FilterParams f;
237 f.addCalo(calo);
238 f.addLayer(layer);
239 TH2D* hAll = new TH2D("hAll", globalTitle + " vs. sampling", 5, 0, 5, 100, lo, hi);
240 hAll->GetXaxis()->SetTitle("sampling");
241 hAll->GetYaxis()->SetTitle("correction");
242 for (unsigned int i = 0; i < nChannels(); i++) {
243 if ((i+1) % 10000 == 0) cout << "Cell # " << i+1 << "/" << nChannels() << endl;
244 const History* history = pass(i, f);
245 if (!history || history->nData() == 0) continue;
246 for (unsigned short k = 0; k < 5; k++) {
248 double xi = history->data(0)->xi(k, type, gain, xip);
249 h1[k]->Fill(xi);
250 hEtaPhi[k]->Fill(history->cellInfo()->eta(), history->cellInfo()->phi(), xi);
251 if (Id::matchCalo(calo, HEC) && history->cellInfo()->region() == 1) // For HEC, take care of larger phi size in region 1
252 hEtaPhi[k]->Fill(history->cellInfo()->eta(), history->cellInfo()->phi() - Geo::phiSize(HEC_A, 1, 1), xi);
253 }
254 }
255 for (unsigned int i = 1; i <= 100; i++)
256 for (unsigned short k = 0; k < 5; k++)
257 hAll->SetBinContent(k+1, i, h1[k]->GetBinContent(i));
258 TFile* ff = TFile::Open(fileName + ".root", "RECREATE");
259 for (unsigned short k = 0; k < 5; k++) {
260 h1[k]->Write();
261 hEtaPhi[k]->Write();
262 }
263 hAll->Write();
264 ff->Close();
265 TCanvas* c1 = new TCanvas("c1", Id::str(calo) + " layer " + Form("%d", layer), 1600, 2000);
266 c1->Divide(2,3);
267 for (unsigned short k = 0; k < 5; k++) {
268 c1->cd(k+1);
269 hEtaPhi[k]->Draw("COLZ");
270 }
271 c1->cd(6);
272 hAll->Draw();
273 c1->Print(fileName + ".png");
274 c1->Print(fileName + ".eps");
275 return true;
276}
277
278
279bool DigitMonitor::residualPlotsRingComp(CaloId calo, unsigned int layer, CaloGain::CaloGain gain, bool xip,
280 double lo, double hi, const TString& fileName)
281{
282 unsigned int nBins = 100;
283 TH1D* h[5][5];
284 TString globalTitle = TString(xip ? "Derivative" : "Offset") + " correction";
285 for (unsigned short l = 0; l < 5; l++)
286 for (unsigned short k = 0; k < 5; k++)
287 h[l][k] = new TH1D(Form("h%d%d", l, k), "", nBins, l < 2 ? -15 : lo, l < 2 ? 15 : hi);
288 FilterParams f;
289 f.addCalo(calo);
290 f.addLayer(layer);
291 TH2D* hAll[5];
292 hAll[0] = new TH2D("h0All", "#xi_{cell}/#delta#xi vs. sampling", 5, 0, 5, nBins, -15, 15);
293 hAll[1] = new TH2D("hDAll", "(#xi_{cell} - #xi_{ring})/#delta#xi vs. sampling", 5, 0, 5, nBins, -15, 15);
294 hAll[2] = new TH2D("hEAll", globalTitle + " per-cell value vs. sampling", 5, 0, 5, nBins, 0, hi);
295 hAll[3] = new TH2D("hRAll", globalTitle + " per-ring value vs. sampling", 5, 0, 5, nBins, lo, hi);
296 hAll[4] = new TH2D("hCAll", globalTitle + " per-cell error vs. sampling", 5, 0, 5, nBins, lo, hi);
297 TH2D* hCor = new TH2D("hCor", globalTitle + " per-ring vs. per-cell values", nBins, lo, hi, nBins, lo, hi);
298 hCor->GetXaxis()->SetTitle("#xi_{cell}");
299 hCor->GetYaxis()->SetTitle("#xi_{ring}");
300 for (unsigned short l = 0; l < 4; l++) {
301 hAll[l]->GetXaxis()->SetTitle("sampling");
302 hAll[l]->GetYaxis()->SetTitle("significance");
303 }
304 for (unsigned int i = 0; i < nChannels(); i++) {
305 if ((i+1) % 10000 == 0) cout << "Cell # " << i+1 << "/" << nChannels() << endl;
306 const History* history = pass(i, f);
307 if (!history || history->nData() == 0) continue;
308 for (unsigned short k = 0; k < 5; k++) {
309 std::unique_ptr<const ShapeErrorData> cellSED = history->shapeErrorData(gain, CellShapeError);
310 std::unique_ptr<const ShapeErrorData> ringSED = history->shapeErrorData(gain, RingShapeError);
311 if (!cellSED || !ringSED) {
312 continue;
313 }
314 double cellVal = (xip ? cellSED->xi()(k) : cellSED->xip()(k));
315 double ringVal = (xip ? ringSED->xi()(k) : ringSED->xip()(k));
316 double cellErr = TMath::Sqrt(xip ? cellSED->xiErr()(k,k) : cellSED->xipErr()(k,k));
317 h[0][k]->Fill(cellVal/cellErr);
318 h[1][k]->Fill((cellVal - ringVal)/cellErr);
319 h[2][k]->Fill(cellVal);
320 h[3][k]->Fill(ringVal);
321 h[4][k]->Fill(cellErr);
322 hCor->Fill(cellVal, ringVal);
323 }
324 }
325 for (unsigned int i = 1; i <= nBins; i++)
326 for (unsigned short l = 0; l < 5; l++)
327 for (unsigned short k = 0; k < 5; k++)
328 hAll[l]->SetBinContent(k+1, i, h[l][k]->GetBinContent(i));
329 TFile* ff = TFile::Open(fileName + ".root", "RECREATE");
330 for (unsigned short l = 0; l < 5; l++) {
331 for (unsigned short k = 0; k < 5; k++)
332 h[l][k]->Write();
333 hAll[l]->Write();
334 }
335 ff->Close();
336 TCanvas* c1 = new TCanvas("c1", Id::str(calo) + " layer " + Form("%d", layer), 1600, 500);
337 c1->Divide(3,1);
338 c1->cd(1); hCor->Draw();
339 TLatex l; l.SetNDC();
340 TString corrStr = Form("Correlation : %.1f%%", 100*hCor->GetCorrelationFactor());
341 l.DrawLatex(0.2, 0.8, corrStr);
342 c1->cd(2); hAll[0]->Draw("COL");
343 c1->cd(3); hAll[1]->Draw("COL");
344 c1->Print(fileName + ".png");
345 c1->Print(fileName + ".eps");
346 return true;
347}
348
349
350bool DigitMonitor::residualPlotsGainComp(CaloId calo, unsigned int layer, bool ring, bool xip,
351 double lo, double hi, const TString& fileName)
352{
353 unsigned int nBins = 100;
354 TH1D* h[5][5];
355 TString globalTitle = TString(xip ? "Derivative" : "Offset") + " correction";
356 for (unsigned short l = 0; l < 5; l++)
357 for (unsigned short k = 0; k < 5; k++)
358 h[l][k] = new TH1D(Form("h%d%d", l, k), "", nBins, l < 2 ? -15 : lo, l < 2 ? 15 : hi);
359 FilterParams f;
360 f.addCalo(calo);
361 f.addLayer(layer);
362 TH2D* hAll[5];
363 hAll[0] = new TH2D("h0All", "#xi_{high}/#delta#xi vs. sampling", 5, 0, 5, nBins, -15, 15);
364 hAll[1] = new TH2D("hDAll", "(#xi_{high} - #xi_{low})/#delta#xi vs. sampling", 5, 0, 5, nBins, -15, 15);
365 hAll[2] = new TH2D("hEAll", globalTitle + " higher-gain value vs. sampling", 5, 0, 5, nBins, 0, hi);
366 hAll[3] = new TH2D("hRAll", globalTitle + " lower-gain value vs. sampling", 5, 0, 5, nBins, lo, hi);
367 hAll[4] = new TH2D("hCAll", globalTitle + " higher-gain error vs. sampling", 5, 0, 5, nBins, lo, hi);
368 TH2D* hCor = new TH2D("hCor", globalTitle + " lower-gain vs. higher-gain values", nBins, lo, hi, nBins, lo, hi);
369 hCor->GetXaxis()->SetTitle("#xi_{high}");
370 hCor->GetYaxis()->SetTitle("#xi_{low}");
371 for (unsigned short l = 0; l < 4; l++) {
372 hAll[l]->GetXaxis()->SetTitle("sampling");
373 hAll[l]->GetYaxis()->SetTitle("significance");
374 }
375
378
379 for (unsigned int i = 0; i < nChannels(); i++) {
380 if ((i+1) % 10000 == 0) cout << "Cell # " << i+1 << "/" << nChannels() << endl;
381 const History* history = pass(i, f);
382 if (!history || history->nData() == 0) continue;
383 for (unsigned short k = 0; k < 5; k++) {
384 std::unique_ptr<const ShapeErrorData> hgSED = history->shapeErrorData(hiGain, ring ? RingShapeError : CellShapeError);
385 std::unique_ptr<const ShapeErrorData> lgSED = history->shapeErrorData(loGain, ring ? RingShapeError : CellShapeError);
386 if (!hgSED || !lgSED) {
387 continue;
388 }
389 double hgVal = (xip ? hgSED->xi()(k) : hgSED->xip()(k));
390 double lgVal = (xip ? lgSED->xi()(k) : lgSED->xip()(k));
391 double hgErr = TMath::Sqrt(xip ? hgSED->xiErr()(k,k) : hgSED->xipErr()(k,k));
392 double error = TMath::Sqrt(xip ? hgSED->xiErr()(k,k) + lgSED->xiErr()(k,k) : hgSED->xipErr()(k,k) + lgSED->xipErr()(k,k));
393 h[0][k]->Fill(hgVal/hgErr);
394 h[1][k]->Fill((hgVal - lgVal)/error);
395 h[2][k]->Fill(hgVal);
396 h[3][k]->Fill(lgVal);
397 h[4][k]->Fill(hgErr);
398 hCor->Fill(hgVal, lgVal);
399 }
400 }
401 for (unsigned int i = 1; i <= nBins; i++)
402 for (unsigned short l = 0; l < 5; l++)
403 for (unsigned short k = 0; k < 5; k++)
404 hAll[l]->SetBinContent(k+1, i, h[l][k]->GetBinContent(i));
405 TFile* ff = TFile::Open(fileName + ".root", "RECREATE");
406 for (unsigned short l = 0; l < 5; l++) {
407 for (unsigned short k = 0; k < 5; k++)
408 h[l][k]->Write();
409 hAll[l]->Write();
410 }
411 ff->Close();
412 TCanvas* c1 = new TCanvas("c1", Id::str(calo) + " layer " + Form("%d", layer), 1600, 500);
413 c1->Divide(3,1);
414 c1->cd(1); hCor->Draw();
415 TLatex l; l.SetNDC();
416 TString corrStr = Form("Correlation : %.1f%%", 100*hCor->GetCorrelationFactor());
417 l.DrawLatex(0.2, 0.8, corrStr);
418 c1->cd(2); hAll[0]->Draw("COL");
419 c1->cd(3); hAll[1]->Draw("COL");
420 c1->Print(fileName + ".png");
421 c1->Print(fileName + ".eps");
422 return true;
423}
424
425
426TH1D* DigitMonitor::shapeErrorDist(unsigned int k, const TString& name, int nBins, double rMin, double rMax, double mean) const
427{
428 TH1D* h = new TH1D(name, "Shape error", nBins, rMin, rMax);
429 h->GetXaxis()->SetTitle("#epsilon^{2}");
430
431 for (unsigned int i = 0; i < nChannels(); i++) {
432 const History* history = cellHistory(i);
433 if (!history) continue;
434 for (unsigned int j = 0; j < history->nData(); j++) {
435 std::unique_ptr<OFC> ofc = history->ofc(j);
436 double delta = history->data(j)->delta(k) - mean;
437 double adcMax = history->data(j)->adcMax();
438 h->Fill((delta*delta - ofc->Gamma()(k, k))/(adcMax*adcMax));
439 }
440 }
441 cout << sqrt(h->GetMean())*100 << " +/- " << sqrt(h->GetMeanError())*100 << endl;
442 return h;
443}
444
445
446double DigitMonitor::residualCorr(unsigned int k1, unsigned int k2) const
447{
448 double sum = 0;
449 int n = 0;
450
451 for (unsigned int i = 0; i < nChannels(); i++) {
452 const History* history = cellHistory(i);
453 if (!history) continue;
454 for (unsigned int j = 0; j < history->nData(); j++) {
455 double residual1 = history->data(j)->delta(k1);
456 double residual2 = history->data(j)->delta(k2);
457 sum += residual1*residual2;
458 n++;
459 }
460 }
461
462 return n==0 ? 0 : sum/n;
463}
464
465
466bool DigitMonitor::residualParams(int lwb, int upb, CovMatrix& k, TVectorD& means) const
467{
468 CovMatrix sum2N(lwb, upb), sumCN(lwb, upb);
469 TVectorD sum(lwb, upb), sumN(lwb, upb);
470 double sumA = 0, n = 0;
471 for (unsigned int i = 0; i < nChannels(); i++) {
472 const History* history = cellHistory(i);
473 if (!history) continue;
474 for (unsigned int j = 0; j < history->nData(); j++) {
475 if (history->data(j)->adcMax() < 0) continue;
476 std::unique_ptr<OFC> ofc = history->ofc(j);
477 if (!ofc || ofc->hasSameRange(lwb, upb)) { cout << "Invalid index bounds!" << endl; return false; }
478 for (int i1 = lwb; i1 <= upb; i1++) {
479 sum(i1) = sum(i1) + history->data(j)->delta(i1);
480 sumN(i1) = sumN(i1) + history->data(j)->delta(i1)/TMath::Power(history->data(j)->adcMax(), 2);
481 for (int i2 = lwb; i2 <= upb; i2++) {
482 sumCN(i1, i2) = sumCN(i1, i2) + ofc->Gamma()(i1, i2)/TMath::Power(history->data(j)->adcMax(), 2);
483 sum2N(i1, i2) = sum2N(i1, i2) + history->data(j)->delta(i1)*history->data(j)->delta(i2)/TMath::Power(history->data(j)->adcMax(),2);
484 }
485 }
486 sumA += 1/TMath::Power(history->data(j)->adcMax(), 2);
487 n++;
488 }
489 }
490
491 if (n > 0) {
492 sumA /= n;
493 }
494 k.ResizeTo(lwb, upb);
495 means.ResizeTo(lwb, upb);
496
497 for (int i1 = lwb; i1 <= upb; i1++) {
498 if (n > 0) {
499 sum(i1) = sum(i1)/n;
500 sumN(i1) = sumN(i1)/n;
501 }
502 means(i1) = sum(i1);
503 }
504
505 for (int i1 = lwb; i1 <= upb; i1++) {
506 for (int i2 = lwb; i2 <= upb; i2++) {
507 if (n != 0) {
508 sumCN(i1, i2) = sumCN(i1, i2)/n;
509 sum2N(i1, i2) = sum2N(i1, i2)/n;
510 }
511 k(i1, i2) = sum2N(i1, i2) - sumN(i1)*sum(i2) - sum(i1)*sumN(i2)
512 + sum(i1)*sum(i2)*sumA - sumCN(i1,i2);
513 k(i1, i2) = (k(i1, i2) >= 0 ? sqrt(k(i1, i2)) : -sqrt(-k(i1, i2)));
514 }
515 }
516 return true;
517}
518
519
520int DigitMonitor::combine(SimpleShape*& pshape, SimpleShape*& pref, const TString& selection, bool timeAligned) const
521{
522 FilterParams f;
523 if (!f.set(selection)) return 0;
524 int n = 0;
525 pshape = pref = nullptr;
526 //don't change existing interface, use these internally
527 std::unique_ptr<SimpleShape> shape;
528 std::unique_ptr<SimpleShape> ref;
529 double maxSum = 0;
530
531 for (unsigned int i = 0; i < nChannels(); i++) {
532 if ((i+1) % 50000 == 0) cout << "Cell # " << i+1 << "/" << nChannels() << endl;
533 const History* history = pass(i, f);
534 if (!history) continue;
535 for (unsigned int j = 0; j < history->nData(); j++) {
536 if (!f.passEvent(*history->data(j))) continue;
537 maxSum += history->data(j)->maxValue();
538 cout << "Adding pulse (" << n+1 << ") at hash " << i << ", index " << j << ", max = " << maxSum/(n+1) << endl;
539 //cout << i << " " << j << endl;
540 auto thisData = std::make_unique<SimpleShape>(*history->data(j));
541 auto thisRef = history->referenceShape(j);
542 if (timeAligned) {
543 if (!SimpleShape::scaleAndShift(thisData, 1, -history->data(j)->ofcTime())) return -1;
544 if (!SimpleShape::scaleAndShift(thisRef, 1, -history->data(j)->ofcTime())) return -1;
545 }
546 if (!SimpleShape::add(shape, *thisData)) return -1;
547 if (!SimpleShape::add(ref, *thisRef)) return -1;
548 n++;
549 }
550 }
551
552 if (n==0) return -1;
553 if (!SimpleShape::scaleAndShift(shape, 1.0/n)) return -1;
554 if (!SimpleShape::scaleAndShift(ref, 1.0/n)) return -1;
555 //hand these back in existing interface
556 pshape = shape.release();
557 pref = ref.release();
558 return n;
559}
560
561std::unique_ptr<Residuals> DigitMonitor::getResiduals(unsigned int hash, CaloGain::CaloGain gain, double absResTrunc, bool adjust, bool zeroTime) const
562{
563 const History* history = cellHistory(hash);
564 if (!history) return nullptr;
565 std::unique_ptr<const History> historyptr;
566 if (adjust) {
567 historyptr = history->adjust();
568 history = historyptr.get();
569 }
570 return history->residuals(gain, absResTrunc, false, zeroTime);
571}
572
573
574bool DigitMonitor::makeResidualCorrections(const TString& outputFile, short resTrunc, short timeTrunc, double absResTrunc,
575 unsigned int minSize, bool weigh, bool adjust, bool zeroTime) const
576{
577 TreeShapeErrorGetter* shapeError = new TreeShapeErrorGetter(outputFile, true);
578 std::vector< std::vector<ResidualCalculator> > ringCalcs;
579 std::vector<ResidualCalculator> gainRingCalcs(Geo::nPhiRings());
580 for (unsigned int g = 0; g < CaloGain::LARNGAIN; g++) ringCalcs.push_back(gainRingCalcs);
581
582 cout << "Processing cells" << endl;
583
584 for (unsigned int i = 0; i < nChannels(); i++) {
585 if ((i+1) % 10000 == 0) cout << "Processing hash = " << (i+1) << endl;
586 for (unsigned int g = 0; g < CaloGain::LARNGAIN; g++) {
587 std::unique_ptr<Residuals> residuals = getResiduals(i, (CaloGain::CaloGain)g, absResTrunc, adjust, zeroTime);
588 if (residuals && (resTrunc > 0 || timeTrunc > 0)) {
589 residuals = residuals->truncate(resTrunc, timeTrunc);
590 }
591
592 if (residuals && residuals->size() < minSize) {
593 residuals = nullptr;
594 }
595
596 ResidualCalculator* resCalc = nullptr;
597 if (residuals) {
598 resCalc = residuals->calculator(weigh);
599 }
600
601 if (!resCalc) {
603 continue;
604 }
605
606 cout << i << " : Phi ring " << cellHistory(i)->cellInfo()->globalPhiRing() << " : adding " << resCalc->size() << endl;
607 shapeError->addCell(*resCalc, (CaloGain::CaloGain)g);
608 short iring = cellHistory(i)->cellInfo()->globalPhiRing();
609 if (iring >= 0) {
610 ringCalcs[g][iring].append(*resCalc);
611 if (iring == 1142) cout << ringCalcs[g][iring].regresser()->mean(0) << endl;
612 }
613 delete resCalc;
614 }
615 }
616
617 cout << "Processing phi-symmetric data" << endl;
618
619 for (unsigned int g = 0; g < CaloGain::LARNGAIN; g++) {
620 for (short i = 0; i < Geo::nPhiRings(); i++) {
621 if (ringCalcs[g][i].size() > 0)
622 cout << "Saving data for gain " << g << ", phi ring " << i << ", n = " << ringCalcs[g][i].size() << ", range = " << ringCalcs[g][i].rangeStr() << endl;
623 shapeError->addRing(ringCalcs[g][i], (CaloGain::CaloGain)g);
624 }
625 }
626
627 cout << "Done!" << endl;
628
629 delete shapeError;
630 return true;
631}
const boost::regex ref(r_ef)
static Double_t * pref
#define max(a, b)
Definition cfImp.cxx:41
const History * pass(unsigned int i, const FilterParams &f) const
virtual const History * cellHistory(unsigned int i) const
double maxValue(bool withErrors=false) const
Definition AbsShape.cxx:30
double eta() const
Definition CellInfo.h:93
short layer() const
Definition CellInfo.h:53
short region() const
Definition CellInfo.h:59
double phi() const
Definition CellInfo.h:94
short globalPhiRing() const
Definition CellInfo.cxx:125
void add(const AbsShapeErrorGetter &getter)
double delta(short sample) const
Definition Data.cxx:430
CaloGain::CaloGain gain() const
Definition Data.h:87
double _adcMax(const DataFuncArgs &) const
Definition Data.h:198
double _minValue(const DataFuncArgs &) const
Definition Data.h:202
double _energy(const DataFuncArgs &) const
Definition Data.h:193
double adcMax() const
Definition Data.h:126
double _maxValue(const DataFuncArgs &) const
Definition Data.h:201
double xi(short sample, ShapeErrorType shapeErrorType=BestShapeError, CaloGain::CaloGain g=CaloGain::UNKNOWNGAIN, bool xip=false) const
Definition Data.cxx:462
double ofcTime() const
Definition Data.h:113
double _noise(const DataFuncArgs &) const
Definition Data.h:200
double _ofcTime(const DataFuncArgs &) const
Definition Data.h:194
TH1D * timeDist(const TString &name, int nBins, double tMin=-25, double tMax=25) const
TH1D * residualDist(unsigned int k, const TString &name, int nBins, double rMin, double rMax, bool norm=false) const
int combine(SimpleShape *&shape, SimpleShape *&ref, const TString &selection="", bool timeAligned=true) const
TH1D * shapeErrorDist(unsigned int k, const TString &name, int nBins, double rMin, double rMax, double mean=0) const
double residualCorr(unsigned int k1, unsigned int k2) const
TH2D * minValueMap(TString name, PartitionId partition) const
TH1D * gainDist(const TString &name) const
bool residualParams(int lwb, int upb, CovMatrix &k, TVectorD &means) const
TH1D * noiseDist(const TString &name, int nBins, double max) const
TH1D * chi2Dist(const TString &name, int nBins, double max, double kFactor=0, double fitMax=-1, int lwb=-1, int upb=-1, unsigned int chi2Pars=DefaultChi2, ShapeErrorType shapeErrorType=BestShapeError, unsigned int nDof=0) const
std::unique_ptr< Residuals > getResiduals(unsigned int hash, CaloGain::CaloGain gain, double absResTrunc=-1, bool adjust=false, bool zeroTime=false) const
bool residualPlotsGainComp(CaloId calo, unsigned int layer, bool ring=false, bool xip=false, double lo=-0.05, double hi=0.05, const TString &fileName="residuals")
static TF1 * fitChi2(TH1D &h, const char *name, double xMin, double xMax, double nDof=-1, double fitMin=Definitions::none, double fitMax=Definitions::none)
TH2D * maxValueMap(TString name, PartitionId partition) const
bool residualPlotsRingComp(CaloId calo, unsigned int layer, CaloGain::CaloGain gain=CaloGain::LARHIGHGAIN, bool xip=false, double lo=-0.05, double hi=0.05, const TString &fileName="residuals")
bool makeResidualCorrections(const TString &outputFile, short resTrunc=-1, short timeTrunc=-1, double absResTrunc=-1, unsigned int minSize=0, bool weigh=false, bool adjust=false, bool zeroTime=false) const
static TF1 * chi2Func(const char *name, double xMin, double xMax)
TH1D * energyDist(const TString &name, int nBins, double eMax=10000) const
TH1D * bestChi2Dist(const TString &name, int nBins, double max, double refErrMin, double refErrMax, unsigned int refErrNBins, int lwb=-1, int upb=-1, unsigned int chi2Pars=DefaultChi2, unsigned int nDof=0) const
TH1D * layerDist(const TString &name) const
TH1D * adcMaxDist(const TString &name, int nBins, double aMin=0, double aMax=4096) const
bool residualPlots(CaloId calo, unsigned int layer, CaloGain::CaloGain gain=CaloGain::LARHIGHGAIN, bool xip=false, bool ring=false, double lo=-0.05, double hi=0.05, const TString &fileName="residuals")
TVectorD means(int lwb, int upb)
static TH2D * etaPhiHist(CaloId calo, short layer, const TString &name, const TString &title)
Definition Geometry.cxx:549
static double phiSize(CaloId calo, short layer, short region, short iPhi=-1)
Definition Geometry.cxx:344
static short nPhiRings()
Definition Geometry.cxx:84
double chi2(int i, int lwb=-1, int upb=-1, int chi2Params=DefaultChi2, ShapeErrorType shapeErrorType=BestShapeError, unsigned int *nDof=0) const
Definition History.cxx:161
std::unique_ptr< const ShapeErrorData > shapeErrorData(CaloGain::CaloGain gain, ShapeErrorType shapeErrorType=BestShapeError, const Residual *res=0) const
Definition History.cxx:293
std::unique_ptr< SimpleShape > referenceShape(unsigned int k, double adcMax=-1, double time=Definitions::none, bool samplesOnly=false) const
Definition History.cxx:527
const CellInfo * cellInfo() const
Definition History.h:59
const Data * data(unsigned int i) const
Definition History.cxx:88
void setShapeErrorGetter(const AbsShapeErrorGetter *err) const
Definition History.h:93
unsigned int nData() const
Definition History.h:54
std::unique_ptr< OFC > ofc(unsigned int i, int lwb=-1, int upb=-1, double time=Definitions::none, bool useCorrs=true) const
Definition History.cxx:201
std::unique_ptr< Residuals > residuals(CaloGain::CaloGain gain=CaloGain::LARNGAIN, double absResTrunc=-1, bool correct=true, bool zeroTime=false) const
Definition History.cxx:593
std::unique_ptr< History > adjust() const
Definition History.cxx:250
static TString str(CaloId id)
Definition CaloId.cxx:15
static bool matchCalo(CaloId id, CaloId idSpec)
Definition CaloId.cxx:188
virtual unsigned int nChannels() const override
TH1D * dist(const DataFuncSet &func, const DataFuncArgs &args, const TString &name, int nBins, double xMin, double xMax, const TString &title="", const TString &xTitle="", const TString &yTitle="", const FilterParams &f=FilterParams()) const
TH2D * partitionMap(const DataFuncSet &func, const DataFuncArgs &args, TString name, PartitionId partition, const TString &title="", CombinationType comb=AverageValue, const FilterParams &f=FilterParams()) const
bool add(unsigned int k, double value, double error)
static bool scaleAndShift(std::unique_ptr< SimpleShape > &s1, double scale, double shift=0)
int addCell(const ResidualCalculator &calc, CaloGain::CaloGain gain)
int addRing(const ResidualCalculator &calc, CaloGain::CaloGain gain)
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="")
const std::string selection
@ LARMEDIUMGAIN
Definition CaloGain.h:18
@ LARLOWGAIN
Definition CaloGain.h:18
@ LARNGAIN
Definition CaloGain.h:19
@ LARHIGHGAIN
Definition CaloGain.h:18