R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandHitCalibrationBar.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
15#include "R3BLogger.h"
17#include "R3BNeulandHitPar.h"
18
19#include "FairLogger.h"
20
21#include "TCanvas.h"
22#include "TDirectory.h"
23#include "TGaxis.h"
24#include "TGraph.h"
25#include "TPad.h"
26
27#include <limits>
28#include <numeric>
29#include <string>
30
31using DPair = std::array<Double_t, 2>;
32
33// FIXME namespace Neuland::Calibration { with c++17
34namespace R3B::Neuland // NOLINT
35{
36 namespace Calibration
37 {
38 bool GetJumps(double* data,
39 unsigned int size,
40 double_t threshold,
41 unsigned int averageSize,
42 std::vector<unsigned int>* output)
43 {
44 if (size < 2 * averageSize)
45 return false;
46
47 const auto scale = 1. / averageSize;
48
49 auto lastAverage = std::accumulate(data, data + averageSize, 0.) * scale;
50 const auto last = data + size;
51 auto currentPos = data + averageSize;
52
53 while (currentPos + averageSize < last)
54 {
55 const auto currentAverage = std::accumulate(currentPos, currentPos + averageSize, 0.) * scale;
56 if (fabs(currentAverage - lastAverage) > threshold)
57 {
58 if (output == nullptr)
59 return true;
60 output->push_back(currentPos - data);
61 }
62 lastAverage = currentAverage;
63 currentPos += averageSize;
64 }
65
66 return false;
67 }
68
69 template <Int_t iterations = 8>
70 const double FastExp(const double val)
71 {
72 auto exp = 1. + val / (1 << iterations);
73 for (auto i = 0; i < iterations; ++i)
74 {
75 exp *= exp;
76 }
77 return exp;
78 }
79
80 const double GetLightAtPMT(const double energyLoss,
81 const double position,
82 const double invLightAttenuationLength,
83 const int side)
84 {
85 return energyLoss *
86 FastExp(-invLightAttenuationLength * (TotalBarLength * 0.5 + (1 - 2 * side) * position));
87 }
88
89 constexpr auto RightSide = 0;
90 constexpr auto LeftSide = 1;
91
92 constexpr auto RightSideHitBit = RightSide;
93 constexpr auto LeftSideHitBit = LeftSide;
94 constexpr auto PosCalibrationBit = 2;
95 constexpr auto TSyncCalibrationBit = 3;
96 constexpr auto EnergyCalibrationBit = 4;
97 constexpr auto ThresholdCalibrationBit = 5;
98 constexpr auto PedestalCalibrationBit = 6;
99 constexpr auto TimeJumpedBit = 20;
100 constexpr auto LogCompleteBit = 30;
101 constexpr auto StrangeBit = 31;
102
103 constexpr auto PositionCalibrationSize = 1024;
104 constexpr auto EnergyCalibrationSize = 4 * 1024;
105 constexpr auto ThresholdCalibrationSize = 1024;
106 constexpr auto CalibrationLogSize = 16 * 1024;
107 constexpr auto LogInitialSize = 128;
108
109 constexpr auto MinEnergyDeposit = 1.; // MeV
110 constexpr auto MaxTSyncError = 0.05; // ns
111 constexpr auto TimeJumpThreshold = 0.075; // ns
112 constexpr auto MaxFastTDiffError = 0.05;
113
114 constexpr auto MaxNumberOfFails = 10U;
115
116 const void SetStatus(Int_t& var, Int_t statusBit) { var |= (1 << statusBit); }
117 const void ClearStatus(Int_t& var, Int_t statusBit) { var &= ~(1 << statusBit); }
118 const bool IsStatus(Int_t var, Int_t statusBit) { return (var & (1 << statusBit)); }
119
120 const char* HitCalibrationBar::CalibrationStatusDescription[] = { "No Data recorded",
121 "Calibration Failed",
122 "Calibration was Successful",
123 "Energy Calibration Failed",
124 "Time Synchronisation Failed",
125 "Time Synchronisation Error large",
126 "Parameters look strange",
127 "Time Jump" };
128 const char* HitCalibrationBar::CalibrationStatusAbbreviation[] = { "XX", "FA", " ", "EC",
129 "TS", "SE", "??", "TJ" };
130
132 : ID(id)
133 , Validity(0)
134 , FailCounter(0)
135 , LastEventNumber(0)
136 , TimeDifference(NaN)
137 , EffectiveSpeed(NaN)
138 , Gain({ NaN, NaN })
139 , InvLightAttenuationLength(NaN)
140 , Pedestal({ 0, 0 })
141 , Saturation({ 0., 0. })
142 , TimeSync({ NaN, NaN })
143 , Threshold({ NaN, NaN })
144 {
145
146 Reset();
147 LastHits.reserve(CalibrationLogSize);
148
149 Log.TimeDifference.Expand(LogInitialSize);
150 Log.TimeDifference.SetTitle("Time Difference Parameter;Event Number;Parameter Value / ns");
151
152 Log.EffectiveSpeed.Expand(LogInitialSize);
153 Log.EffectiveSpeed.SetTitle("Effective Speed Parameter;Event Number;Parameter Value / cm/ns");
154
155 Log.LightAttenuationLength.Expand(LogInitialSize);
156 Log.LightAttenuationLength.SetTitle("Attenuation Length Parameter;Event Number;Parameter Value / cm");
157
158 for (auto side = 0; side < 2; ++side)
159 {
160 Log.Gain[side].Expand(LogInitialSize);
161 Log.Gain[side].SetTitle("Gain Parameter;Event Number;Parameter Value / ns/MeVle");
162
163 Log.TotalHits[side] = TH1F("", "Total Hits; Light at PMT / MeVle; #Hits", 64, 0., 32.);
164 Log.TotalHits[side].SetDirectory(nullptr);
165 Log.TotalHits[side].SetStats(false);
166
167 Log.MissedHits[side] = TH1F("", "Missed Hits; Light at PMT / MeVle; #Hits", 64, 0., 32.);
168 Log.MissedHits[side].SetDirectory(nullptr);
169 Log.MissedHits[side].SetStats(false);
170 Log.MissedHits[side].SetLineColor(kGreen);
171
172 Log.Energy[side] = TH1F("", "Energy Calibration; Gain / ns/MeVle; Counts", 100, 0., 30.);
173 Log.Energy[side].SetDirectory(nullptr);
174 Log.Energy[side].SetStats(false);
175
176 Log.Saturation[side] =
177 TH2F("", "Saturation; Light at PMT / MeVle; QDC / ns", 128, 0., 128., 256, 0., 512.);
178 Log.Saturation[side].SetDirectory(nullptr);
179 Log.Saturation[side].SetStats(false);
180
181 Log.Pedestal[side] = TH1F("", "Pedestal; QDC / ns", 512, 0., 512.);
182 Log.Pedestal[side].SetDirectory(nullptr);
183 Log.Pedestal[side].SetStats(false);
184 }
185
186 FitGraph.Expand(PositionCalibrationSize);
187 }
188
190
192 {
193 if (!std::isnan(par->GetTimeOffset(1)))
194 {
195 TimeDifference = par->GetTDiff();
196 EffectiveSpeed = par->GetEffectiveSpeed();
197 SetStatus(Validity, PosCalibrationBit);
198 }
199
200 if (!std::isnan(par->GetEnergyGain(1)))
201 {
202 Gain[0] = par->GetEnergyGain(1);
203 Gain[1] = par->GetEnergyGain(2);
204
205 InvLightAttenuationLength = 1. / par->GetLightAttenuationLength();
206
208 }
209
210 if (!std::isnan(par->GetPMTThreshold(1)))
211 {
212 Threshold[0] = par->GetPMTThreshold(1);
213 Threshold[1] = par->GetPMTThreshold(2);
215 }
216
217 if (!std::isnan(par->GetPMTSaturation(1)))
218 {
219 Saturation[0] = par->GetPMTSaturation(1);
220 Saturation[1] = par->GetPMTSaturation(2);
221 }
222
223 if (!std::isnan(par->GetPedestal(1)))
224 {
225 Pedestal[0] = par->GetPedestal(1);
226 Pedestal[1] = par->GetPedestal(2);
228 }
229
230 LOG(debug) << "Loaded Parameters for Bar " << ID;
231 return;
232 }
233
234 void HitCalibrationBar::Set(const Int_t side, const Double_t time, const Int_t qdc)
235 {
236 SetStatus(Validity, side);
237 R3BLOG(debug1, fmt::format("side: {}, time: {} qdc: {}", side, time, qdc));
238 CurrentHit.Time[side] = time;
239 CurrentHit.QDC[side] = std::max(qdc - Pedestal[side], 1);
240 Log.Pedestal[side].Fill(qdc);
241
242 if (IsStatus(Validity, LeftSideHitBit) && IsStatus(Validity, RightSideHitBit))
243 {
244 if (fabs(CurrentHit.Time[1] - CurrentHit.Time[0]) > 0.5 * MaxCalTime)
245 {
246 if (CurrentHit.Time[1] < CurrentHit.Time[0])
247 CurrentHit.Time[1] += MaxCalTime;
248 else
249 CurrentHit.Time[0] += MaxCalTime;
250 }
251
252 if (IsStatus(Validity, PosCalibrationBit))
253 {
254 CurrentHit.EntryPosition =
255 (CurrentHit.Time[1] - CurrentHit.Time[0] - TimeDifference) * EffectiveSpeed;
256
257 if (fabs(CurrentHit.EntryPosition) > TotalBarLength)
258 {
259 ++FailCounter;
260 if (FailCounter > 10)
261 {
262 LOG(debug) << "Bar " << ID << ":"
263 << " Positions were too often too strange. Might be a Time-Jump in one PMT. "
264 "Removing Position-Calibration.";
265 FailCounter = 0;
267 }
268 }
269 else
270 {
271 FailCounter = 0;
272 }
273 }
274 }
275 }
276 void HitCalibrationBar::SetGlobalTSync(const Double_t value, const Double_t error)
277 {
278 TimeSync = { value, error };
279 if (std::isnan(value))
281 else
283 }
284
286 {
287 if (LastHits.size() == 0 && !IsStatus(Validity, LogCompleteBit))
288 return CalibrationStatus::noData;
289 if (!IsStatus(Validity, PosCalibrationBit))
290 return CalibrationStatus::completeFail;
291 if (!IsStatus(Validity, EnergyCalibrationBit))
292 return CalibrationStatus::energyGainFail;
293 if (!IsStatus(Validity, TSyncCalibrationBit))
294 return CalibrationStatus::timeSyncFail;
295 if (IsStatus(Validity, TimeJumpedBit))
296 return CalibrationStatus::timeJump;
297 if (TimeSync[1] > MaxTSyncError)
298 return CalibrationStatus::timeSyncError;
299 if (IsStatus(Validity, StrangeBit))
300 return CalibrationStatus::strangeParameters;
301 return CalibrationStatus::success;
302 }
303
304 Bool_t HitCalibrationBar::Add(const Double_t timeOffset,
305 const Double_t entryPosition,
306 const Double_t exitPosition,
307 const Double_t energy,
308 const UInt_t eventNumber)
309 {
310
311 // If we have already a gain calibration, we can log the hits.
312 // This way we might get some information about strange behaviour.
313 if (IsStatus(Validity, EnergyCalibrationBit))
314 {
315 const auto centerPosition = (entryPosition + exitPosition) * 0.5;
316 const DPair lightAtPMT = { GetLightAtPMT(energy, centerPosition, InvLightAttenuationLength, RightSide),
317 GetLightAtPMT(energy, centerPosition, InvLightAttenuationLength, LeftSide) };
318
319 if (!IsStatus(Validity, RightSideHitBit))
320 {
321 Log.MissedHits[RightSide].Fill(lightAtPMT[RightSide]);
322 }
323 if (!IsStatus(Validity, LeftSideHitBit))
324 {
325 Log.MissedHits[LeftSide].Fill(lightAtPMT[LeftSide]);
326 }
327
328 for (auto side = 0; side < 2; ++side)
329 {
330 Log.TotalHits[side].Fill(lightAtPMT[side]);
331 // std::cout << "???????????????lightAtPMT" << lightAtPMT[side] << " \n";
332 }
333 if (IsStatus(Validity, LeftSideHitBit) && IsStatus(Validity, RightSideHitBit))
334 {
335 for (auto side = 0; side < 2; ++side)
336 {
337 Log.Saturation[side].Fill(lightAtPMT[side], CurrentHit.QDC[side]);
338 const auto gain = CurrentHit.QDC[side] /
339 (lightAtPMT[side] * (1. - CurrentHit.QDC[side] * SaturationCoefficient));
340 Log.Energy[side].Fill(gain);
341 }
342 }
343 }
344
345 // We require both sides and reasonable energy deposit.
346 if (!IsStatus(Validity, LeftSideHitBit) || !IsStatus(Validity, RightSideHitBit) ||
347 energy < MinEnergyDeposit)
348 return false;
349
350 // Dont take the bar, if the cosmic was somewhere close to the lightguide:
351 // Strange things happen there ;).
352 const auto maxPos = BarLength * 0.5 - 1.;
353 if (fabs(entryPosition) > maxPos || fabs(exitPosition) > maxPos)
354 return false;
355
356 LastEventNumber = eventNumber;
357
358 for (auto side = 0; side < 2; ++side)
359 CurrentHit.Time[side] -= timeOffset * 0.5;
360
361 CurrentHit.EntryPosition = entryPosition;
362 CurrentHit.ExitPosition = exitPosition;
363 CurrentHit.Energy = energy;
364
365 LastHits.push_back(CurrentHit);
366
367 if (LastHits.size() % PositionCalibrationSize == 0)
368 positionCalibration(LastHits.size() - PositionCalibrationSize, PositionCalibrationSize);
369
370 if (LastHits.size() % EnergyCalibrationSize == 0)
371 energyCalibration(LastHits.size() - EnergyCalibrationSize, EnergyCalibrationSize);
372
373 // If we reached our limit, clear the vector.
374 // The content itself is still there. We will set a (hacky) bit,
375 // so we know when the vector actually contains more data than its size.
376 if (LastHits.size() == CalibrationLogSize)
377 {
378 SetStatus(Validity, LogCompleteBit);
379 LastHits.clear();
380 }
381
382 // We can say only that this was a somewhat good event, if we have a position calibration
383 // and the position of the track matches with the one from the calibration
384 if (IsStatus(Validity, PosCalibrationBit) &&
385 fabs(0.5 * (entryPosition + exitPosition) - GetPosition()) > 25.)
386 return true;
387
388 return true;
389 }
390
392 {
393 ClearStatus(Validity, LeftSideHitBit);
394 ClearStatus(Validity, RightSideHitBit);
395 CurrentHit.EntryPosition = NaN;
396 }
397
398 void HitCalibrationBar::Calibrate(TDirectory* histogramDir)
399 {
400 const auto nHits = (IsStatus(Validity, LogCompleteBit) ? CalibrationLogSize : LastHits.size());
401
402 // Check if we have any Hits at all.
403 // If not, this bar is either not connected or has some problems.
404 // In this case do not bother creating plots and just return.
405 if (nHits == 0)
406 return;
407
408 // Position Calibration
409 if (Log.TimeDifference.GetN() > 0)
410 {
411 TimeDifference = getMean(Log.TimeDifference);
412
413 EffectiveSpeed = getMean(Log.EffectiveSpeed);
414
415 // check for timejumps
416 if (GetJumps(Log.TimeDifference.GetY(), Log.TimeDifference.GetN(), TimeJumpThreshold, 3, nullptr))
417 SetStatus(Validity, TimeJumpedBit);
418
419 SetStatus(Validity, PosCalibrationBit);
420 }
421 else if (nHits > 0.5 * PositionCalibrationSize)
422 {
423 // There were not enough hits for the usual calibration.
424 // If we have at least half the chunk size, we can try to calibrate anyway
425 positionCalibration(0, nHits);
426 }
427 else
428 {
429 // We cannot have a reliable calibration.
430 // Give up at this point.
432 }
433
434 // Gain Calibration
435 if (Log.Gain[0].GetN() > 0 && Log.Gain[1].GetN() > 0)
436 {
437 for (auto side = 0; side < 2; ++side)
438 {
439 Gain[side] = getMean(Log.Gain[side], Gain[side]);
440 Saturation[side] = SaturationCoefficient * Gain[side];
441 }
442
443 InvLightAttenuationLength = 1. / getMean(Log.LightAttenuationLength, 1. / InvLightAttenuationLength);
444
446 }
447 else if (nHits > 0.5 * EnergyCalibrationSize)
448 {
449 // There were not enough hits for the usual calibration.
450 // If we have at least half the chunk size, we can try to calibrate anyway
451 energyCalibration(0, nHits);
452 }
453 else
454 {
455 // We cannot have a reliable calibration.
456 // Give up at this point.
458 }
459
460 // Threshold Calibration
461 thresholdCalibration();
462
463 createHistograms(histogramDir);
464 }
465
467 {
468 R3BNeulandHitModulePar parameter;
469 parameter.SetModuleId(ID);
470
471 parameter.SetTDiff(TimeDifference);
472 parameter.SetTSync(TimeSync[0]);
473 parameter.SetEffectiveSpeed(EffectiveSpeed);
474 parameter.SetLightAttenuationLength(1. / InvLightAttenuationLength);
475
476 for (auto side = 0; side < 2; ++side)
477 {
478 parameter.SetEnergyGain(Gain[side], side + 1);
479 parameter.SetPMTThreshold(Threshold[side], side + 1);
480 parameter.SetPMTSaturation(Saturation[side], side + 1);
481 parameter.SetPedestal(Pedestal[side], side + 1);
482 }
483
484 LOG(debug) << ID << " Bar Parameters: ";
485 LOG(debug) << " Time Offset : " << parameter.GetTimeOffset(1) << " " << parameter.GetTimeOffset(2);
486 LOG(debug) << " Effective Speed : " << parameter.GetEffectiveSpeed();
487 LOG(debug) << " Energy Gain : " << parameter.GetEnergyGain(1) << " " << parameter.GetEnergyGain(2);
488 LOG(debug) << " Attenuationlength : " << parameter.GetLightAttenuationLength();
489 LOG(debug) << " Threshold : " << parameter.GetPMTThreshold(1) << " "
490 << parameter.GetEnergyGain(2);
491 LOG(debug) << " Pedestal : " << parameter.GetPedestal(1) << " " << parameter.GetPedestal(2);
492 LOG(debug) << " Saturation : " << parameter.GetPMTSaturation(1) << " "
493 << parameter.GetPMTSaturation(2);
494
495 return parameter;
496 }
497
499 {
500 return IsStatus(Validity, LeftSideHitBit) && IsStatus(Validity, RightSideHitBit);
501 }
502
503 void HitCalibrationBar::positionCalibration(int firstHit, int nHits)
504 {
505 FitGraph.Set(0);
506 const auto loopSize = (LastHits.size() >= (firstHit + nHits)) ? nHits : (LastHits.size() - firstHit);
507 for (auto index = 0; index < loopSize; ++index)
508 {
509 const auto& hit = LastHits[firstHit + index];
510
511 const auto meanPosition = 0.5 * (hit.EntryPosition + hit.ExitPosition);
512 const auto tdiff = hit.Time[1] - hit.Time[0];
513 FitGraph.SetPoint(index, meanPosition, tdiff);
514 }
515
516 TF1 linearFit("", "pol1", -TotalBarLength * 0.5, TotalBarLength * 0.5);
517 FitGraph.Fit(&linearFit, "NQ");
518 if (linearFit.GetParError(0) > MaxFastTDiffError)
519 cleanupFit(FitGraph, linearFit, 2.5);
520
521 // Use this calibration
522 TimeDifference = linearFit.GetParameter(0);
523 EffectiveSpeed = 1. / linearFit.GetParameter(1);
524
525 // Write parameters to the Log.
526 const auto nPoints = Log.TimeDifference.GetN();
527 Log.TimeDifference.SetPoint(nPoints, LastEventNumber, TimeDifference);
528 Log.TimeDifference.SetPointError(nPoints, 0, linearFit.GetParError(0));
529 Log.EffectiveSpeed.SetPoint(nPoints, LastEventNumber, EffectiveSpeed);
530 Log.EffectiveSpeed.SetPointError(nPoints, 0, linearFit.GetParError(1) * Sqr(EffectiveSpeed));
531
532 SetStatus(Validity, PosCalibrationBit);
533 }
534
535 void HitCalibrationBar::energyCalibration(int firstHit, int nHits)
536 {
537 const auto loopSize = (LastHits.size() >= (firstHit + nHits)) ? nHits : (LastHits.size() - firstHit);
538 if (!IsStatus(Validity, PedestalCalibrationBit))
539 {
540 // This is the first calibration
541 // if we do not have pedestals yet, get them and subtract them from the hits
542
543 pedestalCalibration();
544 for (auto h_index = 0; h_index < loopSize; ++h_index)
545 {
546 auto& hit = LastHits[firstHit + h_index];
547 for (auto side = 0; side < 2; ++side)
548 {
549 hit.QDC[side] = std::max(hit.QDC[side] - Pedestal[side], 1);
550 }
551 }
552 }
553
554 FitGraph.Set(0);
555 for (auto h_index = 0; h_index < loopSize; ++h_index)
556 {
557 const auto& hit = LastHits[firstHit + h_index];
558 const auto centerPosition = (hit.EntryPosition + hit.ExitPosition) * 0.5;
559 FitGraph.SetPoint(h_index, centerPosition, log(hit.QDC[1] * 1. / hit.QDC[0]));
560 }
561
562 TF1 linearFit("", "pol1", -TotalBarLength * 0.5, TotalBarLength * 0.5);
563 // this has usually always outliers, so do not bother doing a normal fit first.
564 cleanupFit(FitGraph, linearFit, 1.);
565
566 InvLightAttenuationLength = 0.5 * linearFit.GetParameter(1);
567
568 const auto logLightAttLenPoints = Log.LightAttenuationLength.GetN();
569 Log.LightAttenuationLength.SetPoint(logLightAttLenPoints, LastEventNumber, 1. / InvLightAttenuationLength);
570 Log.LightAttenuationLength.SetPointError(
571 logLightAttLenPoints, 0, 0.5 * linearFit.GetParError(1) / Sqr(InvLightAttenuationLength));
572
573 for (auto h_index = 0; h_index < loopSize; ++h_index)
574 {
575 const auto& hit = LastHits[firstHit + h_index];
576 const auto centerPosition = (hit.EntryPosition + hit.ExitPosition) * 0.5;
577 for (auto side = 0; side < 2; ++side)
578 {
579 const auto lightAtPMT = GetLightAtPMT(hit.Energy, centerPosition, InvLightAttenuationLength, side);
580 const auto gain = hit.QDC[side] / (lightAtPMT * (1. - hit.QDC[side] * SaturationCoefficient));
581 Log.Energy[side].Fill(gain);
582 }
583 }
584
585 std::array<TH1F, 2> hGainCal = { TH1F("", "", 200, 0., 40.), TH1F("", "", 200, 0., 40.) };
586
587 for (auto h_index = 0; h_index < loopSize; ++h_index)
588 {
589 const auto& hit = LastHits[firstHit + h_index];
590 const auto centerPosition = (hit.EntryPosition + hit.ExitPosition) * 0.5;
591 for (auto side = 0; side < 2; ++side)
592 {
593 const auto lightAtPMT = GetLightAtPMT(hit.Energy, centerPosition, InvLightAttenuationLength, side);
594 const auto distance = 0.5 * TotalBarLength + (1 - 2 * side) * centerPosition;
595 const auto gain = hit.QDC[side] / (lightAtPMT * (1. - hit.QDC[side] * SaturationCoefficient));
596 hGainCal[side].Fill(gain);
597 }
598 }
599
600 TF1 gausFit("", "gaus", 0., 30.);
601 for (auto side = 0; side < 2; ++side)
602 {
603 const auto maxPos = hGainCal[side].GetBinCenter(hGainCal[side].GetMaximumBin());
604 gausFit.SetParameter(0, hGainCal[side].GetMaximum());
605 gausFit.SetParameter(1, maxPos);
606 hGainCal[side].Fit(&gausFit, "NQ", "", maxPos - 1.5, maxPos + 1.5);
607
608 Gain[side] = gausFit.GetParameter(1);
609
610 const auto nPoints = Log.Gain[side].GetN();
611 Log.Gain[side].SetPoint(nPoints, LastEventNumber, Gain[side]);
612 Log.Gain[side].SetPointError(nPoints, 0, gausFit.GetParError(1));
613 }
614
616 }
617
618 void HitCalibrationBar::pedestalCalibration()
619 {
620 for (auto side = 0; side < 2; ++side)
621 for (auto bin = 0; bin < Log.Pedestal[side].GetNbinsX(); ++bin)
622 {
623 if (Log.Pedestal[side].GetBinContent(bin + 1) > 3.)
624 {
625 Pedestal[side] = bin;
626 break;
627 }
628 }
629
631 }
632
633 void HitCalibrationBar::thresholdCalibration()
634 {
635 if (Log.TotalHits[0].GetEntries() < ThresholdCalibrationSize ||
636 Log.TotalHits[1].GetEntries() < ThresholdCalibrationSize)
637 return;
638
639 // We have 50% of the maximum at x=[1]
640 // For the bar at the edges it is hard to check if the cosmic was stopped.
641 // Therefore we usually have [2] != 0.
642 TF1 missFit("missFit", "(1. - [2]) * 0.5 * (1. - TMath::Erf([0] * (x - [1]))) + [2]", 0., 20.);
643
644 for (auto side = 0; side < 2; ++side)
645 {
646 const auto nPoints = Log.TotalHits[side].GetNbinsX();
647 Log.MissRatio[side].Expand(nPoints);
648 Log.MissRatio[side].Set(0);
649 for (auto bin = 1; bin <= nPoints; ++bin)
650 {
651 if (Log.MissedHits[side].GetBinContent(bin) > 0)
652 {
653 const auto nextPoint = Log.MissRatio[side].GetN();
654 const auto binCenter = Log.MissedHits[side].GetBinCenter(bin);
655 const auto ratio =
656 Log.MissedHits[side].GetBinContent(bin) / Log.TotalHits[side].GetBinContent(bin);
657 Log.MissRatio[side].SetPoint(nextPoint, binCenter, ratio);
658 }
659 }
660
661 if (Log.MissRatio[side].GetN() > 5)
662 {
663 missFit.SetParameter(0, 0.25);
664 missFit.SetParameter(1, AvgThreshold);
665 missFit.SetParLimits(1, 0., 10.);
666 missFit.SetParameter(2, 0.);
667 missFit.SetParLimits(2, 0., 0.5);
668 Log.MissRatio[side].Fit(&missFit, "QNE");
669 Threshold[side] = missFit.GetParameter(1);
670 }
671 }
673 }
674
675 Int_t HitCalibrationBar::cleanupFit(TGraph& graph, TF1& fit, Double_t maxDifference) const
676 {
677 std::array<Int_t, 256> remove;
678 auto totalRemoved = 0;
679
680 graph.Fit(&fit, "QN ROB=0.90");
681
682 auto removePosition = 0;
683 auto x = graph.GetX();
684 auto y = graph.GetY();
685 auto nPoints = graph.GetN();
686 for (auto p = 0; p < nPoints; ++p)
687 {
688 if (fabs(y[p] - fit.Eval(x[p])) < maxDifference)
689 continue;
690
691 remove[removePosition] = p;
692 if (++removePosition == remove.size())
693 {
694 removePoints(&remove[0], removePosition, graph);
695 totalRemoved += removePosition;
696 nPoints -= removePosition;
697 p -= removePosition;
698 removePosition = 0;
699
700 x = graph.GetX();
701 y = graph.GetY();
702 }
703 }
704
705 removePoints(&remove[0], removePosition, graph);
706 totalRemoved += removePosition;
707
708 graph.Fit(&fit, "QN");
709
710 return totalRemoved;
711 }
712
713 void HitCalibrationBar::removePoints(Int_t* points, Int_t nPoints, TGraph& graph) const
714 {
715 if (nPoints == 0)
716 return;
717
718 auto offset = 1;
719 const auto totalPoints = graph.GetN() - nPoints;
720 const auto x = graph.GetX();
721 const auto y = graph.GetY();
722
723 for (auto p = points[0]; p < totalPoints; ++p)
724 {
725
726 if (offset < nPoints && p + offset == points[offset])
727 {
728 ++offset;
729 --p;
730 continue;
731 }
732
733 x[p] = x[p + offset];
734 y[p] = y[p + offset];
735 }
736
737 graph.Set(totalPoints);
738 }
739
740 void HitCalibrationBar::createHistograms(TDirectory* histogramDir)
741 {
742 if (!histogramDir)
743 return;
744
745 // First creat all the histograms we do not have yet.
746 const auto nHits = (IsStatus(Validity, LogCompleteBit) ? CalibrationLogSize : LastHits.size());
747
748 TGraph posCal;
749 posCal.SetTitle("Position Calibration;Position / cm;Time Difference / ns");
750 posCal.Expand(nHits);
751
752 const auto loopSize = (LastHits.size() > nHits) ? nHits : LastHits.size();
753 for (auto h_index = 0; h_index < loopSize; ++h_index)
754 {
755 const auto& hit = LastHits[h_index];
756 const auto tdiff = hit.Time[1] - hit.Time[0];
757 posCal.SetPoint(posCal.GetN(), hit.EntryPosition, tdiff);
758 }
759
760 std::array<TF1, 2> fSatCal = { TF1("", "[0]*x/(1+[1]*x)", 0., TotalBarLength),
761 TF1("", "[0]*x/(1+[1]*x)", 0., TotalBarLength) };
762
763 for (auto side = 0; side < 2; ++side)
764 {
765 fSatCal[side].SetParameter(0, Gain[side]);
766 fSatCal[side].FixParameter(1, Saturation[side]);
767 }
768
769 // Go into the directory, create a Canvas, fill it with the plots and write it into the directory.
770 histogramDir->cd();
771 const auto canvname = TString::Format("Bar %d", ID);
772 TCanvas Canvas(canvname, canvname);
773 Canvas.SetWindowSize(1920, 1200);
774 Canvas.Divide(1, 3);
775 auto row = 0;
776 TVirtualPad* currentPad;
777
778 currentPad = Canvas.cd(++row);
779 currentPad->Divide(3, 1);
780
781 currentPad->cd(1);
782 if (posCal.GetN() > 0)
783 posCal.Draw("AP");
784
785 currentPad->cd(2);
786 if (Log.TimeDifference.GetN() > 0)
787 Log.TimeDifference.Draw("A*");
788
789 currentPad->cd(3);
790 if (Log.EffectiveSpeed.GetN() > 0)
791 Log.EffectiveSpeed.Draw("A*");
792
793 currentPad = Canvas.cd(++row);
794 currentPad->Divide(4, 1);
795
796 currentPad->cd(1);
797 Log.Energy[0].Draw("");
798 Log.Energy[1].SetLineColor(kRed);
799 Log.Energy[1].Draw("SAME");
800
801 currentPad->cd(2);
802 if (Log.LightAttenuationLength.GetN() > 0)
803 Log.LightAttenuationLength.Draw("A*");
804
805 currentPad->cd(3);
806 if (Log.Gain[0].GetN() > 0)
807 Log.Gain[0].Draw("A*");
808
809 currentPad->cd(4);
810 if (Log.Gain[1].GetN() > 0)
811 Log.Gain[1].Draw("A*");
812
813 currentPad = Canvas.cd(++row);
814 currentPad->Divide(4, 1);
815
816 TPad* mainPads[2];
817 TPad* overlayPads[2];
818 TGaxis* overlayAxis[2];
819 for (auto side = 0; side < 2; ++side)
820 {
821 // Threshold
822 auto pad = currentPad->cd(1 + side);
823 mainPads[side] = new TPad("", "", 0, 0, 1, 1);
824 mainPads[side]->Draw();
825 mainPads[side]->cd();
826 Log.TotalHits[side].Draw();
827 Log.MissedHits[side].Draw("SAME");
828 mainPads[side]->Update();
829 pad->cd();
830
831 overlayPads[side] = new TPad("", "", 0, 0, 1, 1);
832 overlayPads[side]->SetFillStyle(4000);
833 double ymin = 0.;
834 double ymax = 1.;
835 double dy = (ymax - ymin) / 0.8; // 10 per cent margins top and bottom
836 double xmin = mainPads[side]->GetUxmin();
837 double xmax = mainPads[side]->GetUxmax();
838 double dx = (xmax - xmin) / 0.8; // 10 per cent margins left and right
839 overlayPads[side]->Range(xmin - 0.1 * dx, ymin - 0.1 * dy, xmax + 0.1 * dx, ymax + 0.1 * dy);
840 overlayPads[side]->Draw();
841 overlayPads[side]->cd();
842
843 if (Log.MissRatio[side].GetN() > 0)
844 {
845 Log.MissRatio[side].SetLineColor(kRed);
846 Log.MissRatio[side].Draw("IL");
847 }
848 overlayPads[side]->Update();
849 overlayAxis[side] = new TGaxis(xmax, ymin, xmax, ymax, ymin, ymax, 50510, "+L");
850 overlayAxis[side]->SetLabelColor(kRed);
851 overlayAxis[side]->Draw();
852
853 overlayPads[side]->Update();
854
855 pad->Modified();
856
857 // Saturation
858
859 currentPad->cd(3 + side);
860 Log.Saturation[side].Draw("COLZ");
861 fSatCal[side].Draw("SAME");
862 }
863
864 Canvas.Write();
865
866 for (auto side = 0; side < 2; ++side)
867 {
868 delete mainPads[side];
869 delete overlayPads[side];
870 delete overlayAxis[side];
871 }
872 }
873
874 Double_t HitCalibrationBar::getMean(const TGraphErrors& graph, Double_t expectedValue)
875 {
876 const auto nPoints = graph.GetN();
877 if (nPoints == 0)
878 return NaN;
879
880 if (nPoints == 1)
881 return graph.GetY()[0];
882
883 TF1 constantFit("", "pol0", 0, nPoints);
884
885 FitGraph.Set(0);
886 for (auto point = 0; point < nPoints; ++point)
887 {
888 FitGraph.SetPoint(point, point, graph.GetY()[point]);
889 FitGraph.SetPointError(point, 0, graph.GetEY()[point]);
890 }
891 constantFit.SetParameter(0, expectedValue);
892 FitGraph.Fit(&constantFit, "QN");
893 return constantFit.GetParameter(0);
894 }
895 } // namespace Calibration
896} // namespace R3B::Neuland
#define R3BLOG(severity, x)
Definition R3BLogger.h:35
std::array< Double_t, 2 > DPair
Bool_t Add(const Double_t timeOffset, const Double_t entryPosition, const Double_t exitPosition, const Double_t energy, const UInt_t eventNumber)
void SetGlobalTSync(const Double_t value, const Double_t error)
void Set(const Int_t side, const Double_t time, const Int_t qdc)
Double_t GetPMTSaturation(const Int_t side) const
void SetModuleId(const Int_t i)
void SetPMTSaturation(const Double_t val, const Int_t side)
Int_t GetPedestal(const Int_t side) const
void SetTDiff(const Double_t val)
Double_t GetLightAttenuationLength() const
Double_t GetTimeOffset(const Int_t side) const
void SetTSync(const Double_t val)
void SetPedestal(const Int_t val, const Int_t side)
void SetEnergyGain(const Double_t val, const Int_t side)
void SetPMTThreshold(const Double_t val, const Int_t side)
void SetEffectiveSpeed(const Double_t val)
Double_t GetPMTThreshold(const Int_t side) const
void SetLightAttenuationLength(const Double_t val)
Double_t GetEnergyGain(const Int_t side) const
const double GetLightAtPMT(const double energyLoss, const double position, const double invLightAttenuationLength, const int side)
const bool IsStatus(Int_t var, Int_t statusBit)
bool GetJumps(double *data, unsigned int size, double_t threshold, unsigned int averageSize, std::vector< unsigned int > *output)
const void SetStatus(Int_t &var, Int_t statusBit)
const double FastExp(const double val)
const void ClearStatus(Int_t &var, Int_t statusBit)
Simulation of NeuLAND Bar/Paddle.
constexpr auto SaturationCoefficient
constexpr auto BarLength
constexpr auto TotalBarLength
constexpr auto NaN
constexpr auto AvgThreshold
constexpr T Sqr(const T val)
constexpr auto MaxCalTime