Class Reference for E1039 Core & Analysis Software
PHG4TruthEventAction.cc
Go to the documentation of this file.
1 #include "PHG4TruthEventAction.h"
2 #include "PHG4Particlev2.h"
3 
5 
6 #include "PHG4VtxPoint.h"
7 #include "PHG4Shower.h"
8 
10 #include "PHG4HitContainer.h"
11 #include "PHG4Hit.h"
12 
13 #include <phool/getClass.h>
15 #include <phool/PHNode.h>
16 
17 #include <Geant4/G4Event.hh>
18 #include <Geant4/G4TrajectoryContainer.hh>
19 #include <Geant4/G4VTrajectory.hh>
20 #include <Geant4/G4VTrajectoryPoint.hh>
21 #include <Geant4/G4ParticleTable.hh>
22 #include <Geant4/G4ParticleDefinition.hh>
23 #include <Geant4/globals.hh>
24 
25 #include <map>
26 
27 #include <Eigen/Dense>
28 
29 using namespace std;
30 
31 //___________________________________________________
33  truthInfoList_( 0 ),
34  prev_existing_lower_key( 0 ),
35  prev_existing_upper_key( 0 )
36 {}
37 
38 //___________________________________________________
39 void PHG4TruthEventAction::BeginOfEventAction(const G4Event* evt) {
40 
41  // if we do not find the node we need to make it.
42  if ( !truthInfoList_ ) {
43  std::cout << "PHG4TruthEventAction::EndOfEventAction - unable to find G4TruthInfo node" << std::endl;
44  return;
45  }
46 
47  const PHG4TruthInfoContainer::Map& map = truthInfoList_->GetMap();
48  if (!map.empty()) {
49  prev_existing_lower_key = map.begin()->first;
50  prev_existing_upper_key = map.rbegin()->first;
51  }
52 }
53 
54 //___________________________________________________
55 void PHG4TruthEventAction::EndOfEventAction(const G4Event* evt) {
56 
57  // if we do not find the node we need to make it.
58  if ( !truthInfoList_ ) {
59  std::cout << "PHG4TruthEventAction::EndOfEventAction - unable to find G4TruthInfo node" << std::endl;
60  return;
61  }
62 
63  // construct a list of track ids to preserve in the the output that includes any
64  // track designated in the writeList_ during processing or its ancestry chain
65 
66  std::set<G4int> savelist;
67  std::set<int> savevtxlist;
68 
69  for (std::set<G4int>::const_iterator write_iter = writeList_.begin();
70  write_iter != writeList_.end();
71  ++write_iter) {
72 
73  std::vector<G4int> wrttracks;
74  std::vector<int> wrtvtx;
75 
76  // usertrackid
77  G4int mytrkid = *write_iter;
78  PHG4Particle *particle = truthInfoList_->GetParticle(mytrkid);
79 
80  // if track is already in save list, nothing needs to be done
81  if (savelist.find(mytrkid) != savelist.end()) {
82  continue;
83  } else {
84  wrttracks.push_back(mytrkid);
85  wrtvtx.push_back(particle->get_vtx_id());
86  }
87 
88  // now crawl up the truth info and add parents until we hit
89  // a track which is already being saved
90  G4int parentid = particle->get_parent_id();
91  while (savelist.find(parentid) == savelist.end() && parentid != 0) {
92  particle = truthInfoList_->GetParticle(parentid);
93  wrttracks.push_back(parentid);
94  wrtvtx.push_back(particle->get_vtx_id());
95  parentid = particle->get_parent_id();
96  }
97 
98  // now fill the new tracks into the save list
99  // while emptying the write lists
100  while (wrttracks.begin() != wrttracks.end()) {
101  savelist.insert(wrttracks.back());
102  wrttracks.pop_back();
103  }
104 
105  while (wrtvtx.begin() != wrtvtx.end()) {
106  savevtxlist.insert(wrtvtx.back());
107  wrtvtx.pop_back();
108  }
109  }
110 
111  // the save lists are filled now, except primary track which never
112  // made it into any active volume and their vertex
113 
114  // loop over particles in truth list and remove them if they are not
115  // in the save list and are not primary particles (parent_id == 0)
116 
117  int removed[4] = {0};
118  PHG4TruthInfoContainer::Range truth_range = truthInfoList_->GetParticleRange();
119  PHG4TruthInfoContainer::Iterator truthiter = truth_range.first;
120  while (truthiter != truth_range.second) {
121 
122  removed[0]++;
123  int trackid = truthiter->first;
124  if (savelist.find(trackid) == savelist.end()) {
125  // track not in save list
126 
127  // check if particle below offset
128  // primary particles had parentid = 0
129  // for embedding: particles from initial sim do not have their keep flag set, we want to keep particles with trkid <= trackidoffset
130  // tracks from a previous geant pass will not be recorded as leaving
131  // hit in the sim, so we exclude this range from the removal
132  // for regular sims, that range is zero to zero
133  if (((trackid < prev_existing_lower_key)||(trackid > prev_existing_upper_key)) && ((truthiter->second)->get_parent_id() != 0)) {
134  truthInfoList_->delete_particle(truthiter++);
135  removed[1]++;
136  } else {
137  // save vertex id for primary particle which leaves no hit
138  // in active area
139  savevtxlist.insert((truthiter->second)->get_vtx_id());
140  ++truthiter;
141  }
142 
143  } else {
144  // track was in save list, move on
145  ++truthiter;
146  }
147  }
148 
149  PHG4TruthInfoContainer::VtxRange vtxrange = truthInfoList_->GetVtxRange();
150  PHG4TruthInfoContainer::VtxIterator vtxiter = vtxrange.first;
151  while (vtxiter != vtxrange.second)
152  {
153  removed[2]++;
154  if (savevtxlist.find(vtxiter->first) == savevtxlist.end())
155  {
156  truthInfoList_->delete_vtx(vtxiter++);
157  removed[3]++;
158  }
159  else
160  {
161  ++vtxiter;
162  }
163  }
164 
165  // loop over all input particles and fish out the ones which have the embed flag set
166  // and store their geant track ids in truthinfo container
167  G4PrimaryVertex *pvtx = evt->GetPrimaryVertex();
168  while (pvtx) {
169 
170  G4PrimaryParticle *part = pvtx->GetPrimary();
171  while (part) {
172 
173  PHG4UserPrimaryParticleInformation *userdata = dynamic_cast<PHG4UserPrimaryParticleInformation *> (part->GetUserInformation());
174  if (userdata) {
175  if (userdata->get_embed()) {
176  truthInfoList_->AddEmbededTrkId(userdata->get_user_track_id(),userdata->get_embed());
177  truthInfoList_->AddEmbededVtxId(userdata->get_user_vtx_id(),userdata->get_embed());
178  }
179  }
180  part = part->GetNext();
181  }
182  pvtx = pvtx->GetNext();
183  }
184 
185  PruneShowers();
186  ProcessShowers();
187 
188  return;
189 }
190 
191 //___________________________________________________
193  writeList_.insert(trackid);
194 }
195 
196 //___________________________________________________
198 
199  //now look for the map and grab a pointer to it.
200  truthInfoList_ = findNode::getClass<PHG4TruthInfoContainer>( topNode , "G4TruthInfo" );
201 
202  // if we do not find the node we need to make it.
203  if ( !truthInfoList_ ) {
204  std::cout << "PHG4TruthEventAction::SetInterfacePointers - unable to find G4TruthInfo" << std::endl;
205  }
206 
207  SearchNode(topNode);
208 }
209 
210 void PHG4TruthEventAction::SearchNode(PHCompositeNode* top) {
211 
212  // fill a lookup map between the g4hit container ids and the containers themselves
213  // without knowing what the container names are in advance, only that they
214  // begin G4HIT_*
215 
216  PHNodeIterator nodeiter(top);
217  PHPointerListIterator<PHNode> iter(nodeiter.ls());
218  PHNode *thisNode;
219  while ((thisNode = iter())) {
220 
221  if (thisNode->getType() == "PHCompositeNode") {
222  SearchNode(static_cast<PHCompositeNode*>(thisNode) );
223  } else if (thisNode->getType() == "PHIODataNode") {
224  if (thisNode->getName().find("G4HIT_") == 0) {
225  PHIODataNode<PHG4HitContainer> *DNode = static_cast<PHIODataNode<PHG4HitContainer>*>(thisNode);
226  if (DNode) {
227  PHG4HitContainer* object = dynamic_cast<PHG4HitContainer*>(DNode->getData());
228  if (object) {
229  hitmap_[object->GetID()] = object;
230  }
231  }
232  }
233  }
234  }
235 }
236 
238  writeList_.clear();
239  return 0;
240 }
241 
242 void PHG4TruthEventAction::PruneShowers() {
243 
244  PHG4TruthInfoContainer::ShowerRange range = truthInfoList_->GetShowerRange();
245  for (PHG4TruthInfoContainer::ShowerIterator iter = range.first;
246  iter != range.second;
247  ++iter) {
248  PHG4Shower* shower = iter->second;
249 
250  std::set<int> remove_ids;
251  for (PHG4Shower::ParticleIdIter jter = shower->begin_g4particle_id();
252  jter != shower->end_g4particle_id();
253  ++jter) {
254  int g4particle_id = *jter;
255  PHG4Particle* particle = truthInfoList_->GetParticle(g4particle_id);
256  if (!particle) {
257  remove_ids.insert(g4particle_id);
258  continue;
259  }
260  }
261 
262  for (std::set<int>::iterator jter = remove_ids.begin();
263  jter != remove_ids.end();
264  ++jter) {
265  shower->remove_g4particle_id(*jter);
266  }
267 
268  std::set<int> remove_more_ids;
269  for (std::map<int,std::set<PHG4HitDefs::keytype> >::iterator jter = shower->begin_g4hit_id();
270  jter != shower->end_g4hit_id();
271  ++jter) {
272  int g4hitmap_id = jter->first;
273  std::map<int,PHG4HitContainer*>::iterator mapiter = hitmap_.find(g4hitmap_id);
274  if (mapiter == hitmap_.end()) {
275  continue;
276  }
277 
278  // get the g4hits from this particle in this volume
279  for (std::set<PHG4HitDefs::keytype>::iterator kter = jter->second.begin();
280  kter != jter->second.end();
281  ) {
282  PHG4HitDefs::keytype g4hit_id = *kter;
283 
284  PHG4Hit* g4hit = mapiter->second->findHit(g4hit_id);
285  if (!g4hit) {
286  // some zero edep g4hits have been removed already
287  jter->second.erase(kter++);
288  continue;
289  } else {
290  ++kter;
291  }
292  }
293 
294  if (jter->second.empty()) {
295  remove_more_ids.insert(g4hitmap_id);
296  }
297  }
298 
299  for (std::set<int>::iterator jter = remove_more_ids.begin();
300  jter != remove_more_ids.end();
301  ++jter) {
302  shower->remove_g4hit_volume(*jter);
303  }
304  }
305 
306  range = truthInfoList_->GetShowerRange();
307  for (PHG4TruthInfoContainer::ShowerIterator iter = range.first;
308  iter != range.second;
309  ) {
310  PHG4Shower* shower = iter->second;
311 
312  if (shower->empty_g4particle_id() && shower->empty_g4hit_id()) {
313  truthInfoList_->delete_shower(iter++);
314  continue;
315  }
316 
317  ++iter;
318  }
319 
320 }
321 
322 void PHG4TruthEventAction::ProcessShowers() {
323 
324  PHG4TruthInfoContainer::ShowerRange range = truthInfoList_->GetShowerRange();
325  for (PHG4TruthInfoContainer::ShowerIterator iter = range.first;
326  iter != range.second;
327  ++iter) {
328  PHG4Shower* shower = iter->second;
329 
330  // Data structures to hold weighted pca
331  std::vector<std::vector<float> > points;
332  std::vector<float> weights;
333  float sumw = 0.0;
334  float sumw2 = 0.0;
335 
336  for (std::map<int,std::set<PHG4HitDefs::keytype> >::iterator iter = shower->begin_g4hit_id();
337  iter != shower->end_g4hit_id();
338  ++iter) {
339  int g4hitmap_id = iter->first;
340  std::map<int,PHG4HitContainer*>::iterator mapiter = hitmap_.find(g4hitmap_id);
341  if (mapiter == hitmap_.end()) {
342  continue;
343  }
344 
345  PHG4HitContainer* hits = mapiter->second;
346 
347  unsigned int nhits = 0;
348  float edep = 0.0;
349  float eion = 0.0;
350  float light_yield = 0.0;
351  float edep_e = 0.0;
352  float edep_h = 0.0;
353 
354  // get the g4hits from this particle in this volume
355  for (std::set<PHG4HitDefs::keytype>::iterator kter = iter->second.begin();
356  kter != iter->second.end();
357  ++kter) {
358  PHG4HitDefs::keytype g4hit_id = *kter;
359 
360  PHG4Hit* g4hit = hits->findHit(g4hit_id);
361  if (!g4hit) {
362  cout << PHWHERE << " missing g4hit" << endl;
363  continue;
364  }
365 
366  PHG4Particle* particle = truthInfoList_->GetParticle(g4hit->get_trkid());
367  if (!particle) {
368  cout << PHWHERE << " missing g4particle for track "
369  << g4hit->get_trkid() << endl;
370  continue;
371  }
372 
373  PHG4VtxPoint* vtx = truthInfoList_->GetVtx( particle->get_vtx_id() );
374  if (!vtx) {
375  cout << PHWHERE << " missing g4vertex" << endl;
376  continue;
377  }
378 
379  // shower location and shape info
380 
381  if (!isnan(g4hit->get_x(0)) &&
382  !isnan(g4hit->get_y(0)) &&
383  !isnan(g4hit->get_z(0))) {
384 
385  std::vector<float> entry(3);
386  entry[0] = g4hit->get_x(0);
387  entry[1] = g4hit->get_y(0);
388  entry[2] = g4hit->get_z(0);
389 
390  points.push_back(entry);
391  float w = g4hit->get_edep();
392  weights.push_back(w);
393  sumw += w;
394  sumw2 += w*w;
395  }
396 
397  if (!isnan(g4hit->get_x(1)) &&
398  !isnan(g4hit->get_y(1)) &&
399  !isnan(g4hit->get_z(1))) {
400 
401  std::vector<float> entry(3);
402  entry[0] = g4hit->get_x(1);
403  entry[1] = g4hit->get_y(1);
404  entry[2] = g4hit->get_z(1);
405 
406  points.push_back(entry);
407  float w = g4hit->get_edep();
408  weights.push_back(w);
409  sumw += w;
410  sumw2 += w*w;
411  }
412 
413  // e/h ratio
414 
415  if (!isnan(g4hit->get_edep())) {
416  if (abs(particle->get_pid()) == 11) {
417  edep_e += g4hit->get_edep();
418  } else {
419  edep_h += g4hit->get_edep();
420  }
421  }
422 
423  // summary info
424 
425  if (g4hit) ++nhits;
426  if (!isnan(g4hit->get_edep())) edep += g4hit->get_edep();
427  if (!isnan(g4hit->get_eion())) eion += g4hit->get_eion();
428  if (!isnan(g4hit->get_light_yield())) light_yield += g4hit->get_light_yield();
429  } // g4hit loop
430 
431  // summary info
432 
433  if (nhits) shower->set_nhits(g4hitmap_id,nhits);
434  if (edep != 0.0) shower->set_edep(g4hitmap_id,edep);
435  if (eion != 0.0) shower->set_eion(g4hitmap_id,eion);
436  if (light_yield != 0.0) shower->set_light_yield(g4hitmap_id,light_yield);
437  if (edep_h != 0.0) shower->set_eh_ratio(g4hitmap_id,edep_e/edep_h);
438  } // volume loop
439 
440  // fill Eigen matrices to compute wPCA
441  // resizing these non-destructively is expensive
442  // so I fill vectors and then copy
443  Eigen::Matrix<double, Eigen::Dynamic, 3> X(points.size(),3);
444  Eigen::Matrix<double, Eigen::Dynamic, 1> W(weights.size(),1);
445 
446  for (unsigned int i=0; i<points.size(); ++i) {
447  for (unsigned int j=0; j<3; ++j) {
448  X(i,j) = points[i][j];
449  }
450  W(i,0) = weights[i];
451  }
452 
453  // mean value of shower
454  double prefactor = 1.0 / sumw;
455  Eigen::Matrix<double, 1, 3> mean = prefactor * W.transpose() * X;
456 
457  // compute residual relative to the mean
458  for (unsigned int i=0; i<points.size(); ++i) {
459  for (unsigned int j=0; j<3; ++j) X(i,j) = points[i][j] - mean(0,j);
460  }
461 
462  // weighted covariance matrix
463  prefactor = sumw / (pow(sumw,2) - sumw2); // effectivelly 1/(N-1) when w_i = 1.0
464  Eigen::Matrix<double, 3, 3> covar = prefactor * (X.transpose() * W.asDiagonal() * X);
465 
466  shower->set_x(mean(0,0));
467  shower->set_y(mean(0,1));
468  shower->set_z(mean(0,2));
469 
470  for (unsigned int i = 0; i < 3; ++i) {
471  for (unsigned int j = 0; j <= i; ++j) {
472  shower->set_covar(i,j,covar(i,j));
473  }
474  }
475 
476  // shower->identify();
477  }
478 }
T * getData()
Definition: PHDataNode.h:21
PHG4Hit * findHit(PHG4HitDefs::keytype key)
virtual float get_z(const int i) const
Definition: PHG4Hit.h:23
virtual float get_eion() const
Definition: PHG4Hit.h:32
virtual float get_edep() const
Definition: PHG4Hit.h:31
virtual float get_light_yield() const
Definition: PHG4Hit.h:33
virtual float get_y(const int i) const
Definition: PHG4Hit.h:22
virtual float get_x(const int i) const
Definition: PHG4Hit.h:21
virtual int get_trkid() const
Definition: PHG4Hit.h:41
virtual int get_pid() const
Definition: PHG4Particle.h:14
virtual int get_parent_id() const
Definition: PHG4Particle.h:23
virtual int get_vtx_id() const
Definition: PHG4Particle.h:22
virtual ParticleIdIter begin_g4particle_id()
Definition: PHG4Shower.h:81
ParticleIdSet::iterator ParticleIdIter
Definition: PHG4Shower.h:17
virtual void set_y(float y)
Definition: PHG4Shower.h:52
virtual void set_eh_ratio(int volume, float eh_ratio)
Definition: PHG4Shower.h:76
virtual void set_eion(int volume, float eion)
Definition: PHG4Shower.h:70
virtual bool empty_g4particle_id() const
Definition: PHG4Shower.h:78
virtual void set_x(float x)
Definition: PHG4Shower.h:49
virtual bool empty_g4hit_id() const
Definition: PHG4Shower.h:98
virtual void set_z(float x)
Definition: PHG4Shower.h:55
virtual size_t remove_g4hit_volume(int volume)
Definition: PHG4Shower.h:108
virtual size_t remove_g4particle_id(int id)
Definition: PHG4Shower.h:85
virtual void set_light_yield(int volume, float light_yield)
Definition: PHG4Shower.h:73
virtual void set_covar(unsigned int i, unsigned int j, float entry)
Definition: PHG4Shower.h:61
virtual HitIdIter begin_g4hit_id()
Definition: PHG4Shower.h:101
virtual ParticleIdIter end_g4particle_id()
Definition: PHG4Shower.h:83
virtual void set_nhits(int volume, unsigned int nhits)
Definition: PHG4Shower.h:64
virtual void set_edep(int volume, float edep)
Definition: PHG4Shower.h:67
virtual HitIdIter end_g4hit_id()
Definition: PHG4Shower.h:105
void BeginOfEventAction(const G4Event *)
int ResetEvent(PHCompositeNode *)
PHG4TruthEventAction(void)
constructor
void EndOfEventAction(const G4Event *)
void AddTrackidToWritelist(const G4int trackid)
add id into track list
void SetInterfacePointers(PHCompositeNode *)
get relevant nodes from top node passed as argument
VtxRange GetVtxRange()
Get a range of iterators covering the entire vertex container.
std::pair< ShowerIterator, ShowerIterator > ShowerRange
void AddEmbededTrkId(const int id, const int flag)
void delete_shower(ShowerIterator piter)
Range GetParticleRange()
Get a range of iterators covering the entire container.
std::pair< Iterator, Iterator > Range
PHG4VtxPoint * GetVtx(const int vtxid)
std::map< int, PHG4Particle * > Map
void delete_particle(Iterator piter)
std::pair< VtxIterator, VtxIterator > VtxRange
PHG4Particle * GetParticle(const int particleid)
ShowerRange GetShowerRange()
Get a range of iterators covering the entire container.
void delete_vtx(VtxIterator viter)
ShowerMap::iterator ShowerIterator
void AddEmbededVtxId(const int id, const int flag)
const Map & GetMap() const
Get the Particle Map storage.
Definition: PHNode.h:15
const std::string getType() const
Definition: PHNode.h:31
const std::string getName() const
Definition: PHNode.h:32
unsigned int keytype
Definition: PHG4HitDefs.h:8
#define PHWHERE
Definition: phool.h:23