R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandHitCalibrationBar.h
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
14#ifndef R3BNEULANDHITCALIBRATIONBAR_H
15#define R3BNEULANDHITCALIBRATIONBAR_H
16
17#include <array>
18#include <vector>
19
20#include "TF1.h"
21#include "TGraphErrors.h"
22#include "TH1F.h"
23#include "TH2F.h"
24
26
27class TGraph;
29class TDirectory;
30
31namespace R3B::Neuland // NOLINT
32{
33 namespace Calibration
34 {
35 class TSyncer;
36
38 {
39 using DPair = std::array<Double_t, 2>;
40 using IPair = std::array<Int_t, 2>;
41
42 public:
54
55 HitCalibrationBar(const Int_t id = 0);
57
58 void Update(const R3BNeulandHitModulePar* par);
59
60 void Set(const Int_t side, const Double_t time, const Int_t qdc);
61 Bool_t Add(const Double_t timeOffset,
62 const Double_t entryPosition,
63 const Double_t exitPosition,
64 const Double_t energy,
65 const UInt_t eventNumber);
66 void Reset();
67 void Calibrate(TDirectory* histogramDir = nullptr);
70 Bool_t IsValid() const;
71 inline Double_t GetTime() const { return 0.5 * (CurrentHit.Time[0] + CurrentHit.Time[1]); }
72 inline Double_t GetPosition() const { return CurrentHit.EntryPosition; }
73
74 void SetGlobalTSync(const Double_t value, const Double_t error);
75
76 inline static const char* const GetCalibrationStatusDescription(const CalibrationStatus status)
77 {
78 return CalibrationStatusDescription[static_cast<int>(status)];
79 }
80 inline static const char* const GetCalibrationStatusAbbreviation(const CalibrationStatus status)
81 {
82 return CalibrationStatusAbbreviation[static_cast<int>(status)];
83 }
84
85 private:
86 static const char* CalibrationStatusDescription[];
87 static const char* CalibrationStatusAbbreviation[];
88
89 Int_t ID;
90 Int_t Validity;
91 UInt_t FailCounter;
92 UInt_t LastEventNumber;
93
94 struct Hit
95 {
96 DPair Time;
97 IPair QDC;
98 Double_t EntryPosition;
99 Double_t ExitPosition;
100 Double_t Energy;
101 } CurrentHit;
102 std::vector<Hit> LastHits;
103
104 Double_t TimeDifference;
105 Double_t EffectiveSpeed;
106 DPair Gain;
107 Double_t InvLightAttenuationLength;
108 IPair Pedestal;
109 DPair Saturation;
110 DPair Threshold;
111 DPair TimeSync;
112
113 struct
114 {
115 TGraphErrors TimeDifference;
116 TGraphErrors EffectiveSpeed;
117 std::array<TGraphErrors, 2> Gain;
119 std::array<TH1F, 2> TotalHits;
120 std::array<TH1F, 2> MissedHits;
121 std::array<TGraph, 2> MissRatio;
122 std::array<TH1F, 2> Energy;
123 std::array<TH2F, 2> Saturation;
124 std::array<TH1F, 2> Pedestal;
125 } Log;
126
127 // Some members so we do not have to create and delete them all the time
128 TGraphErrors FitGraph;
129
130 void positionCalibration(int firstHit, int nHits);
131 void energyCalibration(int firstHit, int nHits);
132 void pedestalCalibration();
133 void thresholdCalibration();
134
135 Int_t cleanupFit(TGraph& graph, TF1& fit, Double_t maxDifference) const;
136
144 void removePoints(Int_t* points, Int_t nPoints, TGraph& graph) const;
145
151 void createHistograms(TDirectory* histogramDir);
152
153 Double_t getMean(const TGraphErrors& graph, Double_t expectedValue = 0.);
154 };
155 } // namespace Calibration
156} // namespace R3B::Neuland
157#endif
std::array< Double_t, 2 > DPair
R3B::Neuland::Calibration::HitCalibrationBar::CalibrationStatus CalibrationStatus
static const char *const GetCalibrationStatusDescription(const CalibrationStatus status)
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)
static const char *const GetCalibrationStatusAbbreviation(const CalibrationStatus status)
std::array< Double_t, 2 > DPair
Simulation of NeuLAND Bar/Paddle.