cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.3
A generic central exclusive processes event generator
Loading...
Searching...
No Matches
ROOTTreeInfo.h
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2019-2024 Laurent Forthomme
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#ifndef CepGen_EventInterfaces_ROOTTreeInfo_h
20#define CepGen_EventInterfaces_ROOTTreeInfo_h
21
22#include <TFile.h>
23#include <TTree.h>
24
25#include "CepGen/Event/Event.h"
26
27namespace ROOT {
29 class CepGenRun {
30 public:
31 static constexpr const char* TREE_NAME = "run";
32
33 static CepGenRun load(TFile*, const std::string& run_tree = TREE_NAME);
34 static CepGenRun load(const std::string&, const std::string& run_tree = TREE_NAME);
35
36 double sqrt_s{-1.};
37 double xsect{-1.};
38 double errxsect{-1.};
39 unsigned int num_events{0};
40 unsigned int litigious_events{0};
41 std::string process_name;
42 std::string process_parameters;
43
44 explicit CepGenRun();
45 void clear();
46 void create();
47 inline TTree* tree() { return tree_.get(); }
48 void fill();
49 void attach(TFile* file, const std::string& run_tree = TREE_NAME);
51 void attach(const std::string& filename, const std::string& run_tree = TREE_NAME) {
52 attach(TFile::Open(filename.data()), run_tree);
53 }
54
55 private:
56 std::shared_ptr<TTree> tree_;
57 };
58
61 public:
62 inline CepGenEvent() { clear(); }
63
64 // book a sufficiently large number to allow the large multiplicity of excited proton fragmentation products
65 static constexpr size_t MAX_PART = 5000;
66 static constexpr const char* TREE_NAME = "events";
67
68 static CepGenEvent load(TFile*, const std::string& events_tree = TREE_NAME);
69 static CepGenEvent load(const std::string&, const std::string& events_tree = TREE_NAME);
70
72 float gen_time{-1.};
73 float tot_time{-1.};
74 float weight{-1.};
75 int np{0};
76 double pt[MAX_PART];
77 double eta[MAX_PART];
78 double phi[MAX_PART];
80 double E[MAX_PART];
81 double m[MAX_PART];
82 double charge[MAX_PART];
89
90 void clear();
91 inline TTree* tree() { return tree_.get(); }
92 void create();
93
94 void attach();
95 void attach(TFile* f, const std::string& events_tree = TREE_NAME);
97 void attach(const std::string& filename, const std::string& events_tree = TREE_NAME);
98
99 //--- direct cepgen::Event I/O helpers
100
101 void fill(const cepgen::Event&, bool compress = false);
102 bool next(cepgen::Event&);
103
104 private:
105 std::shared_ptr<TTree> tree_;
106 std::unique_ptr<TFile> file_;
107 bool tree_attached_{false};
108 unsigned long long num_read_events_{0ull};
109 };
110} // namespace ROOT
111
112#endif
All useful information about a generated event.
double m[MAX_PART]
Particles mass, in GeV/c .
static CepGenEvent load(TFile *, const std::string &events_tree=TREE_NAME)
double rapidity[MAX_PART]
Particles rapidity.
double charge[MAX_PART]
Particles charges, in e.
bool next(cepgen::Event &)
Read the next event in the file.
double eta[MAX_PART]
Particles pseudo-rapidity.
int parent1[MAX_PART]
First particles mother.
void fill(const cepgen::Event &, bool compress=false)
Fill the tree with a new event.
cepgen::Event::EventMetadata metadata
int status[MAX_PART]
Integer status code.
double phi[MAX_PART]
Particles azimuthal angle.
TTree * tree()
Retrieve the ROOT tree.
float weight
Event weight.
double E[MAX_PART]
Particles energy, in GeV.
int np
Number of particles in the event.
float gen_time
Event generation time.
float tot_time
Total event generation time.
int parent2[MAX_PART]
Last particles mother.
int role[MAX_PART]
Particles role in the event.
static constexpr size_t MAX_PART
Maximal number of particles in event.
void attach()
Attach the event tree reader to a tree.
int pdg_id[MAX_PART]
Integer particles PDG id.
void clear()
Reinitialise the event content.
void create()
Populate the tree and all associated branches.
static constexpr const char * TREE_NAME
Output tree name.
int stable[MAX_PART]
Whether the particle must decay or not.
double pt[MAX_PART]
Particles transverse momentum.
All useful information about a generation run.
std::string process_name
Unique name of the process generated in this run.
unsigned int litigious_events
Number of litigious events in run.
double xsect
Process cross section, in pb.
double errxsect
Uncertainty on process cross section, in pb.
TTree * tree()
Retrieve the ROOT tree.
void attach(TFile *file, const std::string &run_tree=TREE_NAME)
Attach the run tree reader to a given tree Attach the run tree reader to a given file.
void clear()
Reinitialise the run tree.
std::string process_parameters
Serialised process parameters.
unsigned int num_events
Number of events generated in run.
static CepGenRun load(TFile *, const std::string &run_tree=TREE_NAME)
void create()
Populate the run tree.
static constexpr const char * TREE_NAME
Output tree name.
void fill()
Fill the run tree.
void attach(const std::string &filename, const std::string &run_tree=TREE_NAME)
double sqrt_s
Centre of mass energy for beam particles.
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
Collection of key -> value pairs storing event metadata.
Definition Event.h:123