forked from sPHENIX-Collaboration/coresoftware
/
G4CellNtuple.cc
150 lines (137 loc) · 4.51 KB
/
G4CellNtuple.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
#include "G4CellNtuple.h"
#include <g4detectors/PHG4Cell.h>
#include <g4detectors/PHG4CellContainer.h>
#include <g4detectors/PHG4CellDefs.h> // for get_phibin
#include <g4detectors/PHG4CylinderCellGeom.h>
#include <g4detectors/PHG4CylinderCellGeomContainer.h>
#include <fun4all/Fun4AllHistoManager.h>
#include <fun4all/SubsysReco.h> // for SubsysReco
#include <phool/getClass.h>
#include <TFile.h>
#include <TH1.h>
#include <TNtuple.h>
#include <cassert>
#include <cmath> // for isfinite
#include <iostream> // for operator<<
#include <sstream>
#include <utility> // for pair
using namespace std;
G4CellNtuple::G4CellNtuple(const std::string &name, const std::string &filename)
: SubsysReco(name)
, nblocks(0)
, hm(nullptr)
, _filename(filename)
, ntup(nullptr)
, outfile(nullptr)
{
}
G4CellNtuple::~G4CellNtuple()
{
// delete ntup;
delete hm;
}
int G4CellNtuple::Init(PHCompositeNode *)
{
hm = new Fun4AllHistoManager(Name());
outfile = new TFile(_filename.c_str(), "RECREATE");
ntup = new TNtuple("cellntup", "G4Cells", "detid:layer:phi:eta:edep");
// ntup->SetDirectory(0);
TH1 *h1 = new TH1F("edep1GeV", "edep 0-1GeV", 1000, 0, 1);
eloss.push_back(h1);
h1 = new TH1F("edep100GeV", "edep 0-100GeV", 1000, 0, 100);
eloss.push_back(h1);
return 0;
}
int G4CellNtuple::process_event(PHCompositeNode *topNode)
{
ostringstream nodename;
ostringstream geonodename;
set<string>::const_iterator iter;
vector<TH1 *>::const_iterator eiter;
for (iter = _node_postfix.begin(); iter != _node_postfix.end(); ++iter)
{
int detid = (_detid.find(*iter))->second;
nodename.str("");
nodename << "G4CELL_" << *iter;
geonodename.str("");
geonodename << "CYLINDERCELLGEOM_" << *iter;
PHG4CylinderCellGeomContainer *cellgeos = findNode::getClass<PHG4CylinderCellGeomContainer>(topNode, geonodename.str());
if (!cellgeos)
{
cout << "no geometry node " << geonodename.str() << " for " << *iter << endl;
continue;
}
PHG4CellContainer *cells = findNode::getClass<PHG4CellContainer>(topNode, nodename.str());
if (cells)
{
int previouslayer = -1;
PHG4CylinderCellGeom *cellgeom = nullptr;
double esum = 0;
// double numcells = cells->size();
// ncells[i]->Fill(numcells);
// cout << "number of cells: " << cells->size() << endl;
PHG4CellContainer::ConstRange cell_range = cells->getCells();
for (PHG4CellContainer::ConstIterator cell_iter = cell_range.first; cell_iter != cell_range.second; cell_iter++)
{
double edep = cell_iter->second->get_edep();
if (!isfinite(edep))
{
cout << "invalid edep: " << edep << endl;
}
esum += cell_iter->second->get_edep();
int phibin = ~0x0;
int etabin = ~0x0;
double phi = NAN;
double eta = NAN;
int layer = cell_iter->second->get_layer();
// to search the map fewer times, cache the geom object until the layer changes
if (layer != previouslayer)
{
cellgeom = cellgeos->GetLayerCellGeom(layer);
previouslayer = layer;
}
if (cell_iter->second->has_binning(PHG4CellDefs::etaphibinning))
{
phibin = PHG4CellDefs::EtaPhiBinning::get_phibin(cell_iter->second->get_cellid());
etabin = PHG4CellDefs::EtaPhiBinning::get_etabin(cell_iter->second->get_cellid());
phi = cellgeom->get_phicenter(phibin);
eta = cellgeom->get_etacenter(etabin);
}
else if (cell_iter->second->has_binning(PHG4CellDefs::sizebinning))
{
phibin = PHG4CellDefs::SizeBinning::get_phibin(cell_iter->second->get_cellid());
etabin = PHG4CellDefs::SizeBinning::get_zbin(cell_iter->second->get_cellid());
phi = cellgeom->get_phicenter(phibin);
eta = cellgeom->get_zcenter(etabin);
}
assert(cellgeom != nullptr);
ntup->Fill(detid,
cell_iter->second->get_layer(),
phi,
eta,
cell_iter->second->get_edep());
}
for (eiter = eloss.begin(); eiter != eloss.end(); ++eiter)
{
(*eiter)->Fill(esum);
}
}
}
return 0;
}
int G4CellNtuple::End(PHCompositeNode */*topNode*/)
{
outfile->cd();
ntup->Write();
outfile->Write();
outfile->Close();
delete outfile;
hm->dumpHistos(_filename, "UPDATE");
return 0;
}
void G4CellNtuple::AddNode(const std::string &name, const int detid)
{
_node_postfix.insert(name);
_detid[name] = detid;
return;
}