GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
CrossTalkHelper.cxx
Go to the documentation of this file.
1#include "CrossTalkHelper.hh"
2
3bool pileup_reject = false;
4
5const double gg_time_low = -200.;
6const double gg_time_high = 300.;
7
8// Beta gamma histograms
9const double gb_time_low = -250.;
10const double gb_time_high = 350.;
11
12// default k-value this might have to be adjusted for each experiment!!!
13const Short_t defaultKValue = 379;
14
16{
17 return ((one->GetDetector() == two->GetDetector()) && (std::fabs(one->GetTime() - two->GetTime()) < 300.));
18}
19
21{
22
23 return ((two->GetTime() - one->GetTime()) >= gg_time_low) && ((two->GetTime() - one->GetTime()) <= gg_time_high);
24}
25
26void CrossTalkHelper::CreateHistograms(unsigned int slot)
27{
28 for(int det_num = 1; det_num <= 16; ++det_num) {
29 std::string aedet_str = Form("aEdet%d", det_num);
30 std::string gedet_str = Form("gEdet%d", det_num);
31 std::string ae2det_str = Form("aE2det%d", det_num);
32 fH1[slot][aedet_str] = new TH1D(Form("aEdet%d", det_num), Form("Addback detector %d", det_num), 1500, 0, 1500);
33 fH1[slot][gedet_str] = new TH1D(Form("geEdet%d", det_num), Form("Singles detector %d", det_num), 1500, 0, 1500);
34 fH1[slot][ae2det_str] = new TH1D(Form("aE2det%d", det_num), Form("Addback with 2 hits, detector %d", det_num), 1500, 0, 1500);
35 for(int crys_1 = 0; crys_1 < 4; ++crys_1) {
36 for(int crys_2 = crys_1 + 1; crys_2 < 4; ++crys_2) {
37 std::string name_str = Form("det_%d_%d_%d", det_num, crys_1, crys_2);
38 const char* hist_name = name_str.c_str();
39 std::cout << "Creating histogram: " << hist_name;
40 fH2[slot][name_str] = new TH2I(hist_name, hist_name, 1500, 0, 1500, 1500, 0, 1500);
41 std::cout << " at address: " << fH2[slot][hist_name] << std::endl;
42 }
43 }
44 }
45
46 fH1[slot]["aMult"] = new TH1D("aMult", "addback multilpicity", 20, 0, 20);
47 fH2[slot]["gE_chan"] = new TH2D("gE_chan", "gE_chan", 65, 0, 65, 1500, 0, 1500);
48 fH1[slot]["aE"] = new TH1D("aE", "Summed Addback", 1500, 0, 1500);
49 fH1[slot]["gE"] = new TH1D("gE", "Summed Singles", 1500, 0, 1500);
50 fH1[slot]["gEnoCT"] = new TH1D("gEnoCT", "Singles, no CT correction", 1500, 0, 1500);
51}
52
53void CrossTalkHelper::Exec(unsigned int slot, TGriffin& grif, TGriffinBgo& grifBgo)
54{
55 // find the multiplicity in each clover over the entire event
56 // we do this because we want to force a multiplicity of 2
57 std::array<Int_t, 17> detMultiplicity = {0};
58 for(auto gr1 = 0; gr1 < grif.GetSuppressedMultiplicity(&grifBgo); ++gr1) {
59 ++(detMultiplicity[grif.GetSuppressedHit(gr1)->GetDetector()]);
60 }
61
62 for(auto gr1 = 0; gr1 < grif.GetSuppressedMultiplicity(&grifBgo); ++gr1) {
63 auto grif1 = grif.GetSuppressedHit(gr1);
64 if(pileup_reject && (grif1->GetKValue() != defaultKValue)) continue; // This pileup number might have to change for other expmnts
65 fH1[slot].at(Form("gEdet%d", grif1->GetDetector()))->Fill(grif1->GetEnergy());
66 fH2[slot].at("gE_chan")->Fill(grif1->GetArrayNumber(), grif1->GetEnergy());
67 fH1[slot].at("gE")->Fill(grif1->GetEnergy());
68 fH1[slot].at("gEnoCT")->Fill(grif1->GetNoCTEnergy());
69 for(auto gr2 = gr1 + 1; gr2 < grif.GetSuppressedMultiplicity(&grifBgo); ++gr2) {
70 auto grif2 = grif.GetSuppressedHit(gr2);
71 if(pileup_reject && grif2->GetKValue() != defaultKValue) continue; // This pileup number might have to change for other expmnts
72 if((detMultiplicity[grif1->GetDetector()] == 2) && Addback(grif1, grif2)) {
73 TGriffinHit *low_crys_hit, *high_crys_hit;
74 if(grif1->GetCrystal() < grif2->GetCrystal()) {
75 low_crys_hit = grif1;
76 high_crys_hit = grif2;
77 } else {
78 low_crys_hit = grif2;
79 high_crys_hit = grif1;
80 }
81 if(low_crys_hit->GetCrystal() != high_crys_hit->GetCrystal()) {
82 fH2[slot].at(Form("det_%d_%d_%d", low_crys_hit->GetDetector(), low_crys_hit->GetCrystal(), high_crys_hit->GetCrystal()))->Fill(low_crys_hit->GetNoCTEnergy(), high_crys_hit->GetNoCTEnergy());
83 }
84 }
85 }
86 }
87
88 for(auto gr1 = 0; gr1 < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++gr1) {
89 auto grif1 = grif.GetSuppressedAddbackHit(gr1);
90 if(pileup_reject && (grif1->GetKValue() != defaultKValue)) continue; // This pileup number might have to change for other expmnts
91 fH1[slot].at("aE")->Fill(grif1->GetEnergy());
92 fH1[slot].at(Form("aEdet%d", grif1->GetDetector()))->Fill(grif1->GetEnergy());
93 fH1[slot].at("aMult")->Fill(grif.GetNSuppressedAddbackFrags(gr1));
94 if(grif.GetNSuppressedAddbackFrags(gr1) == 2) fH1[slot].at(Form("aE2det%d", grif1->GetDetector()))->Fill(grif1->GetEnergy());
95 }
96}
bool pileup_reject
const Short_t defaultKValue
const double gb_time_low
const double gg_time_high
bool Addback(TGriffinHit *one, TGriffinHit *two)
bool PromptCoincidence(TGriffinHit *one, TGriffinHit *two)
const double gg_time_low
const double gb_time_high
void Exec(unsigned int slot, TGriffin &grif, TGriffinBgo &grifBgo)
void CreateHistograms(unsigned int slot) override
Virtual helper function that the user uses to create their histograms.
virtual Int_t GetCrystal() const
!
virtual Int_t GetDetector() const
!
virtual Double_t GetTime(const ETimeFlag &correct_flag=ETimeFlag::kAll, Option_t *opt="") const
Returns a time value to the nearest nanosecond!
std::vector< TGRSIMap< std::string, TH2 * > > fH2
Definition TGRSIHelper.h:49
std::vector< TGRSIMap< std::string, TH1 * > > fH1
Definition TGRSIHelper.h:48
Double_t GetNoCTEnergy(Option_t *opt="") const
UShort_t GetNSuppressedAddbackFrags(const size_t &idx)
Definition TGriffin.h:129
TGriffinHit * GetSuppressedAddbackHit(const int &i)
Definition TGriffin.h:122
Short_t GetSuppressedMultiplicity(const TBgo *bgo)
Definition TGriffin.h:111
Short_t GetSuppressedAddbackMultiplicity(const TBgo *bgo)
Definition TGriffin.h:119
TGriffinHit * GetSuppressedHit(const int &i)
!
Definition TGriffin.h:108