GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
ExampleEventSelector.C
Go to the documentation of this file.
1#define ExampleEventSelector_cxx
2// The class definition in ExampleEventSelector.h has been generated automatically
4
6{
7 // Define Histograms
8 fH1["gE"] = new TH1D("gE", "#gamma Singles", 12000, 0, 3000);
9 fH1["gE_b"] = new TH1D("gE_b", "#gamma Singles in rough #beta coincidence", 12000, 0, 3000);
10 fH1["aE"] = new TH1D("aE", "Addback Singles", 12000, 0, 3000);
11 fH2["ggE"] = new TH2F("ggE", "#gamma #gamma Coincidence", 8192, 0, 4096., 8192, 0, 4096.);
12 fH1["gsE"] = new TH1D("gsE", "suppressed #gamma Singles", 12000, 0, 3000);
13 fH1["gsE_b"] = new TH1D("gsE_b", "suppressed #gamma Singles in rough #beta coincidence", 12000, 0, 3000);
14 fH1["asE"] = new TH1D("asE", "suppressed Addback Singles", 12000, 0, 3000);
15
16 // Send histograms to Output list to be added and written.
17 for(auto it : fH1) {
18 GetOutputList()->Add(it.second);
19 }
20 for(auto it : fH2) {
21 GetOutputList()->Add(it.second);
22 }
23 for(auto it : fHSparse) {
24 GetOutputList()->Add(it.second);
25 }
26}
27
29{
30 // Check if hits are less then 300 ns apart.
31 return std::fabs(g->GetTime() - s->GetTime()) < 300.;
32}
33
35{
36 // Check if hits are less then 500 ns apart.
37 return std::fabs(g1->GetTime() - g2->GetTime()) < 500.;
38}
39
41{
42 // Loop over all Griffin Hits
43 for(auto i = 0; i < fGrif->GetMultiplicity(); ++i) {
44 auto grif1 = fGrif->GetGriffinHit(i);
45 fH1.at("gE")->Fill(grif1->GetEnergy());
46 // second loop over all Griffin Hits
47 for(auto j = 0; j < fGrif->GetMultiplicity(); ++j) {
48 if(i == j) continue;
49 auto grif2 = fGrif->GetGriffinHit(j);
50 fH2.at("ggE")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
51 }
52 // Loop over all sceptar hits
53 for(auto j = 0; j < fScep->GetMultiplicity(); ++j) {
54 if(PromptCoincidence(grif1, fScep->GetSceptarHit(j))) {
55 fH1.at("gE_b")->Fill(grif1->GetEnergy());
56 }
57 }
58 }
59 // Loop over all addback Griffin Hits
60 for(auto i = 0; i < fGrif->GetAddbackMultiplicity(); ++i) {
61 auto grif1 = fGrif->GetAddbackHit(i);
62 fH1.at("aE")->Fill(grif1->GetEnergy());
63 }
64 // Loop over all suppressed Griffin Hits
65 for(auto i = 0; i < fGrif->GetSuppressedMultiplicity(fGriffinBgo); ++i) {
66 auto grif1 = fGrif->GetSuppressedHit(i);
67 fH1.at("gsE")->Fill(grif1->GetEnergy());
68 // Loop over all sceptar hits
69 for(auto j = 0; j < fScep->GetMultiplicity(); ++j) {
70 if(PromptCoincidence(grif1, fScep->GetSceptarHit(j))) {
71 fH1.at("gsE_b")->Fill(grif1->GetEnergy());
72 }
73 }
74 }
75 // Loop over all suppressed addback Griffin Hits
76 for(auto i = 0; i < fGrif->GetSuppressedAddbackMultiplicity(fGriffinBgo); ++i) {
77 auto grif1 = fGrif->GetSuppressedAddbackHit(i);
78 fH1.at("asE")->Fill(grif1->GetEnergy());
79 }
80}
bool PromptCoincidence(TGriffinHit *one, TGriffinHit *two)
bool PromptCoincidence(TGriffinHit *g, TSceptarHit *s)
const Double_t s
virtual Double_t GetTime(const ETimeFlag &correct_flag=ETimeFlag::kAll, Option_t *opt="") const
Returns a time value to the nearest nanosecond!
virtual Short_t GetMultiplicity() const
Definition TDetector.h:73
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()
TGriffinHit * GetAddbackHit(const int &i)
Definition TGriffin.h:87
Short_t GetAddbackMultiplicity()
Definition TGriffin.h:84
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
Short_t GetMultiplicity() const override
Definition TGriffin.h:52
TGriffinHit * GetSuppressedHit(const int &i)
!
Definition TGriffin.h:108
TGriffinHit * GetGriffinHit(const int &i)
!
Definition TGriffin.h:47
TSceptarHit * GetSceptarHit(const int &i) const
Definition TSceptar.h:36