Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4HitNtuple.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4HitNtuple.cc
1 #include "G4HitNtuple.h"
2 
3 #include <g4main/PHG4Hit.h>
5 
7 #include <fun4all/SubsysReco.h> // for SubsysReco
8 
9 #include <phool/getClass.h>
10 
11 #include <TFile.h>
12 #include <TH1.h>
13 #include <TNtuple.h>
14 
15 #include <sstream>
16 #include <utility> // for pair
17 
18 using namespace std;
19 
21  : SubsysReco(name)
22  , nblocks(0)
23  , hm(nullptr)
24  , _filename(filename)
25  , ntup(nullptr)
26  , outfile(nullptr)
27 {
28 }
29 
31 {
32  // delete ntup;
33  delete hm;
34 }
35 
37 {
38  hm = new Fun4AllHistoManager(Name());
39  outfile = new TFile(_filename.c_str(), "RECREATE");
40  ntup = new TNtuple("hitntup", "G4Hits", "detid:lyr:slat:x0:y0:z0:x1:y1:z1:edep");
41  // ntup->SetDirectory(0);
42  TH1 *h1 = new TH1F("edep1GeV", "edep 0-1GeV", 1000, 0, 1);
43  eloss.push_back(h1);
44  h1 = new TH1F("edep100GeV", "edep 0-100GeV", 1000, 0, 100);
45  eloss.push_back(h1);
46  return 0;
47 }
48 
50 {
51  ostringstream nodename;
52  set<string>::const_iterator iter;
53  vector<TH1 *>::const_iterator eiter;
54  for (iter = _node_postfix.begin(); iter != _node_postfix.end(); ++iter)
55  {
56  int detid = (_detid.find(*iter))->second;
57  nodename.str("");
58  nodename << "G4HIT_" << *iter;
59  PHG4HitContainer *hits = findNode::getClass<PHG4HitContainer>(topNode, nodename.str());
60  if (hits)
61  {
62  double esum = 0;
63  // double numhits = hits->size();
64  // nhits[i]->Fill(numhits);
65  PHG4HitContainer::ConstRange hit_range = hits->getHits();
66  for (PHG4HitContainer::ConstIterator hit_iter = hit_range.first; hit_iter != hit_range.second; hit_iter++)
67 
68  {
69  esum += hit_iter->second->get_edep();
70  ntup->Fill(detid,
71  hit_iter->second->get_layer(),
72  hit_iter->second->get_scint_id(),
73  hit_iter->second->get_x(0),
74  hit_iter->second->get_y(0),
75  hit_iter->second->get_z(0),
76  hit_iter->second->get_x(1),
77  hit_iter->second->get_y(1),
78  hit_iter->second->get_z(1),
79  hit_iter->second->get_edep());
80  }
81  for (eiter = eloss.begin(); eiter != eloss.end(); ++eiter)
82  {
83  (*eiter)->Fill(esum);
84  }
85  }
86  }
87  return 0;
88 }
89 
91 {
92  outfile->cd();
93  ntup->Write();
94  outfile->Write();
95  outfile->Close();
96  delete outfile;
97  hm->dumpHistos(_filename, "UPDATE");
98  return 0;
99 }
100 
101 void G4HitNtuple::AddNode(const std::string &name, const int detid)
102 {
103  _node_postfix.insert(name);
104  _detid[name] = detid;
105  return;
106 }