GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TacOffsetHelper.cxx
Go to the documentation of this file.
1#include "TacOffsetHelper.hh"
2
3void TacOffsetHelper::CreateHistograms(unsigned int slot)
4{
5 // find the offsets from the current calibration
6 int channelNumber = 0;
7 int currentIndex = 0;
8 while(currentIndex < 8) {
9 auto* channel = TChannel::GetChannelByNumber(++channelNumber);
10 if(channel == nullptr) {
11 break;
12 }
13 if(channel->GetClassType() != TTAC::Class()) {
14 continue;
15 }
16 fOffset[currentIndex] = static_cast<double>(channel->GetTimeOffset());
17 if(slot == 0) {
18 std::cout << "Current TAC offset in the calfile: " << fOffset[currentIndex] << " for channel #" << channelNumber << std::endl;
19 }
20 ++currentIndex;
21 }
22
23 //TAC offset histograms
24 for(int i = 0; i < fOffset.size(); ++i) {
25 fH1[slot][Form("TacOffset_%d", i)] = new TH1D(Form("TacOffset_%d", i), Form("Time difference between TAC and LaBr %d; time (ns); counts/ns", i), 10000, -5000., 5000.);
26 fH1[slot][Form("TacOffsetCorrected_%d", i)] = new TH1D(Form("TacOffsetCorrected_%d", i), Form("Time difference between TAC and LaBr %d, corrected by TAC offset; time (ns); counts/ns", i), 10000, -5000., 5000.);
27 fH1[slot][Form("TimeDiff_%d", i)] = new TH1D(Form("TimeDiff_%d", i), Form("Time difference for LaBr %d - LaBr with TAC coincidence; time (ns); counts/ns", i), 10000, -5000., 5000.);
28 fH1[slot][Form("TimeDiffNoTac_%d", i)] = new TH1D(Form("TimeDiffNoTac_%d", i), Form("Time difference for LaBr %d - LaBr without TAC coincidence; time (ns); counts/ns", i), 10000, -5000., 5000.);
29 fH1[slot][Form("TimeStampDiff_%d", i)] = new TH1D(Form("TimeStampDiff_%d", i), Form("Timestamp difference for LaBr %d - LaBr with TAC coincidence; time (ns); counts/ns", i), 10000, -5000., 5000.);
30 fH1[slot][Form("TimeStampDiffNoTac_%d", i)] = new TH1D(Form("TimeStampDiffNoTac_%d", i), Form("Timestamp difference for LaBr %d - LaBr without TAC coincidence; time (ns); counts/ns", i), 10000, -5000., 5000.);
31 }
32 fH1[slot]["TimeStampDiffGriffin"] = new TH1D("TimeStampDiffGriffin", "Timestamp difference for HPGe - LaBr, with TAC coincidence; time (ns); counts/ns", 10000, -5000., 5000.);
33 fH1[slot]["TimeDiffGriffin"] = new TH1D("TimeDiffGriffin", "Time difference for HPGe - LaBr, with TAC coincidence; time (ns); counts/ns", 10000, -5000., 5000.);
34}
35
36void TacOffsetHelper::Exec(unsigned int slot, TGriffin& grif, TTAC& tac, TLaBr& labr) // NOLINT
37{
38 // we use .at() here instead of [] so that we get meaningful error message if a histogram we try to fill wasn't created
39 // e.g. because of a typo
40
41 Short_t nofLaBr = labr.GetMultiplicity();
42 Short_t nofTac = tac.GetMultiplicity();
43 Short_t nofGrif = grif.GetMultiplicity();
44
45 // ignoring any events without exactly two LaBr
46 if(nofLaBr != 2) { return; }
47
48 auto* labr0 = labr.GetLaBrHit(0);
49 auto* labr1 = labr.GetLaBrHit(1);
50
51 // ignoring events where both LaBr hits are from the same detector
52 if(labr0->GetDetector() == labr1->GetDetector()) { return; }
53
54 // these histograms get filled with labr0 refering to the earlier hit and labr1 to the later hit,
55 // not to the lower and higher detector number like the other histograms!
56 fH1[slot].at(Form("TimeDiffNoTac_%d", labr1->GetDetector() - 1))->Fill(labr1->GetTime() - labr0->GetTime());
57 fH1[slot].at(Form("TimeStampDiffNoTac_%d", labr1->GetDetector() - 1))->Fill(static_cast<Double_t>(labr1->GetTimeStampNs() - labr0->GetTimeStampNs()));
58
59 // ensure that labr0 always refers to the LaBr with the lower detector number
60 if(labr0->GetDetector() > labr1->GetDetector()) {
61 std::swap(labr0, labr1);
62 }
63
64 if(nofTac == 1) {
65 auto* tac0 = tac.GetTACHit(0);
66 if(labr0->GetDetector() == tac0->GetDetector()) {
67 fH1[slot].at(Form("TacOffset_%d", tac0->GetDetector() - 1))->Fill(labr0->GetTime() - tac0->GetTime() + fOffset[tac0->GetDetector() - 1] * 10.);
68 fH1[slot].at(Form("TacOffsetCorrected_%d", tac0->GetDetector() - 1))->Fill(labr0->GetTime() - tac0->GetTime());
69 fH1[slot].at(Form("TimeDiff_%d", labr1->GetDetector() - 1))->Fill(labr1->GetTime() - labr1->GetTime());
70 fH1[slot].at(Form("TimeStampDiff_%d", labr1->GetDetector() - 1))->Fill(static_cast<Double_t>(labr1->GetTimeStampNs() - labr0->GetTimeStampNs()));
71 for(int i = 0; i < nofGrif; ++i) {
72 auto* grif1 = grif.GetGriffinHit(i);
73 fH1[slot].at("TimeDiffGriffin")->Fill(grif1->GetTime() - labr0->GetTime());
74 fH1[slot].at("TimeStampDiffGriffin")->Fill(static_cast<Double_t>(grif1->GetTimeStampNs() - labr0->GetTimeStampNs()));
75 }
76 }
77 }
78}
79
80void TacOffsetHelper::EndOfSort(std::shared_ptr<std::map<std::string, TList>>& list)
81{
82}
static TChannel * GetChannelByNumber(int temp_num)
Definition TChannel.cxx:491
virtual Short_t GetMultiplicity() const
Definition TDetector.h:70
std::vector< TGRSIMap< std::string, TH1 * > > fH1
Definition TGRSIHelper.h:44
TGriffinHit * GetGriffinHit(const int &i)
!
Definition TGriffin.cxx:243
Definition TLaBr.h:28
TLaBrHit * GetLaBrHit(const int &i) const
Definition TLaBr.h:49
Definition TTAC.h:27
TTACHit * GetTACHit(const int &i) const
Definition TTAC.h:36
void Exec(unsigned int slot, TGriffin &grif, TTAC &tac, TLaBr &labr)
std::array< double, 8 > fOffset
void EndOfSort(std::shared_ptr< std::map< std::string, TList > > &list) override
This method gets called at the end of Finalize()
void CreateHistograms(unsigned int slot) override
Virtual helper function that the user uses to create their histograms.