Class Reference for E1039 Core & Analysis Software
AnaSimHit.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <memory>
4 #include <g4main/PHG4Particle.h>
9 #include <phool/getClass.h>
10 #include <UtilAna/UtilSQHit.h>
11 #include "AnaSimHit.h"
12 using namespace std;
13 
15  : SubsysReco("AnaSimHit")
16  , m_evt (0)
17  , m_evt_true(0)
18  , m_hit_vec (0)
19  , m_g4true (0)
20 {
21  ;
22 }
23 
25 {
27 }
28 
30 {
31  //m_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
32  //m_evt_true = findNode::getClass<SQMCEvent>(topNode, "SQMCEvent");
33  m_hit_vec = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
34  m_g4true = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
36 }
37 
39 {
40  static unsigned int n_evt = 0;
41  cout << "AnaSimHit: Event " << ++n_evt << "\n";
42  if (!m_hit_vec || !m_g4true) return Fun4AllReturnCodes::EVENT_OK;
43 
44  shared_ptr<SQHitVector> hv_h1t(UtilSQHit::FindHits(m_hit_vec, "H1T"));
45  shared_ptr<SQHitVector> hv_h1b(UtilSQHit::FindHits(m_hit_vec, "H1B"));
46  shared_ptr<SQHitVector> hv_h2t(UtilSQHit::FindHits(m_hit_vec, "H2T"));
47  shared_ptr<SQHitVector> hv_h2b(UtilSQHit::FindHits(m_hit_vec, "H2B"));
48  shared_ptr<SQHitVector> hv_h3t(UtilSQHit::FindHits(m_hit_vec, "H3T"));
49  shared_ptr<SQHitVector> hv_h3b(UtilSQHit::FindHits(m_hit_vec, "H3B"));
50  shared_ptr<SQHitVector> hv_h4t(UtilSQHit::FindHits(m_hit_vec, "H4T"));
51  shared_ptr<SQHitVector> hv_h4b(UtilSQHit::FindHits(m_hit_vec, "H4B"));
52  PrintHits("H1T", hv_h1t.get());
53  PrintHits("H1B", hv_h1b.get());
54  PrintHits("H2T", hv_h2t.get());
55  PrintHits("H2B", hv_h2b.get());
56  PrintHits("H3T", hv_h3t.get());
57  PrintHits("H3B", hv_h3b.get());
58  PrintHits("H4T", hv_h4t.get());
59  PrintHits("H4B", hv_h4b.get());
60 
61  int n_h1_pri = 0;
62  int n_h1_sec = 0;
63  CountHits(hv_h1t.get(), n_h1_pri, n_h1_sec);
64  CountHits(hv_h1b.get(), n_h1_pri, n_h1_sec);
65  int n_h2_pri = 0;
66  int n_h2_sec = 0;
67  CountHits(hv_h2t.get(), n_h2_pri, n_h2_sec);
68  CountHits(hv_h2b.get(), n_h2_pri, n_h2_sec);
69  int n_h3_pri = 0;
70  int n_h3_sec = 0;
71  CountHits(hv_h3t.get(), n_h3_pri, n_h3_sec);
72  CountHits(hv_h3b.get(), n_h3_pri, n_h3_sec);
73  int n_h4_pri = 0;
74  int n_h4_sec = 0;
75  CountHits(hv_h4t.get(), n_h4_pri, n_h4_sec);
76  CountHits(hv_h4b.get(), n_h4_pri, n_h4_sec);
77  cout << " Primary : " << n_h1_pri << " " << n_h2_pri << " " << n_h3_pri << " " << n_h4_pri << "\n"
78  << " Secondary: " << n_h1_sec << " " << n_h2_sec << " " << n_h3_sec << " " << n_h4_sec << "\n";
79 
81 }
82 
84 {
86 }
87 
88 int AnaSimHit::GetParticleID(const int track_id)
89 {
90  for (auto it = m_g4true->GetParticleRange().first; it != m_g4true->GetParticleRange().second; ++it) {
91  PHG4Particle* par = it->second;
92  if (par->get_track_id() == track_id) return par->get_pid();
93  }
94  cout << "!!ERROR!! AnaSimHit::GetParticleID(): Failed for track_id = " << track_id << ". Abort.\n";
95  exit(1);
96 }
97 
98 void AnaSimHit::PrintHits(const std::string name, SQHitVector* hv, std::ostream& os)
99 {
100  os << " " << name;
101  for (SQHitVector::ConstIter it = hv->begin(); it != hv->end(); it++) {
102  int trk_id = (*it)->get_track_id();
103  os << " " << trk_id << ":" << GetParticleID(trk_id);
104  }
105  os << endl;
106 }
107 
108 void AnaSimHit::CountHits(SQHitVector* hv, int& n_pri, int& n_sec)
109 {
110  for (SQHitVector::ConstIter it = hv->begin(); it != hv->end(); it++) {
111  if ((*it)->get_track_id() >= 0) n_pri++;
112  else n_sec++;
113  }
114 }
int process_event(PHCompositeNode *topNode)
Definition: AnaSimHit.cc:38
int InitRun(PHCompositeNode *topNode)
Definition: AnaSimHit.cc:29
int Init(PHCompositeNode *topNode)
Definition: AnaSimHit.cc:24
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: AnaSimHit.cc:83
virtual int get_track_id() const
Definition: PHG4Particle.h:21
virtual int get_pid() const
Definition: PHG4Particle.h:14
Range GetParticleRange()
Get a range of iterators covering the entire container.
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
std::vector< SQHit * >::const_iterator ConstIter
Definition: SQHitVector.h:37
virtual ConstIter end() const =0
virtual ConstIter begin() const =0
SQHitVector * FindHits(const SQHitVector *vec_in, const std::string det_name, const bool in_time=false)
Extract a set of hits that are of the given detector (det_name).
Definition: UtilSQHit.cc:16