File indexing completed on 2025-08-05 08:18:11
0001 #include "PHG4TruthEventAction.h"
0002
0003 #include "PHG4Hit.h"
0004 #include "PHG4HitContainer.h"
0005 #include "PHG4HitDefs.h" // for keytype
0006 #include "PHG4Particle.h" // for PHG4Particle
0007 #include "PHG4Shower.h"
0008 #include "PHG4TruthInfoContainer.h"
0009 #include "PHG4UserPrimaryParticleInformation.h"
0010
0011 #include <phool/PHCompositeNode.h>
0012 #include <phool/PHIODataNode.h> // for PHIODataNode
0013 #include <phool/PHNode.h>
0014 #include <phool/PHNodeIterator.h> // for PHNodeIterator
0015 #include <phool/PHPointerListIterator.h>
0016 #include <phool/getClass.h>
0017 #include <phool/phool.h> // for PHWHERE
0018
0019 #include <Geant4/G4Event.hh>
0020 #include <Geant4/G4PrimaryParticle.hh> // for G4PrimaryPart...
0021 #include <Geant4/G4PrimaryVertex.hh> // for G4PrimaryVertex
0022 #include <Geant4/G4VUserPrimaryParticleInformation.hh> // for G4VUserPrimar...
0023
0024
0025 #pragma GCC diagnostic push
0026 #pragma GCC diagnostic ignored "-Wshadow"
0027 #include <Eigen/Dense>
0028 #pragma GCC diagnostic pop
0029
0030 #include <cmath> // for isnan
0031 #include <cstdlib> // for abs
0032 #include <iostream> // for operator<<, endl
0033 #include <iterator> // for reverse_iterator
0034 #include <map>
0035 #include <string> // for operator==
0036 #include <utility> // for swap, pair
0037 #include <vector> // for vector, vecto...
0038
0039 class PHG4VtxPoint;
0040
0041 using namespace std;
0042
0043
0044 PHG4TruthEventAction::PHG4TruthEventAction()
0045 : m_TruthInfoContainer(nullptr)
0046 , m_LowerKeyPrevExist(0)
0047 , m_UpperKeyPrevExist(0)
0048 {
0049 }
0050
0051
0052 void PHG4TruthEventAction::BeginOfEventAction(const G4Event* )
0053 {
0054
0055 if (!m_TruthInfoContainer)
0056 {
0057 std::cout << "PHG4TruthEventAction::EndOfEventAction - unable to find G4TruthInfo node" << std::endl;
0058 return;
0059 }
0060
0061 const PHG4TruthInfoContainer::Map& map = m_TruthInfoContainer->GetMap();
0062 if (!map.empty())
0063 {
0064 m_LowerKeyPrevExist = map.begin()->first;
0065 m_UpperKeyPrevExist = map.rbegin()->first;
0066 }
0067 }
0068
0069
0070 void PHG4TruthEventAction::EndOfEventAction(const G4Event* evt)
0071 {
0072
0073 if (!m_TruthInfoContainer)
0074 {
0075 std::cout << "PHG4TruthEventAction::EndOfEventAction - unable to find G4TruthInfo node" << std::endl;
0076 return;
0077 }
0078
0079
0080
0081 PruneShowers();
0082 ProcessShowers();
0083
0084
0085
0086 std::set<int> savelist;
0087 std::set<int> savevtxlist;
0088
0089 for (int mytrkid : m_WriteSet)
0090 {
0091 std::vector<int> wrttracks;
0092 std::vector<int> wrtvtx;
0093
0094
0095 PHG4Particle* particle = m_TruthInfoContainer->GetParticle(mytrkid);
0096
0097
0098 if (savelist.find(mytrkid) != savelist.end())
0099 {
0100 continue;
0101 }
0102 else
0103 {
0104 wrttracks.push_back(mytrkid);
0105 wrtvtx.push_back(particle->get_vtx_id());
0106 }
0107
0108
0109
0110 int parentid = particle->get_parent_id();
0111 while (savelist.find(parentid) == savelist.end() && parentid != 0)
0112 {
0113 particle = m_TruthInfoContainer->GetParticle(parentid);
0114 wrttracks.push_back(parentid);
0115 wrtvtx.push_back(particle->get_vtx_id());
0116 parentid = particle->get_parent_id();
0117 }
0118
0119
0120
0121 while (wrttracks.begin() != wrttracks.end())
0122 {
0123 savelist.insert(wrttracks.back());
0124 wrttracks.pop_back();
0125 }
0126
0127 while (wrtvtx.begin() != wrtvtx.end())
0128 {
0129 savevtxlist.insert(wrtvtx.back());
0130 wrtvtx.pop_back();
0131 }
0132 }
0133
0134
0135
0136
0137
0138
0139
0140 int removed[4] = {0};
0141 PHG4TruthInfoContainer::Range truth_range = m_TruthInfoContainer->GetParticleRange();
0142 PHG4TruthInfoContainer::Iterator truthiter = truth_range.first;
0143 while (truthiter != truth_range.second)
0144 {
0145 removed[0]++;
0146 int trackid = truthiter->first;
0147 if (savelist.find(trackid) == savelist.end())
0148 {
0149
0150
0151
0152
0153
0154
0155
0156
0157 if (((trackid < m_LowerKeyPrevExist) || (trackid > m_UpperKeyPrevExist)) && ((truthiter->second)->get_parent_id() != 0))
0158 {
0159 m_TruthInfoContainer->delete_particle(truthiter++);
0160 removed[1]++;
0161 }
0162 else
0163 {
0164
0165
0166 savevtxlist.insert((truthiter->second)->get_vtx_id());
0167 ++truthiter;
0168 }
0169 }
0170 else
0171 {
0172
0173 ++truthiter;
0174 }
0175 }
0176
0177 PHG4TruthInfoContainer::VtxRange vtxrange = m_TruthInfoContainer->GetVtxRange();
0178 PHG4TruthInfoContainer::VtxIterator vtxiter = vtxrange.first;
0179 while (vtxiter != vtxrange.second)
0180 {
0181 removed[2]++;
0182 if (savevtxlist.find(vtxiter->first) == savevtxlist.end())
0183 {
0184 m_TruthInfoContainer->delete_vtx(vtxiter++);
0185 removed[3]++;
0186 }
0187 else
0188 {
0189 ++vtxiter;
0190 }
0191 }
0192
0193
0194
0195 G4PrimaryVertex* pvtx = evt->GetPrimaryVertex();
0196 while (pvtx)
0197 {
0198 G4PrimaryParticle* part = pvtx->GetPrimary();
0199 while (part)
0200 {
0201 PHG4UserPrimaryParticleInformation* userdata = dynamic_cast<PHG4UserPrimaryParticleInformation*>(part->GetUserInformation());
0202 if (userdata)
0203 {
0204 if (userdata->get_embed())
0205 {
0206 m_TruthInfoContainer->AddEmbededTrkId(userdata->get_user_track_id(), userdata->get_embed());
0207 m_TruthInfoContainer->AddEmbededVtxId(userdata->get_user_vtx_id(), userdata->get_embed());
0208 }
0209 }
0210 part = part->GetNext();
0211 }
0212 pvtx = pvtx->GetNext();
0213 }
0214
0215 return;
0216 }
0217
0218
0219 void PHG4TruthEventAction::AddTrackidToWritelist(const int trackid)
0220 {
0221 m_WriteSet.insert(trackid);
0222 }
0223
0224
0225 void PHG4TruthEventAction::SetInterfacePointers(PHCompositeNode* topNode)
0226 {
0227
0228 m_TruthInfoContainer = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0229
0230
0231 if (!m_TruthInfoContainer)
0232 {
0233 std::cout << "PHG4TruthEventAction::SetInterfacePointers - unable to find G4TruthInfo" << std::endl;
0234 }
0235
0236 SearchNode(topNode);
0237 }
0238
0239
0240 void PHG4TruthEventAction::SearchNode(PHCompositeNode* top)
0241 {
0242
0243
0244
0245
0246 PHNodeIterator nodeiter(top);
0247 PHPointerListIterator<PHNode> iter(nodeiter.ls());
0248 PHNode* thisNode;
0249 while ((thisNode = iter()))
0250 {
0251 if (thisNode->getType() == "PHCompositeNode")
0252 {
0253 SearchNode(static_cast<PHCompositeNode*>(thisNode));
0254 }
0255 else if (thisNode->getType() == "PHIODataNode")
0256 {
0257 if (thisNode->getName().find("G4HIT_") == 0)
0258 {
0259 PHIODataNode<PHG4HitContainer>* DNode = static_cast<PHIODataNode<PHG4HitContainer>*>(thisNode);
0260 if (DNode)
0261 {
0262 PHG4HitContainer* object = dynamic_cast<PHG4HitContainer*>(DNode->getData());
0263 if (object)
0264 {
0265 m_HitContainerMap[object->GetID()] = object;
0266 }
0267 }
0268 }
0269 }
0270 }
0271 }
0272
0273 int PHG4TruthEventAction::ResetEvent(PHCompositeNode* )
0274 {
0275 m_WriteSet.clear();
0276 return 0;
0277 }
0278
0279 void PHG4TruthEventAction::PruneShowers()
0280 {
0281 PHG4TruthInfoContainer::ShowerRange range = m_TruthInfoContainer->GetShowerRange();
0282 for (PHG4TruthInfoContainer::ShowerIterator iter = range.first;
0283 iter != range.second;
0284 ++iter)
0285 {
0286 PHG4Shower* shower = iter->second;
0287
0288 std::set<int> remove_ids;
0289 for (PHG4Shower::ParticleIdIter jter = shower->begin_g4particle_id();
0290 jter != shower->end_g4particle_id();
0291 ++jter)
0292 {
0293 int g4particle_id = *jter;
0294 PHG4Particle* particle = m_TruthInfoContainer->GetParticle(g4particle_id);
0295 if (!particle)
0296 {
0297 remove_ids.insert(g4particle_id);
0298 continue;
0299 }
0300 }
0301
0302 for (int remove_id : remove_ids)
0303 {
0304 shower->remove_g4particle_id(remove_id);
0305 }
0306
0307 std::set<int> remove_more_ids;
0308 for (std::map<int, std::set<PHG4HitDefs::keytype> >::iterator jter = shower->begin_g4hit_id();
0309 jter != shower->end_g4hit_id();
0310 ++jter)
0311 {
0312 int g4hitmap_id = jter->first;
0313 std::map<int, PHG4HitContainer*>::iterator mapiter = m_HitContainerMap.find(g4hitmap_id);
0314 if (mapiter == m_HitContainerMap.end())
0315 {
0316 continue;
0317 }
0318
0319
0320 for (std::set<PHG4HitDefs::keytype>::iterator kter = jter->second.begin();
0321 kter != jter->second.end();)
0322 {
0323 PHG4HitDefs::keytype g4hit_id = *kter;
0324
0325 PHG4Hit* g4hit = mapiter->second->findHit(g4hit_id);
0326 if (!g4hit)
0327 {
0328
0329 jter->second.erase(kter++);
0330 continue;
0331 }
0332 else
0333 {
0334 ++kter;
0335 }
0336 }
0337
0338 if (jter->second.empty())
0339 {
0340 remove_more_ids.insert(g4hitmap_id);
0341 }
0342 }
0343
0344 for (int remove_more_id : remove_more_ids)
0345 {
0346 shower->remove_g4hit_volume(remove_more_id);
0347 }
0348 }
0349
0350 range = m_TruthInfoContainer->GetShowerRange();
0351 for (PHG4TruthInfoContainer::ShowerIterator iter = range.first;
0352 iter != range.second;)
0353 {
0354 PHG4Shower* shower = iter->second;
0355
0356 if (shower->empty_g4particle_id() && shower->empty_g4hit_id())
0357 {
0358 if (shower->get_edep() == 0)
0359 {
0360 m_TruthInfoContainer->delete_shower(iter++);
0361 continue;
0362 }
0363 }
0364
0365 ++iter;
0366 }
0367 }
0368
0369 void PHG4TruthEventAction::ProcessShowers()
0370 {
0371 PHG4TruthInfoContainer::ShowerRange range = m_TruthInfoContainer->GetShowerRange();
0372 for (PHG4TruthInfoContainer::ShowerIterator shwiter = range.first;
0373 shwiter != range.second;
0374 ++shwiter)
0375 {
0376 PHG4Shower* shower = shwiter->second;
0377
0378
0379 std::vector<std::vector<float> > points;
0380 std::vector<float> weights;
0381 float sumw = 0.0;
0382 float sumw2 = 0.0;
0383
0384 for (std::map<int, std::set<PHG4HitDefs::keytype> >::iterator iter = shower->begin_g4hit_id();
0385 iter != shower->end_g4hit_id();
0386 ++iter)
0387 {
0388 int g4hitmap_id = iter->first;
0389 std::map<int, PHG4HitContainer*>::iterator mapiter = m_HitContainerMap.find(g4hitmap_id);
0390 if (mapiter == m_HitContainerMap.end())
0391 {
0392 continue;
0393 }
0394
0395 PHG4HitContainer* hits = mapiter->second;
0396
0397 unsigned int nhits = 0;
0398 float edep = 0.0;
0399 float eion = 0.0;
0400 float light_yield = 0.0;
0401 float edep_e = 0.0;
0402 float edep_h = 0.0;
0403
0404
0405 for (unsigned long long g4hit_id : iter->second)
0406 {
0407 PHG4Hit* g4hit = hits->findHit(g4hit_id);
0408 if (!g4hit)
0409 {
0410 cout << PHWHERE << " missing g4hit" << endl;
0411 continue;
0412 }
0413
0414 PHG4Particle* particle = m_TruthInfoContainer->GetParticle(g4hit->get_trkid());
0415 if (!particle)
0416 {
0417 cout << PHWHERE << " missing g4particle for track "
0418 << g4hit->get_trkid() << endl;
0419 continue;
0420 }
0421
0422 PHG4VtxPoint* vtx = m_TruthInfoContainer->GetVtx(particle->get_vtx_id());
0423 if (!vtx)
0424 {
0425 cout << PHWHERE << " missing g4vertex" << endl;
0426 continue;
0427 }
0428
0429
0430
0431 if (!isnan(g4hit->get_x(0)) &&
0432 !isnan(g4hit->get_y(0)) &&
0433 !isnan(g4hit->get_z(0)))
0434 {
0435 std::vector<float> entry(3);
0436 entry[0] = g4hit->get_x(0);
0437 entry[1] = g4hit->get_y(0);
0438 entry[2] = g4hit->get_z(0);
0439
0440 points.push_back(entry);
0441 float w = g4hit->get_edep();
0442 weights.push_back(w);
0443 sumw += w;
0444 sumw2 += w * w;
0445 }
0446
0447 if (!isnan(g4hit->get_x(1)) &&
0448 !isnan(g4hit->get_y(1)) &&
0449 !isnan(g4hit->get_z(1)))
0450 {
0451 std::vector<float> entry(3);
0452 entry[0] = g4hit->get_x(1);
0453 entry[1] = g4hit->get_y(1);
0454 entry[2] = g4hit->get_z(1);
0455
0456 points.push_back(entry);
0457 float w = g4hit->get_edep();
0458 weights.push_back(w);
0459 sumw += w;
0460 sumw2 += w * w;
0461 }
0462
0463
0464
0465 if (!isnan(g4hit->get_edep()))
0466 {
0467 if (abs(particle->get_pid()) == 11)
0468 {
0469 edep_e += g4hit->get_edep();
0470 }
0471 else
0472 {
0473 edep_h += g4hit->get_edep();
0474 }
0475 }
0476
0477
0478
0479 ++nhits;
0480 if (!isnan(g4hit->get_edep()))
0481 {
0482 edep += g4hit->get_edep();
0483 }
0484 if (!isnan(g4hit->get_eion()))
0485 {
0486 eion += g4hit->get_eion();
0487 }
0488 if (!isnan(g4hit->get_light_yield()))
0489 {
0490 light_yield += g4hit->get_light_yield();
0491 }
0492 }
0493
0494
0495
0496 if (nhits)
0497 {
0498 shower->set_nhits(g4hitmap_id, nhits);
0499 }
0500 if (edep != 0.0)
0501 {
0502 shower->set_edep(g4hitmap_id, edep);
0503 }
0504 if (eion != 0.0)
0505 {
0506 shower->set_eion(g4hitmap_id, eion);
0507 }
0508 if (light_yield != 0.0)
0509 {
0510 shower->set_light_yield(g4hitmap_id, light_yield);
0511 }
0512 if (edep_h != 0.0)
0513 {
0514 shower->set_eh_ratio(g4hitmap_id, edep_e / edep_h);
0515 }
0516 }
0517
0518
0519
0520
0521 Eigen::Matrix<double, Eigen::Dynamic, 3> X(points.size(), 3);
0522 Eigen::Matrix<double, Eigen::Dynamic, 1> W(weights.size(), 1);
0523
0524 for (unsigned int i = 0; i < points.size(); ++i)
0525 {
0526 for (unsigned int j = 0; j < 3; ++j)
0527 {
0528 X(i, j) = points[i][j];
0529 }
0530 W(i, 0) = weights[i];
0531 }
0532
0533
0534 double prefactor = 1.0 / sumw;
0535 Eigen::Matrix<double, 1, 3> mean = prefactor * W.transpose() * X;
0536
0537
0538 for (unsigned int i = 0; i < points.size(); ++i)
0539 {
0540 for (unsigned int j = 0; j < 3; ++j)
0541 {
0542 X(i, j) = points[i][j] - mean(0, j);
0543 }
0544 }
0545
0546
0547 prefactor = sumw / (sumw * sumw - sumw2);
0548 Eigen::Matrix<double, 3, 3> covar = prefactor * (X.transpose() * W.asDiagonal() * X);
0549
0550 shower->set_x(mean(0, 0));
0551 shower->set_y(mean(0, 1));
0552 shower->set_z(mean(0, 2));
0553
0554 for (unsigned int i = 0; i < 3; ++i)
0555 {
0556 for (unsigned int j = 0; j <= i; ++j)
0557 {
0558 shower->set_covar(i, j, covar(i, j));
0559 }
0560 }
0561
0562
0563 }
0564 }