GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
ExampleTreeSelector.h
Go to the documentation of this file.
1#ifndef ExampleTreeSelector_h
2#define ExampleTreeSelector_h
3
4/** \addtogroup Selectors
5 * @{
6 */
7
8////////////////////////////////////////////////////////////////////////////////
9/// \class ExampleTreeSelector
10///
11/// This selector shows how to create a tree with selected events (beta-tagged
12/// with gamma multiplicities of at least three), and selected information
13/// (suppressed addback energies, and beta-gamma timing differences).
14///
15////////////////////////////////////////////////////////////////////////////////
16
17#include "TChain.h"
18#include "TFile.h"
19
20#include "TH1.h"
21#include "TH2.h"
22#include "THnSparse.h"
23
24// Header file for the classes stored in the TTree if any.
25#include "TGRSISelector.h"
26#include "TGriffin.h"
27#include "TGriffinBgo.h"
28#include "TSceptar.h"
29#include "TZeroDegree.h"
30
31// Fixed size dimensions of array or collections stored in the TTree if any.
32
33class ExampleTreeSelector : public TGRSISelector { // Must be same name as .C and .h
34public:
35 explicit ExampleTreeSelector(TTree* /*tree*/ = nullptr) : TGRSISelector(), fGrif(nullptr), fScep(nullptr), fZds(nullptr), fGriffinBgo(nullptr)
36 {
37 SetOutputPrefix("ExampleTree"); // Changes prefix of output file
38 }
39 // These functions are expected to exist
40 virtual ~ExampleTreeSelector() = default;
41 virtual Int_t Version() const { return 1; }
42 void CreateHistograms();
43 void FillHistograms();
44 void InitializeBranches(TTree* tree);
45
46private:
47 // branches of input tree
48 TGriffin* fGrif; ///< pointer for griffin branch
49 TSceptar* fScep; ///< pointer for sceptar branch
50 TGriffinBgo* fGriffinBgo; ///< pointer for griffin BGO branch
51 TZeroDegree* fZds; ///< pointer for ZDS branch
52
53 // branches of output trees
54 std::vector<double> fSuppressedAddback; ///< vector of suppressed addback energies
55 std::vector<double> fBetaGammaTiming; ///< vector of beta-gamma timing
56
57 ClassDef(ExampleTreeSelector, 2); // Makes ROOT happier
58};
59/*! @} */
60#endif
61
62#ifdef ExampleTreeSelector_cxx
64{
65 if(!tree) return;
66 if(tree->SetBranchAddress("TGriffin", &fGrif) == TTree::kMissingBranch) {
67 fGrif = new TGriffin;
68 }
69 if(tree->SetBranchAddress("TSceptar", &fScep) == TTree::kMissingBranch) {
70 fScep = new TSceptar;
71 }
72 if(tree->SetBranchAddress("TZeroDegree", &fZds) == TTree::kMissingBranch) {
73 fZds = new TZeroDegree;
74 }
75 if(tree->SetBranchAddress("TGriffinBgo", &fGriffinBgo) == TTree::kMissingBranch) {
77 }
78}
79#endif // #ifdef ExampleTreeSelector_cxx
ClassDef(ExampleTreeSelector, 2)
TSceptar * fScep
pointer for sceptar branch
std::vector< double > fSuppressedAddback
vector of suppressed addback energies
TGriffinBgo * fGriffinBgo
pointer for griffin BGO branch
void InitializeBranches(TTree *tree)
ExampleTreeSelector(TTree *=nullptr)
std::vector< double > fBetaGammaTiming
vector of beta-gamma timing
TZeroDegree * fZds
pointer for ZDS branch
virtual ~ExampleTreeSelector()=default
TGriffin * fGrif
pointer for griffin branch
virtual Int_t Version() const
void SetOutputPrefix(const char *prefix)