GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
CrossTalk.C
Go to the documentation of this file.
1#define CrossTalk_cxx
2// The class definition in CrossTalk.h has been generated automatically
3#include "CrossTalk.h"
4
5#include "TH1.h"
6
7bool pileup_reject = true;
8
9const double gg_time_low = -200.;
10const double gg_time_high = 300.;
11
12// Beta gamma histograms
13const double gb_time_low = -250.;
14const double gb_time_high = 350.;
15
16// default k-value this might have to be adjusted for each experiment!!!
17const Short_t defaultKValue = 379;
18
20{
21 return ((one.GetDetector() == two.GetDetector()) && (std::fabs(one.GetTime() - two.GetTime()) < 300.));
22}
23
25{
26
27 return ((two->GetTime() - one->GetTime()) >= gg_time_low) && ((two->GetTime() - one->GetTime()) <= gg_time_high);
28}
29
31{
32 fH2.clear();
33 for(int det_num = 1; det_num <= 16; ++det_num) {
34 std::string aedet_str = Form("aEdet%d", det_num);
35 std::string gedet_str = Form("gEdet%d", det_num);
36 std::string ae2det_str = Form("aE2det%d", det_num);
37 fH1[aedet_str] = new TH1D(Form("aEdet%d", det_num), Form("Addback detector %d", det_num), 1500, 0, 1500);
38 fH1[gedet_str] = new TH1D(Form("geEdet%d", det_num), Form("Singles detector %d", det_num), 1500, 0, 1500);
39 fH1[ae2det_str] = new TH1D(Form("aE2det%d", det_num), Form("Addback with 2 hits, detector %d", det_num), 1500, 0, 1500);
40 for(int crys_1 = 0; crys_1 < 4; ++crys_1) {
41 for(int crys_2 = crys_1 + 1; crys_2 < 4; ++crys_2) {
42 std::string name_str = Form("det_%d_%d_%d", det_num, crys_1, crys_2);
43 const char* hist_name = name_str.c_str();
44 std::cout << "Creating histogram: " << hist_name;
45 fH2[name_str] = new TH2I(hist_name, hist_name, 1500, 0, 1500, 1500, 0, 1500);
46 std::cout << " at address: " << fH2[hist_name] << std::endl;
47 }
48 }
49 }
50
51 fH1["aMult"] = new TH1D("aMult", "addback multilpicity", 20, 0, 20);
52 fH2["gE_chan"] = new TH2D("gE_chan", "gE_chan", 65, 0, 65, 1500, 0, 1500);
53 fH1["aE"] = new TH1D("aE", "Summed Addback", 1500, 0, 1500);
54 fH1["gE"] = new TH1D("gE", "Summed Singles", 1500, 0, 1500);
55 fH1["gEnoCT"] = new TH1D("gEnoCT", "Singles, no CT correction", 1500, 0, 1500);
56
57 for(auto it : fH1) {
58 GetOutputList()->Add(it.second);
59 }
60 for(auto it : fH2) {
61 GetOutputList()->Add(it.second);
62 }
63 for(auto it : fHSparse) {
64 GetOutputList()->Add(it.second);
65 }
66}
67
69{
70 // find the multiplicity in each clover over the entire event
71 // we do this because we want to force a multiplicity of 2
72 Int_t det_multiplicity[17] = {0};
73 for(auto gr1 = 0; gr1 < fGrif->GetSuppressedMultiplicity(fGriffinBgo); ++gr1) {
74 ++(det_multiplicity[fGrif->GetSuppressedHit(gr1)->GetDetector()]);
75 }
76
77 for(auto gr1 = 0; gr1 < fGrif->GetSuppressedMultiplicity(fGriffinBgo); ++gr1) {
78 if(pileup_reject && (fGrif->GetSuppressedHit(gr1)->GetKValue() != defaultKValue)) continue; // This pileup number might have to change for other expmnts
79 fH1[Form("gEdet%d", fGrif->GetSuppressedHit(gr1)->GetDetector())]->Fill(fGrif->GetSuppressedHit(gr1)->GetEnergy());
80 fH2["gE_chan"]->Fill(fGrif->GetSuppressedHit(gr1)->GetArrayNumber(), fGrif->GetSuppressedHit(gr1)->GetEnergy());
81 fH1["gE"]->Fill(fGrif->GetSuppressedHit(gr1)->GetEnergy());
82 fH1["gEnoCT"]->Fill(fGrif->GetSuppressedHit(gr1)->GetNoCTEnergy());
83 for(auto gr2 = gr1 + 1; gr2 < fGrif->GetSuppressedMultiplicity(fGriffinBgo); ++gr2) {
84 if(pileup_reject && fGrif->GetSuppressedHit(gr2)->GetKValue() != defaultKValue) continue; // This pileup number might have to change for other expmnts
85 if((det_multiplicity[fGrif->GetSuppressedHit(gr1)->GetDetector()] == 2) && Addback(*(fGrif->GetSuppressedHit(gr1)), *(fGrif->GetSuppressedHit(gr2)))) {
86 TGriffinHit *low_crys_hit, *high_crys_hit;
88 low_crys_hit = fGrif->GetSuppressedHit(gr1);
89 high_crys_hit = fGrif->GetSuppressedHit(gr2);
90 } else {
91 low_crys_hit = fGrif->GetSuppressedHit(gr2);
92 high_crys_hit = fGrif->GetSuppressedHit(gr1);
93 }
94 if(low_crys_hit->GetCrystal() != high_crys_hit->GetCrystal()) {
95 fH2[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());
96 }
97 }
98 }
99 }
100
101 for(auto gr1 = 0; gr1 < fGrif->GetSuppressedAddbackMultiplicity(fGriffinBgo); ++gr1) {
103 continue; // This pileup number might have to change for other expmnts
104 fH1["aE"]->Fill(fGrif->GetSuppressedAddbackHit(gr1)->GetEnergy());
105 fH1[Form("aEdet%d", fGrif->GetSuppressedAddbackHit(gr1)->GetDetector())]->Fill(fGrif->GetSuppressedAddbackHit(gr1)->GetEnergy());
106 fH1["aMult"]->Fill(fGrif->GetNSuppressedAddbackFrags(gr1));
107 if(fGrif->GetNSuppressedAddbackFrags(gr1) == 2)
108 fH1[Form("aE2det%d", fGrif->GetSuppressedAddbackHit(gr1)->GetDetector())]->Fill(fGrif->GetSuppressedAddbackHit(gr1)->GetEnergy());
109 }
110}
bool pileup_reject
Definition CrossTalk.C:7
bool Addback(TGriffinHit &one, TGriffinHit &two)
Definition CrossTalk.C:19
const Short_t defaultKValue
Definition CrossTalk.C:17
const double gb_time_low
Definition CrossTalk.C:13
const double gg_time_high
Definition CrossTalk.C:10
bool PromptCoincidence(TGriffinHit *one, TGriffinHit *two)
Definition CrossTalk.C:24
const double gg_time_low
Definition CrossTalk.C:9
const double gb_time_high
Definition CrossTalk.C:14
bool pileup_reject
const Short_t defaultKValue
bool Addback(TGriffinHit *one, TGriffinHit *two)
TGriffinBgo * fGriffinBgo
Definition CrossTalk.h:30
void FillHistograms()
Definition CrossTalk.C:68
TGriffin * fGrif
Definition CrossTalk.h:28
void CreateHistograms()
Definition CrossTalk.C:30
virtual double GetEnergy(Option_t *opt="") const
virtual Short_t GetKValue() const
!
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!
TGRSIMap< std::string, TH1 * > fH1
TGRSIMap< std::string, TH2 * > fH2
TGRSIMap< std::string, THnSparseF * > fHSparse
TList * GetOutputList() const override
this does the same as TSelector::GetOutputList()
UShort_t GetArrayNumber() const override
!
Definition TGriffinHit.h:63
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