19#include "FairLogger.h"
22#include "TDirectory.h"
31using DPair = std::array<Double_t, 2>;
41 unsigned int averageSize,
42 std::vector<unsigned int>* output)
44 if (size < 2 * averageSize)
47 const auto scale = 1. / averageSize;
49 auto lastAverage = std::accumulate(data, data + averageSize, 0.) * scale;
50 const auto last = data + size;
51 auto currentPos = data + averageSize;
53 while (currentPos + averageSize < last)
55 const auto currentAverage = std::accumulate(currentPos, currentPos + averageSize, 0.) * scale;
56 if (fabs(currentAverage - lastAverage) > threshold)
58 if (output ==
nullptr)
60 output->push_back(currentPos - data);
62 lastAverage = currentAverage;
63 currentPos += averageSize;
69 template <Int_t iterations = 8>
72 auto exp = 1. + val / (1 << iterations);
73 for (
auto i = 0; i < iterations; ++i)
81 const double position,
82 const double invLightAttenuationLength,
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)); }
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",
128 const char* HitCalibrationBar::CalibrationStatusAbbreviation[] = {
"XX",
"FA",
" ",
"EC",
129 "TS",
"SE",
"??",
"TJ" };
136 , TimeDifference(
NaN)
137 , EffectiveSpeed(
NaN)
139 , InvLightAttenuationLength(
NaN)
141 , Saturation({ 0., 0. })
150 Log.TimeDifference.SetTitle(
"Time Difference Parameter;Event Number;Parameter Value / ns");
153 Log.EffectiveSpeed.SetTitle(
"Effective Speed Parameter;Event Number;Parameter Value / cm/ns");
156 Log.LightAttenuationLength.SetTitle(
"Attenuation Length Parameter;Event Number;Parameter Value / cm");
158 for (
auto side = 0; side < 2; ++side)
161 Log.Gain[side].SetTitle(
"Gain Parameter;Event Number;Parameter Value / ns/MeVle");
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);
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);
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);
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);
181 Log.Pedestal[side] = TH1F(
"",
"Pedestal; QDC / ns", 512, 0., 512.);
182 Log.Pedestal[side].SetDirectory(
nullptr);
183 Log.Pedestal[side].SetStats(
false);
230 LOG(debug) <<
"Loaded Parameters for Bar " << ID;
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);
244 if (fabs(CurrentHit.Time[1] - CurrentHit.Time[0]) > 0.5 *
MaxCalTime)
246 if (CurrentHit.Time[1] < CurrentHit.Time[0])
254 CurrentHit.EntryPosition =
255 (CurrentHit.Time[1] - CurrentHit.Time[0] - TimeDifference) * EffectiveSpeed;
260 if (FailCounter > 10)
262 LOG(debug) <<
"Bar " << ID <<
":"
263 <<
" Positions were too often too strange. Might be a Time-Jump in one PMT. "
264 "Removing Position-Calibration.";
278 TimeSync = { value, error };
279 if (std::isnan(value))
288 return CalibrationStatus::noData;
290 return CalibrationStatus::completeFail;
292 return CalibrationStatus::energyGainFail;
294 return CalibrationStatus::timeSyncFail;
296 return CalibrationStatus::timeJump;
298 return CalibrationStatus::timeSyncError;
300 return CalibrationStatus::strangeParameters;
301 return CalibrationStatus::success;
305 const Double_t entryPosition,
306 const Double_t exitPosition,
307 const Double_t energy,
308 const UInt_t eventNumber)
315 const auto centerPosition = (entryPosition + exitPosition) * 0.5;
328 for (
auto side = 0; side < 2; ++side)
330 Log.TotalHits[side].Fill(lightAtPMT[side]);
335 for (
auto side = 0; side < 2; ++side)
337 Log.Saturation[side].Fill(lightAtPMT[side], CurrentHit.QDC[side]);
338 const auto gain = CurrentHit.QDC[side] /
340 Log.Energy[side].Fill(gain);
352 const auto maxPos =
BarLength * 0.5 - 1.;
353 if (fabs(entryPosition) > maxPos || fabs(exitPosition) > maxPos)
356 LastEventNumber = eventNumber;
358 for (
auto side = 0; side < 2; ++side)
359 CurrentHit.Time[side] -= timeOffset * 0.5;
361 CurrentHit.EntryPosition = entryPosition;
362 CurrentHit.ExitPosition = exitPosition;
363 CurrentHit.Energy = energy;
365 LastHits.push_back(CurrentHit);
385 fabs(0.5 * (entryPosition + exitPosition) -
GetPosition()) > 25.)
395 CurrentHit.EntryPosition =
NaN;
409 if (Log.TimeDifference.GetN() > 0)
411 TimeDifference = getMean(Log.TimeDifference);
413 EffectiveSpeed = getMean(Log.EffectiveSpeed);
425 positionCalibration(0, nHits);
435 if (Log.Gain[0].GetN() > 0 && Log.Gain[1].GetN() > 0)
437 for (
auto side = 0; side < 2; ++side)
439 Gain[side] = getMean(Log.Gain[side], Gain[side]);
443 InvLightAttenuationLength = 1. / getMean(Log.LightAttenuationLength, 1. / InvLightAttenuationLength);
451 energyCalibration(0, nHits);
461 thresholdCalibration();
463 createHistograms(histogramDir);
476 for (
auto side = 0; side < 2; ++side)
484 LOG(debug) << ID <<
" Bar Parameters: ";
503 void HitCalibrationBar::positionCalibration(
int firstHit,
int nHits)
506 const auto loopSize = (LastHits.size() >= (firstHit + nHits)) ? nHits : (LastHits.size() - firstHit);
507 for (
auto index = 0; index < loopSize; ++index)
509 const auto& hit = LastHits[firstHit + index];
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);
517 FitGraph.Fit(&linearFit,
"NQ");
519 cleanupFit(FitGraph, linearFit, 2.5);
522 TimeDifference = linearFit.GetParameter(0);
523 EffectiveSpeed = 1. / linearFit.GetParameter(1);
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));
535 void HitCalibrationBar::energyCalibration(
int firstHit,
int nHits)
537 const auto loopSize = (LastHits.size() >= (firstHit + nHits)) ? nHits : (LastHits.size() - firstHit);
543 pedestalCalibration();
544 for (
auto h_index = 0; h_index < loopSize; ++h_index)
546 auto& hit = LastHits[firstHit + h_index];
547 for (
auto side = 0; side < 2; ++side)
549 hit.QDC[side] = std::max(hit.QDC[side] - Pedestal[side], 1);
555 for (
auto h_index = 0; h_index < loopSize; ++h_index)
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]));
564 cleanupFit(FitGraph, linearFit, 1.);
566 InvLightAttenuationLength = 0.5 * linearFit.GetParameter(1);
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));
573 for (
auto h_index = 0; h_index < loopSize; ++h_index)
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)
579 const auto lightAtPMT =
GetLightAtPMT(hit.Energy, centerPosition, InvLightAttenuationLength, side);
581 Log.Energy[side].Fill(gain);
585 std::array<TH1F, 2> hGainCal = { TH1F(
"",
"", 200, 0., 40.), TH1F(
"",
"", 200, 0., 40.) };
587 for (
auto h_index = 0; h_index < loopSize; ++h_index)
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)
593 const auto lightAtPMT =
GetLightAtPMT(hit.Energy, centerPosition, InvLightAttenuationLength, side);
596 hGainCal[side].Fill(gain);
600 TF1 gausFit(
"",
"gaus", 0., 30.);
601 for (
auto side = 0; side < 2; ++side)
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);
608 Gain[side] = gausFit.GetParameter(1);
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));
618 void HitCalibrationBar::pedestalCalibration()
620 for (
auto side = 0; side < 2; ++side)
621 for (
auto bin = 0; bin < Log.Pedestal[side].GetNbinsX(); ++bin)
623 if (Log.Pedestal[side].GetBinContent(bin + 1) > 3.)
625 Pedestal[side] = bin;
633 void HitCalibrationBar::thresholdCalibration()
642 TF1 missFit(
"missFit",
"(1. - [2]) * 0.5 * (1. - TMath::Erf([0] * (x - [1]))) + [2]", 0., 20.);
644 for (
auto side = 0; side < 2; ++side)
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)
651 if (Log.MissedHits[side].GetBinContent(bin) > 0)
653 const auto nextPoint = Log.MissRatio[side].GetN();
654 const auto binCenter = Log.MissedHits[side].GetBinCenter(bin);
656 Log.MissedHits[side].GetBinContent(bin) / Log.TotalHits[side].GetBinContent(bin);
657 Log.MissRatio[side].SetPoint(nextPoint, binCenter, ratio);
661 if (Log.MissRatio[side].GetN() > 5)
663 missFit.SetParameter(0, 0.25);
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);
675 Int_t HitCalibrationBar::cleanupFit(TGraph& graph, TF1& fit, Double_t maxDifference)
const
677 std::array<Int_t, 256> remove;
678 auto totalRemoved = 0;
680 graph.Fit(&fit,
"QN ROB=0.90");
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)
688 if (fabs(y[p] - fit.Eval(x[p])) < maxDifference)
691 remove[removePosition] = p;
692 if (++removePosition == remove.size())
694 removePoints(&remove[0], removePosition, graph);
695 totalRemoved += removePosition;
696 nPoints -= removePosition;
705 removePoints(&remove[0], removePosition, graph);
706 totalRemoved += removePosition;
708 graph.Fit(&fit,
"QN");
713 void HitCalibrationBar::removePoints(Int_t* points, Int_t nPoints, TGraph& graph)
const
719 const auto totalPoints = graph.GetN() - nPoints;
720 const auto x = graph.GetX();
721 const auto y = graph.GetY();
723 for (
auto p = points[0]; p < totalPoints; ++p)
726 if (offset < nPoints && p + offset == points[offset])
733 x[p] = x[p + offset];
734 y[p] = y[p + offset];
737 graph.Set(totalPoints);
740 void HitCalibrationBar::createHistograms(TDirectory* histogramDir)
749 posCal.SetTitle(
"Position Calibration;Position / cm;Time Difference / ns");
750 posCal.Expand(nHits);
752 const auto loopSize = (LastHits.size() > nHits) ? nHits : LastHits.size();
753 for (
auto h_index = 0; h_index < loopSize; ++h_index)
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);
760 std::array<TF1, 2> fSatCal = { TF1(
"",
"[0]*x/(1+[1]*x)", 0.,
TotalBarLength),
763 for (
auto side = 0; side < 2; ++side)
765 fSatCal[side].SetParameter(0, Gain[side]);
766 fSatCal[side].FixParameter(1, Saturation[side]);
771 const auto canvname = TString::Format(
"Bar %d", ID);
772 TCanvas Canvas(canvname, canvname);
773 Canvas.SetWindowSize(1920, 1200);
776 TVirtualPad* currentPad;
778 currentPad = Canvas.cd(++row);
779 currentPad->Divide(3, 1);
782 if (posCal.GetN() > 0)
786 if (Log.TimeDifference.GetN() > 0)
787 Log.TimeDifference.Draw(
"A*");
790 if (Log.EffectiveSpeed.GetN() > 0)
791 Log.EffectiveSpeed.Draw(
"A*");
793 currentPad = Canvas.cd(++row);
794 currentPad->Divide(4, 1);
797 Log.Energy[0].Draw(
"");
798 Log.Energy[1].SetLineColor(kRed);
799 Log.Energy[1].Draw(
"SAME");
802 if (Log.LightAttenuationLength.GetN() > 0)
803 Log.LightAttenuationLength.Draw(
"A*");
806 if (Log.Gain[0].GetN() > 0)
807 Log.Gain[0].Draw(
"A*");
810 if (Log.Gain[1].GetN() > 0)
811 Log.Gain[1].Draw(
"A*");
813 currentPad = Canvas.cd(++row);
814 currentPad->Divide(4, 1);
817 TPad* overlayPads[2];
818 TGaxis* overlayAxis[2];
819 for (
auto side = 0; side < 2; ++side)
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();
831 overlayPads[side] =
new TPad(
"",
"", 0, 0, 1, 1);
832 overlayPads[side]->SetFillStyle(4000);
835 double dy = (ymax - ymin) / 0.8;
836 double xmin = mainPads[side]->GetUxmin();
837 double xmax = mainPads[side]->GetUxmax();
838 double dx = (xmax - xmin) / 0.8;
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();
843 if (Log.MissRatio[side].GetN() > 0)
845 Log.MissRatio[side].SetLineColor(kRed);
846 Log.MissRatio[side].Draw(
"IL");
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();
853 overlayPads[side]->Update();
859 currentPad->cd(3 + side);
860 Log.Saturation[side].Draw(
"COLZ");
861 fSatCal[side].Draw(
"SAME");
866 for (
auto side = 0; side < 2; ++side)
868 delete mainPads[side];
869 delete overlayPads[side];
870 delete overlayAxis[side];
874 Double_t HitCalibrationBar::getMean(
const TGraphErrors& graph, Double_t expectedValue)
876 const auto nPoints = graph.GetN();
881 return graph.GetY()[0];
883 TF1 constantFit(
"",
"pol0", 0, nPoints);
886 for (
auto point = 0; point < nPoints; ++point)
888 FitGraph.SetPoint(point, point, graph.GetY()[point]);
889 FitGraph.SetPointError(point, 0, graph.GetEY()[point]);
891 constantFit.SetParameter(0, expectedValue);
892 FitGraph.Fit(&constantFit,
"QN");
893 return constantFit.GetParameter(0);
#define R3BLOG(severity, x)
std::array< Double_t, 2 > DPair
HitCalibrationBar(const Int_t id=0)
R3BNeulandHitModulePar GetParameters()
Bool_t Add(const Double_t timeOffset, const Double_t entryPosition, const Double_t exitPosition, const Double_t energy, const UInt_t eventNumber)
void Update(const R3BNeulandHitModulePar *par)
void SetGlobalTSync(const Double_t value, const Double_t error)
CalibrationStatus GetCalibrationStatus() const
void Calibrate(TDirectory *histogramDir=nullptr)
void Set(const Int_t side, const Double_t time, const Int_t qdc)
Double_t GetPosition() const
Double_t GetPMTSaturation(const Int_t side) const
void SetModuleId(const Int_t i)
void SetPMTSaturation(const Double_t val, const Int_t side)
Double_t GetEffectiveSpeed() const
Int_t GetPedestal(const Int_t side) const
void SetTDiff(const Double_t val)
Double_t GetLightAttenuationLength() const
Double_t GetTDiff() 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
constexpr auto LogCompleteBit
constexpr auto CalibrationLogSize
constexpr auto PedestalCalibrationBit
constexpr auto MaxTSyncError
constexpr auto RightSideHitBit
const double GetLightAtPMT(const double energyLoss, const double position, const double invLightAttenuationLength, const int side)
constexpr auto LogInitialSize
constexpr auto MaxFastTDiffError
constexpr auto ThresholdCalibrationBit
constexpr auto MaxNumberOfFails
constexpr auto TimeJumpThreshold
constexpr auto EnergyCalibrationBit
const bool IsStatus(Int_t var, Int_t statusBit)
constexpr auto TimeJumpedBit
bool GetJumps(double *data, unsigned int size, double_t threshold, unsigned int averageSize, std::vector< unsigned int > *output)
constexpr auto TSyncCalibrationBit
constexpr auto PosCalibrationBit
const void SetStatus(Int_t &var, Int_t statusBit)
constexpr auto StrangeBit
constexpr auto LeftSideHitBit
constexpr auto MinEnergyDeposit
const double FastExp(const double val)
const void ClearStatus(Int_t &var, Int_t statusBit)
constexpr auto ThresholdCalibrationSize
constexpr auto EnergyCalibrationSize
constexpr auto PositionCalibrationSize
Simulation of NeuLAND Bar/Paddle.
constexpr auto SaturationCoefficient
constexpr auto TotalBarLength
constexpr auto AvgThreshold
constexpr T Sqr(const T val)
constexpr auto MaxCalTime