File indexing completed on 2025-08-06 08:17:36
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "DecayFinder.h"
0012
0013 #include "DecayFinderContainerBase.h" // for DecayFinderContainerBase::Iter
0014 #include "DecayFinderContainer_v1.h" // for DecayFinderContainer_v1
0015
0016 #include <g4main/PHG4TruthInfoContainer.h>
0017
0018 #include <phhepmc/PHHepMCGenEvent.h>
0019 #include <phhepmc/PHHepMCGenEventMap.h>
0020
0021 #include <fun4all/Fun4AllReturnCodes.h>
0022
0023 #include <phool/PHCompositeNode.h>
0024 #include <phool/PHIODataNode.h> // for PHIODataNode
0025 #include <phool/PHNodeIterator.h> // for PHNodeIterator
0026 #include <phool/getClass.h>
0027
0028 #pragma GCC diagnostic push
0029 #pragma GCC diagnostic ignored "-Wdeprecated-declarations"
0030 #include <HepMC/GenEvent.h>
0031 #include <HepMC/GenVertex.h> // for GenVertex::particle_iterator
0032 #pragma GCC diagnostic pop
0033
0034 #include <HepMC/GenParticle.h>
0035 #include <HepMC/IteratorRange.h>
0036 #include <HepMC/SimpleVector.h>
0037
0038 #include <TDatabasePDG.h>
0039
0040 #include <algorithm> // for find, for_each, sort, transform
0041 #include <cctype> // for toupper
0042 #include <cstdlib> // for abs, size_t
0043 #include <iostream> // for operator<<, endl, basic_ostream
0044 #include <iterator> // for end, begin
0045 #include <map> // for map, map<>::mapped_type, _Rb...
0046 #include <memory> // for allocator_traits<>::value_type
0047
0048 int listOfResonantPIDs[] = {111, 113, 213, 333, 310, 311, 313, 323, 413, 423, 513, 523, 441, 443, 100443, 9000111, 9000211, 100111, 100211, 10111,
0049 10211, 9010111, 9010211, 10113, 10213, 20113, 20213, 9000113, 9000213, 100113, 100213, 9010113, 9010213, 9020113, 9020213,
0050 30113, 30213, 9030113, 9030213, 9040113, 9040213, 115, 215, 10115, 10215, 9000115, 9000215, 9010115, 9010215, 117, 217,
0051 9000117, 9000217, 9010117, 9010217, 119, 219, 221, 331, 9000221, 9010221, 100221, 10221, 9020221, 100331, 9030221, 10331,
0052 9040221, 9050221, 9060221, 9070221, 9080221, 223, 10223, 20223, 10333, 20333, 1000223, 9000223, 9010223, 30223, 100333, 225,
0053 9000225, 335, 9010225, 9020225, 10225, 9030225, 10335, 9040225, 9050225, 9060225, 9070225, 9080225, 9090225, 227, 337, 229,
0054 9000229, 9010229, 9000311, 9000321, 10311, 10321, 100311, 100321, 9010311, 9010321, 9020311, 9020321, 10313, 10323, 20313,
0055 20323, 100313, 100323, 9000313, 9000323, 30313, 30323, 315, 325, 9000315, 9000325, 10315, 10325, 20315, 20325, 9010315,
0056 9010325, 9020315, 9020325, 317, 327, 9010317, 9010327, 319, 329, 9000319, 9000329, 10411, 10421, 10413, 10423, 20413, 20423,
0057 415, 425, 431, 10431, 433, 10433, 20433, 435, 10511, 10521, 10513, 10523, 20513, 20523, 515, 525, 10531, 533, 10533, 20533, 535,
0058 10541, 543, 10543, 20543, 545, 10441, 100441, 10443, 20443, 30443, 9000443, 9010443, 9020443, 445, 100445, 551, 10551, 100551,
0059 110551, 200551, 210551, 10553, 20553, 30553, 110553, 120553, 130553, 210553, 220553, 9000553, 9010553, 555, 10555, 20555, 100555,
0060 110555, 120555, 200555, 557, 100557, 2224, 2214, 2114, 1114, 3212, 3224, 3214, 3114, 3324, 3314, 4222, 4212, 4112, 4224, 4214,
0061 4114, 4232, 4132, 4322, 4312, 4324, 4314, 4332, 4334, 4412, 4422, 4414, 4424, 4432, 4434, 4444};
0062
0063 DecayFinder::DecayFinder()
0064 : SubsysReco("DECAYFINDER")
0065 , m_save_dst(false)
0066 {
0067 }
0068
0069 DecayFinder::DecayFinder(const std::string& name)
0070 : SubsysReco(name)
0071 , m_save_dst(false)
0072 {
0073 }
0074
0075 int DecayFinder::Init(PHCompositeNode* topNode)
0076 {
0077 if (Verbosity() >= VERBOSITY_SOME)
0078 {
0079 std::cout << "DecayFinder name: " << Name() << std::endl;
0080 std::cout << "Decay descriptor: " << m_decayDescriptor << std::endl;
0081 }
0082
0083 int canSearchDecay = parseDecayDescriptor();
0084
0085 if (m_save_dst)
0086 {
0087 createDecayNode(topNode);
0088 }
0089
0090 return canSearchDecay;
0091 }
0092
0093 int DecayFinder::process_event(PHCompositeNode* topNode)
0094 {
0095 bool decayFound = findDecay(topNode);
0096
0097 if (decayFound && m_save_dst && Verbosity() >= VERBOSITY_MORE)
0098 {
0099 printNode(topNode);
0100 }
0101
0102 if (m_triggerOnDecay && !decayFound)
0103 {
0104 if (Verbosity() >= VERBOSITY_MORE)
0105 {
0106 std::cout << "The decay, " << m_decayDescriptor << " was not found or reconstructable in this event, skipping" << std::endl;
0107 }
0108 return Fun4AllReturnCodes::ABORTEVENT;
0109 }
0110 else
0111 {
0112 return Fun4AllReturnCodes::EVENT_OK;
0113 }
0114 }
0115
0116 int DecayFinder::End(PHCompositeNode* )
0117 {
0118 if (Verbosity() >= VERBOSITY_SOME)
0119 {
0120 printInfo();
0121 }
0122
0123 return 0;
0124 }
0125
0126
0127
0128
0129
0130
0131
0132
0133 int DecayFinder::parseDecayDescriptor()
0134 {
0135 bool ddCanBeParsed = true;
0136
0137 size_t daughterLocator;
0138
0139 std::string mother;
0140 std::string intermediate;
0141 std::string daughter;
0142
0143 int mother_charge = 0;
0144 std::vector<int> intermediates_charge;
0145 std::vector<int> daughters_charge;
0146
0147 std::string decayArrow = "->";
0148 std::string chargeIndicator = "^";
0149 std::string startIntermediate = "{";
0150 std::string endIntermediate = "}";
0151
0152
0153 std::string specialTracks[] = {"e", "mu", "pi", "K"};
0154
0155 std::string manipulateDecayDescriptor = m_decayDescriptor;
0156
0157
0158 size_t pos;
0159 while ((pos = manipulateDecayDescriptor.find(' ')) != std::string::npos)
0160 {
0161 manipulateDecayDescriptor.replace(pos, 1, "");
0162 }
0163
0164
0165 std::string checkForCC = manipulateDecayDescriptor.substr(0, 1) + manipulateDecayDescriptor.substr(manipulateDecayDescriptor.size() - 3, 3);
0166 std::for_each(checkForCC.begin(), checkForCC.end(), [](char& c)
0167 { c = ::toupper(c); });
0168
0169
0170 if (checkForCC == "[]CC")
0171 {
0172 manipulateDecayDescriptor = manipulateDecayDescriptor.substr(1, manipulateDecayDescriptor.size() - 4);
0173 m_getChargeConjugate = true;
0174 }
0175
0176
0177 size_t findMotherEndPoint = manipulateDecayDescriptor.find(decayArrow);
0178 mother = manipulateDecayDescriptor.substr(0, findMotherEndPoint);
0179 if (findParticle(mother))
0180 {
0181 m_mother_ID = abs(get_pdgcode(mother));
0182 }
0183 else
0184 {
0185 ddCanBeParsed = false;
0186 }
0187 manipulateDecayDescriptor.erase(0, findMotherEndPoint + decayArrow.length());
0188
0189
0190 while ((pos = manipulateDecayDescriptor.find(startIntermediate)) != std::string::npos)
0191 {
0192 size_t findIntermediateStartPoint = manipulateDecayDescriptor.find(startIntermediate, pos);
0193 size_t findIntermediateEndPoint = manipulateDecayDescriptor.find(endIntermediate, pos);
0194 std::string intermediateDecay = manipulateDecayDescriptor.substr(pos + 1, findIntermediateEndPoint - (pos + 1));
0195
0196 intermediate = intermediateDecay.substr(0, intermediateDecay.find(decayArrow));
0197 if (findParticle(intermediate))
0198 {
0199 m_intermediates_ID.push_back(abs(get_pdgcode(intermediate)));
0200 }
0201 else
0202 {
0203 ddCanBeParsed = false;
0204 }
0205
0206
0207 int nDaughters = 0;
0208 intermediateDecay.erase(0, intermediateDecay.find(decayArrow) + decayArrow.length());
0209 while ((daughterLocator = intermediateDecay.find(chargeIndicator)) != std::string::npos)
0210 {
0211 daughter = intermediateDecay.substr(0, daughterLocator);
0212 std::string daughterChargeString = intermediateDecay.substr(daughterLocator + 1, 1);
0213 if (std::find(std::begin(specialTracks), std::end(specialTracks), daughter) != std::end(specialTracks))
0214 {
0215 daughter += daughterChargeString;
0216 }
0217 if (findParticle(daughter))
0218 {
0219 m_daughters_ID.push_back(abs(get_pdgcode(daughter)));
0220 daughters_charge.push_back(get_charge(daughter));
0221 }
0222 else
0223 {
0224 ddCanBeParsed = false;
0225 }
0226 intermediateDecay.erase(0, daughterLocator + 2);
0227 ++nDaughters;
0228 }
0229 manipulateDecayDescriptor.erase(findIntermediateStartPoint, findIntermediateEndPoint + 1 - findIntermediateStartPoint);
0230 m_nTracksFromIntermediates.push_back(nDaughters);
0231 m_nTracksFromMother += 1;
0232 }
0233
0234
0235 while ((daughterLocator = manipulateDecayDescriptor.find(chargeIndicator)) != std::string::npos)
0236 {
0237 daughter = manipulateDecayDescriptor.substr(0, daughterLocator);
0238 std::string daughterChargeString = manipulateDecayDescriptor.substr(daughterLocator + 1, 1);
0239 if (std::find(std::begin(specialTracks), std::end(specialTracks), daughter) != std::end(specialTracks))
0240 {
0241 daughter += daughterChargeString;
0242 }
0243 if (findParticle(daughter))
0244 {
0245 m_daughters_ID.push_back(abs(get_pdgcode(daughter)));
0246 daughters_charge.push_back(get_charge(daughter));
0247 }
0248 else
0249 {
0250 ddCanBeParsed = false;
0251 }
0252 manipulateDecayDescriptor.erase(0, daughterLocator + 2);
0253 m_nTracksFromMother += 1;
0254 }
0255
0256 unsigned int trackStart = 0;
0257 unsigned int trackEnd = 0;
0258 for (unsigned int i = 0; i < m_intermediates_ID.size(); ++i)
0259 {
0260 trackStart = trackEnd;
0261 trackEnd = m_nTracksFromIntermediates[i] + trackStart;
0262
0263 int vtxCharge = 0;
0264
0265 for (unsigned int j = trackStart; j < trackEnd; ++j)
0266 {
0267 vtxCharge += daughters_charge[j];
0268 }
0269
0270 intermediates_charge.push_back(vtxCharge);
0271 }
0272
0273 for (unsigned int i = 0; i < m_daughters_ID.size(); ++i)
0274 {
0275 mother_charge += daughters_charge[i];
0276 }
0277
0278 m_mother_ID = mother_charge == 0 ? m_mother_ID : mother_charge * m_mother_ID;
0279 for (unsigned int i = 0; i < m_intermediates_ID.size(); ++i)
0280 {
0281 m_intermediates_ID[i] = intermediates_charge[i] == 0 ? m_intermediates_ID[i] : intermediates_charge[i] * m_intermediates_ID[i];
0282 m_motherDecayProducts.push_back(m_intermediates_ID[i]);
0283 }
0284 for (unsigned int i = 0; i < m_daughters_ID.size(); ++i)
0285 {
0286 m_daughters_ID[i] = daughters_charge[i] == 0 ? m_daughters_ID[i] : daughters_charge[i] * m_daughters_ID[i];
0287 if (i >= trackEnd)
0288 {
0289 m_motherDecayProducts.push_back(m_daughters_ID[i]);
0290 }
0291 }
0292
0293 if (ddCanBeParsed)
0294 {
0295 if (Verbosity() >= VERBOSITY_MORE)
0296 {
0297 std::cout << "Your decay descriptor can be parsed" << std::endl;
0298 }
0299 return 0;
0300 }
0301 else
0302 {
0303 if (Verbosity() >= VERBOSITY_SOME)
0304 {
0305 std::cout << "Your decay descriptor cannot be parsed, " << Name() << " will not be registered" << std::endl;
0306 }
0307 return Fun4AllReturnCodes::DONOTREGISTERSUBSYSTEM;
0308 }
0309 }
0310
0311
0312
0313
0314
0315
0316
0317 bool DecayFinder::findDecay(PHCompositeNode* topNode)
0318 {
0319 bool decayWasFound = false;
0320 bool reconstructableDecayWasFound = false;
0321 bool aTrackFailedPT = false;
0322 bool aTrackFailedETA = false;
0323 bool aMotherHasPhoton = false;
0324 bool aMotherHasPi0 = false;
0325 std::vector<int> correctMotherProducts;
0326 std::vector<int> positive_motherDecayProducts;
0327
0328 int n = sizeof(listOfResonantPIDs) / sizeof(listOfResonantPIDs[0]);
0329 for (int i : m_intermediates_ID)
0330 {
0331 n = deleteElement(listOfResonantPIDs, n, i);
0332 }
0333
0334 n = deleteElement(listOfResonantPIDs, n, m_mother_ID);
0335
0336 for (int m_motherDecayProduct : m_motherDecayProducts)
0337 {
0338 positive_motherDecayProducts.push_back(std::abs(m_motherDecayProduct));
0339 }
0340
0341 m_truthinfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0342 if (!m_truthinfo)
0343 {
0344 std::cout << "DecayFinder: Missing node G4TruthInfo" << std::endl;
0345 }
0346
0347 m_geneventmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
0348 if (!m_geneventmap)
0349 {
0350 std::cout << "DecayFinder: Missing node PHHepMCGenEventMap" << std::endl;
0351 }
0352
0353 if (!m_truthinfo && !m_geneventmap)
0354 {
0355 std::cout << "You have neither the PHHepMCGenEventMap or G4TruthInfo nodes" << std::endl;
0356 std::cout << "DecayFinder will crash, exiting now!" << std::endl;
0357 exit(1);
0358 }
0359
0360 if (m_truthinfo && !m_geneventmap)
0361 {
0362 if (Verbosity() >= VERBOSITY_SOME)
0363 {
0364 std::cout << "You have a valid G4TruthInfo but no PHHepMCGenEventMap container" << std::endl;
0365 std::cout << "You are probably using a particle gun" << std::endl;
0366 }
0367
0368 PHG4TruthInfoContainer::ConstRange range = m_truthinfo->GetParticleRange();
0369 for (PHG4TruthInfoContainer::ConstIterator iter = range.first; iter != range.second; ++iter)
0370 {
0371 PHG4Particle* g4particle = iter->second;
0372 int this_pid = m_getChargeConjugate ? abs(g4particle->get_pid()) : g4particle->get_pid();
0373 if (this_pid == m_mother_ID)
0374 {
0375 if (Verbosity() >= VERBOSITY_MAX)
0376 {
0377 std::cout << "parent->pdg_id(): " << g4particle->get_pid() << std::endl;
0378 }
0379
0380 bool breakOut = false;
0381 correctMotherProducts.clear();
0382 decayChain.clear();
0383 decayChain.push_back(std::make_pair(std::make_pair(g4particle->get_primary_id(), g4particle->get_barcode()), g4particle->get_pid()));
0384
0385 searchGeant4Record(g4particle->get_barcode(), g4particle->get_pid(), positive_motherDecayProducts, breakOut, aMotherHasPhoton, aMotherHasPi0, aTrackFailedPT, aTrackFailedETA, correctMotherProducts);
0386
0387 if (breakOut)
0388 {
0389 continue;
0390 }
0391
0392 decayWasFound = compareDecays(m_motherDecayProducts, correctMotherProducts);
0393
0394 if (decayWasFound)
0395 {
0396 m_counter += 1;
0397 if (aTrackFailedPT && !aTrackFailedETA)
0398 {
0399 m_nCandFail_pT += 1;
0400 }
0401 else if (!aTrackFailedPT && aTrackFailedETA)
0402 {
0403 m_nCandFail_eta += 1;
0404 }
0405 else if (aTrackFailedPT && aTrackFailedETA)
0406 {
0407 m_nCandFail_pT_and_eta += 1;
0408 }
0409 else
0410 {
0411 m_nCandReconstructable += 1;
0412 reconstructableDecayWasFound = true;
0413 if (m_save_dst)
0414 {
0415 fillDecayNode(topNode, decayChain);
0416 }
0417 if (aMotherHasPhoton && !aMotherHasPi0)
0418 {
0419 m_nCandHas_Photon += 1;
0420 }
0421 else if (!aMotherHasPhoton && aMotherHasPi0)
0422 {
0423 m_nCandHas_Pi0 += 1;
0424 }
0425 else if (aMotherHasPhoton && aMotherHasPi0)
0426 {
0427 m_nCandHas_Photon_and_Pi0 += 1;
0428 }
0429 else
0430 {
0431 m_nCandHas_noPhoton_and_noPi0 += 1;
0432 }
0433 }
0434 }
0435 }
0436 }
0437 return reconstructableDecayWasFound;
0438 }
0439
0440 for (PHHepMCGenEventMap::ConstReverseIter riter = m_geneventmap->rbegin();
0441 riter != m_geneventmap->rend(); ++riter)
0442 {
0443 m_genevt = riter->second;
0444 if (!m_genevt)
0445 {
0446 std::cout << "DecayFinder: Missing node PHHepMCGenEvent" << std::endl;
0447 return false;
0448 }
0449
0450 HepMC::GenEvent* theEvent = m_genevt->getEvent();
0451
0452 for (HepMC::GenEvent::particle_const_iterator p = theEvent->particles_begin(); p != theEvent->particles_end(); ++p)
0453 {
0454 int this_pid = m_getChargeConjugate ? abs((*p)->pdg_id()) : (*p)->pdg_id();
0455 if (this_pid == m_mother_ID)
0456 {
0457 if (Verbosity() >= VERBOSITY_MAX)
0458 {
0459 std::cout << "parent->pdg_id(): " << (*p)->pdg_id() << std::endl;
0460 }
0461
0462 bool breakOut = false;
0463 correctMotherProducts.clear();
0464 decayChain.clear();
0465 decayChain.push_back(std::make_pair(std::make_pair(m_genevt->get_embedding_id(), (*p)->barcode()), (*p)->pdg_id()));
0466
0467
0468 if (!(*p)->end_vertex())
0469 {
0470 searchGeant4Record((*p)->barcode(), (*p)->pdg_id(), positive_motherDecayProducts, breakOut, aMotherHasPhoton, aMotherHasPi0, aTrackFailedPT, aTrackFailedETA, correctMotherProducts);
0471 }
0472 else
0473 {
0474 searchHepMCRecord((*p), positive_motherDecayProducts, breakOut, aMotherHasPhoton, aMotherHasPi0, aTrackFailedPT, aTrackFailedETA, correctMotherProducts);
0475 }
0476
0477 if (breakOut)
0478 {
0479 continue;
0480 }
0481
0482 decayWasFound = compareDecays(m_motherDecayProducts, correctMotherProducts);
0483
0484 if (decayWasFound)
0485 {
0486 m_counter += 1;
0487 if (aTrackFailedPT && !aTrackFailedETA)
0488 {
0489 m_nCandFail_pT += 1;
0490 }
0491 else if (!aTrackFailedPT && aTrackFailedETA)
0492 {
0493 m_nCandFail_eta += 1;
0494 }
0495 else if (aTrackFailedPT && aTrackFailedETA)
0496 {
0497 m_nCandFail_pT_and_eta += 1;
0498 }
0499 else
0500 {
0501 m_nCandReconstructable += 1;
0502 reconstructableDecayWasFound = true;
0503 if (m_save_dst)
0504 {
0505 fillDecayNode(topNode, decayChain);
0506 }
0507 if (aMotherHasPhoton && !aMotherHasPi0)
0508 {
0509 m_nCandHas_Photon += 1;
0510 }
0511 else if (!aMotherHasPhoton && aMotherHasPi0)
0512 {
0513 m_nCandHas_Pi0 += 1;
0514 }
0515 else if (aMotherHasPhoton && aMotherHasPi0)
0516 {
0517 m_nCandHas_Photon_and_Pi0 += 1;
0518 }
0519 else
0520 {
0521 m_nCandHas_noPhoton_and_noPi0 += 1;
0522 }
0523 }
0524 }
0525 }
0526 }
0527 }
0528
0529 return reconstructableDecayWasFound;
0530 }
0531
0532
0533
0534
0535
0536 void DecayFinder::searchHepMCRecord(HepMC::GenParticle* particle, std::vector<int> decayProducts,
0537 bool& breakLoop, bool& hasPhoton, bool& hasPi0, bool& failedPT, bool& failedETA,
0538 std::vector<int>& actualDecayProducts)
0539 {
0540 for (HepMC::GenVertex::particle_iterator children = particle->end_vertex()->particles_begin(HepMC::children);
0541 children != particle->end_vertex()->particles_end(HepMC::children); ++children)
0542 {
0543 if (Verbosity() >= VERBOSITY_MAX)
0544 {
0545 std::cout << "--children->pdg_id(): " << (*children)->pdg_id() << std::endl;
0546 }
0547
0548 if ((*children)->pdg_id() == 22)
0549 {
0550 hasPhoton = true;
0551 }
0552 if ((*children)->pdg_id() == 111)
0553 {
0554 hasPi0 = true;
0555 }
0556 if ((!m_allowPhotons && (*children)->pdg_id() == 22) || (!m_allowPi0 && (*children)->pdg_id() == 111))
0557 {
0558 breakLoop = true;
0559 break;
0560 }
0561
0562
0563 if (std::find(decayProducts.begin(), decayProducts.end(),
0564 std::abs((*children)->pdg_id())) != decayProducts.end())
0565 {
0566 if (Verbosity() >= VERBOSITY_MAX)
0567 {
0568 std::cout << "This is a child you were looking for" << std::endl;
0569 }
0570
0571 std::vector<int> positive_intermediates_ID;
0572 for (int i : m_intermediates_ID)
0573 {
0574 positive_intermediates_ID.push_back(abs(i));
0575 }
0576
0577 if (!(*children)->end_vertex() && std::find(positive_intermediates_ID.begin(), positive_intermediates_ID.end(),
0578 abs((*children)->pdg_id())) != positive_intermediates_ID.end())
0579 {
0580 std::vector<int> requiredIntermediateDecayProducts;
0581 std::vector<int> positive_requiredIntermediateDecayProducts;
0582 std::vector<int> actualIntermediateDecayProducts;
0583
0584 auto it = std::find(positive_intermediates_ID.begin(), positive_intermediates_ID.end(), abs((*children)->pdg_id()));
0585 int index = std::distance(positive_intermediates_ID.begin(), it);
0586
0587 unsigned int trackStart = 0, trackStop = 0;
0588 if (index == 0)
0589 {
0590 trackStop = m_nTracksFromIntermediates[0];
0591 }
0592 else
0593 {
0594 for (int i = 0; i < index; ++i)
0595 {
0596 trackStart += m_nTracksFromIntermediates[i];
0597 }
0598 trackStop = trackStart + m_nTracksFromIntermediates[index];
0599 }
0600 for (unsigned int i = trackStart; i < trackStop; ++i)
0601 {
0602 requiredIntermediateDecayProducts.push_back(m_daughters_ID[i]);
0603 positive_requiredIntermediateDecayProducts.push_back(abs(m_daughters_ID[i]));
0604 }
0605
0606 searchGeant4Record((*children)->barcode(), (*children)->pdg_id(), positive_requiredIntermediateDecayProducts,
0607 breakLoop, hasPhoton, hasPi0, failedPT, failedETA, actualIntermediateDecayProducts);
0608
0609 bool needThisParticle = compareDecays(requiredIntermediateDecayProducts, actualIntermediateDecayProducts);
0610 if (needThisParticle)
0611 {
0612 actualDecayProducts.push_back((*children)->pdg_id());
0613 decayChain.push_back(std::make_pair(std::make_pair(m_genevt->get_embedding_id(), (*children)->barcode()), (*children)->pdg_id()));
0614 }
0615 else
0616 {
0617 decayChain.erase(decayChain.end() - (m_intermediate_product_counter + 1), decayChain.end());
0618 decayChain.pop_back();
0619 }
0620
0621 m_intermediate_product_counter = 0;
0622 }
0623 else
0624 {
0625 bool needThisParticle = checkIfCorrectHepMCParticle((*children), failedPT, failedETA);
0626 if (needThisParticle)
0627 {
0628 actualDecayProducts.push_back((*children)->pdg_id());
0629 decayChain.push_back(std::make_pair(std::make_pair(m_genevt->get_embedding_id(), (*children)->barcode()), (*children)->pdg_id()));
0630 }
0631 }
0632 }
0633 else if (std::find(std::begin(listOfResonantPIDs), std::end(listOfResonantPIDs),
0634 std::abs((*children)->pdg_id())) != std::end(listOfResonantPIDs))
0635 {
0636 if (Verbosity() >= VERBOSITY_MAX)
0637 {
0638 std::cout << "This is a resonance to investigate further" << std::endl;
0639 }
0640 if (!(*children)->end_vertex())
0641 {
0642 searchGeant4Record((*children)->barcode(), (*children)->pdg_id(), decayProducts,
0643 breakLoop, hasPhoton, hasPi0, failedPT, failedETA, actualDecayProducts);
0644 }
0645 else
0646 {
0647 for (HepMC::GenVertex::particle_iterator grandchildren = (*children)->end_vertex()->particles_begin(HepMC::children);
0648 grandchildren != (*children)->end_vertex()->particles_end(HepMC::children); ++grandchildren)
0649 {
0650 bool needThisParticle = checkIfCorrectHepMCParticle((*grandchildren), failedPT, failedETA);
0651 if (needThisParticle)
0652 {
0653 actualDecayProducts.push_back((*grandchildren)->pdg_id());
0654 decayChain.push_back(std::make_pair(std::make_pair(m_genevt->get_embedding_id(), (*grandchildren)->barcode()), (*grandchildren)->pdg_id()));
0655 }
0656 }
0657 }
0658 }
0659 else
0660 {
0661 breakLoop = true;
0662 break;
0663 }
0664 }
0665 }
0666
0667
0668
0669
0670
0671 void DecayFinder::searchGeant4Record(int barcode, int pid, std::vector<int> decayProducts, bool& breakLoop, bool& hasPhoton, bool& hasPi0, bool& failedPT, bool& failedETA, std::vector<int>& actualDecayProducts)
0672 {
0673 PHG4TruthInfoContainer::ConstRange range = m_truthinfo->GetParticleRange();
0674 for (PHG4TruthInfoContainer::ConstIterator iter = range.first; iter != range.second; ++iter)
0675 {
0676 if (decayChain.size() == 100)
0677 {
0678 breakLoop = true;
0679 break;
0680 }
0681
0682 PHG4Particle* g4particle = iter->second;
0683 PHG4Particle* mother = nullptr;
0684 if (g4particle->get_parent_id() != 0)
0685 {
0686 mother = m_truthinfo->GetParticle(g4particle->get_parent_id());
0687 }
0688 else
0689 {
0690 continue;
0691 }
0692 if (mother->get_barcode() == barcode && abs(mother->get_pid()) == abs(pid))
0693 {
0694 int particleID = g4particle->get_pid();
0695 if (Verbosity() >= VERBOSITY_MAX)
0696 {
0697 std::cout << "--children->pdg_id(): " << particleID << std::endl;
0698 }
0699 if (particleID == 22)
0700 {
0701 hasPhoton = true;
0702 }
0703 if (particleID == 111)
0704 {
0705 hasPi0 = true;
0706 }
0707
0708 if ((!m_allowPhotons && particleID == 22) || (!m_allowPi0 && particleID == 111))
0709 {
0710 breakLoop = true;
0711 break;
0712 }
0713
0714 if (std::find(decayProducts.begin(), decayProducts.end(),
0715 std::abs(particleID)) != decayProducts.end())
0716 {
0717 bool needThisParticle = checkIfCorrectGeant4Particle(g4particle, hasPhoton, hasPi0, failedPT, failedETA);
0718 if (needThisParticle)
0719 {
0720 if (Verbosity() >= VERBOSITY_MAX)
0721 {
0722 std::cout << "This is a child you were looking for" << std::endl;
0723 }
0724 actualDecayProducts.push_back(particleID);
0725 int embedding_id = m_geneventmap ? m_genevt->get_embedding_id() : g4particle->get_primary_id();
0726 decayChain.push_back(std::make_pair(std::make_pair(embedding_id, g4particle->get_barcode()), particleID));
0727 }
0728 }
0729 else if ((m_allowPhotons && particleID == 22) || (m_allowPi0 && particleID == 111))
0730 {
0731 continue;
0732 }
0733 else if (std::find(std::begin(listOfResonantPIDs), std::end(listOfResonantPIDs),
0734 std::abs(particleID)) != std::end(listOfResonantPIDs))
0735 {
0736 if (Verbosity() >= VERBOSITY_MAX)
0737 {
0738 std::cout << "This is a resonance to investigate further" << std::endl;
0739 }
0740 searchGeant4Record(g4particle->get_barcode(), g4particle->get_pid(), decayProducts,
0741 breakLoop, hasPhoton, hasPi0, failedPT, failedETA, actualDecayProducts);
0742 }
0743 else
0744 {
0745 breakLoop = true;
0746 break;
0747 }
0748 }
0749 }
0750 }
0751
0752
0753
0754
0755 bool DecayFinder::checkIfCorrectHepMCParticle(HepMC::GenParticle* particle, bool& trackFailedPT, bool& trackFailedETA)
0756 {
0757 bool acceptParticle = false;
0758
0759 std::vector<int> positive_intermediates_ID;
0760 for (int i : m_intermediates_ID)
0761 {
0762 positive_intermediates_ID.push_back(abs(i));
0763 }
0764
0765
0766 if (std::find(positive_intermediates_ID.begin(), positive_intermediates_ID.end(),
0767 abs(particle->pdg_id())) != positive_intermediates_ID.end())
0768 {
0769 std::vector<int> requiredIntermediateDecayProducts;
0770 std::vector<int> actualIntermediateDecayProducts;
0771
0772 auto it = std::find(positive_intermediates_ID.begin(), positive_intermediates_ID.end(), abs(particle->pdg_id()));
0773 int index = std::distance(positive_intermediates_ID.begin(), it);
0774
0775 unsigned int trackStart = 0, trackStop = 0;
0776 if (index == 0)
0777 {
0778 trackStop = m_nTracksFromIntermediates[0];
0779 }
0780 else
0781 {
0782 for (int i = 0; i < index; ++i)
0783 {
0784 trackStart += m_nTracksFromIntermediates[i];
0785 }
0786 trackStop = trackStart + m_nTracksFromIntermediates[index];
0787 }
0788
0789 for (unsigned int i = trackStart; i < trackStop; ++i)
0790 {
0791 requiredIntermediateDecayProducts.push_back(m_daughters_ID[i]);
0792 }
0793
0794 for (HepMC::GenVertex::particle_iterator grandchildren = particle->end_vertex()->particles_begin(HepMC::children);
0795 grandchildren != particle->end_vertex()->particles_end(HepMC::children); ++grandchildren)
0796 {
0797 if (Verbosity() >= VERBOSITY_MAX)
0798 {
0799 std::cout << "----grandchildren->pdg_id(): " << (*grandchildren)->pdg_id() << std::endl;
0800 }
0801
0802 if (std::find(std::begin(listOfResonantPIDs), std::end(listOfResonantPIDs),
0803 std::abs((*grandchildren)->pdg_id())) != std::end(listOfResonantPIDs))
0804 {
0805 for (HepMC::GenVertex::particle_iterator greatgrandchildren = (*grandchildren)->end_vertex()->particles_begin(HepMC::children);
0806 greatgrandchildren != (*grandchildren)->end_vertex()->particles_end(HepMC::children); ++greatgrandchildren)
0807 {
0808 if (Verbosity() >= VERBOSITY_MAX)
0809 {
0810 std::cout << "--------greatgrandchildren->pdg_id(): " << (*greatgrandchildren)->pdg_id() << std::endl;
0811 }
0812
0813 if ((m_allowPhotons && (*greatgrandchildren)->pdg_id() == 22) || (m_allowPi0 && (*greatgrandchildren)->pdg_id() == 111))
0814 {
0815 continue;
0816 }
0817 else if ((!m_allowPhotons && (*greatgrandchildren)->pdg_id() == 22) || (!m_allowPi0 && (*greatgrandchildren)->pdg_id() == 111))
0818 {
0819 break;
0820 }
0821 else
0822 {
0823 actualIntermediateDecayProducts.push_back((*greatgrandchildren)->pdg_id());
0824 decayChain.push_back(std::make_pair(std::make_pair(m_genevt->get_embedding_id(), (*greatgrandchildren)->barcode()), (*greatgrandchildren)->pdg_id()));
0825 ++m_intermediate_product_counter;
0826
0827 HepMC::FourVector myFourVector = (*greatgrandchildren)->momentum();
0828
0829 HepMC::GenVertex* thisVtx = (*greatgrandchildren)->production_vertex();
0830 double vtxPos[3] = {thisVtx->point3d().x(), thisVtx->point3d().y(), thisVtx->point3d().z()};
0831 if (m_recalcualteEtaRange)
0832 {
0833 recalculateEta(myFourVector.py(), vtxPos);
0834 }
0835
0836 if (myFourVector.perp() < m_pt_req)
0837 {
0838 trackFailedPT = true;
0839 }
0840 if (!isInRange(m_eta_low_req, myFourVector.eta(), m_eta_high_req))
0841 {
0842 trackFailedETA = true;
0843 }
0844
0845 if (Verbosity() >= VERBOSITY_MAX)
0846 {
0847 std::cout << "pT = " << myFourVector.perp() << ", eta = " << myFourVector.eta() << std::endl;
0848 std::cout << "The track " << passOrFail(myFourVector.perp() >= m_pt_req) << " the pT requirement, ";
0849 std::cout << "the track " << passOrFail(isInRange(m_eta_low_req, myFourVector.eta(), m_eta_high_req)) << " the eta requirement." << std::endl;
0850 }
0851 }
0852 }
0853 }
0854 else if ((m_allowPhotons && (*grandchildren)->pdg_id() == 22) || (m_allowPi0 && (*grandchildren)->pdg_id() == 111))
0855 {
0856 continue;
0857 }
0858 else if ((!m_allowPhotons && (*grandchildren)->pdg_id() == 22) || (!m_allowPi0 && (*grandchildren)->pdg_id() == 111))
0859 {
0860 break;
0861 }
0862 else
0863 {
0864 actualIntermediateDecayProducts.push_back((*grandchildren)->pdg_id());
0865 decayChain.push_back(std::make_pair(std::make_pair(m_genevt->get_embedding_id(), (*grandchildren)->barcode()), (*grandchildren)->pdg_id()));
0866 ++m_intermediate_product_counter;
0867
0868 HepMC::FourVector myFourVector = (*grandchildren)->momentum();
0869
0870 HepMC::GenVertex* thisVtx = (*grandchildren)->production_vertex();
0871 double vtxPos[3] = {thisVtx->point3d().x(), thisVtx->point3d().y(), thisVtx->point3d().z()};
0872 if (m_recalcualteEtaRange)
0873 {
0874 recalculateEta(myFourVector.py(), vtxPos);
0875 }
0876
0877 if (myFourVector.perp() < m_pt_req)
0878 {
0879 trackFailedPT = true;
0880 }
0881 if (!isInRange(m_eta_low_req, myFourVector.eta(), m_eta_high_req))
0882 {
0883 trackFailedETA = true;
0884 }
0885
0886 if (Verbosity() >= VERBOSITY_MAX)
0887 {
0888 std::cout << "pT = " << myFourVector.perp() << ", eta = " << myFourVector.eta() << std::endl;
0889 std::cout << "The track " << passOrFail(myFourVector.perp() >= m_pt_req) << " the pT requirement, ";
0890 std::cout << "the track " << passOrFail(isInRange(m_eta_low_req, myFourVector.eta(), m_eta_high_req)) << " the eta requirement." << std::endl;
0891 }
0892 }
0893 }
0894
0895 acceptParticle = compareDecays(requiredIntermediateDecayProducts, actualIntermediateDecayProducts);
0896 }
0897 else if ((particle->pdg_id() == 22) || (particle->pdg_id() == 111))
0898 {
0899 return false;
0900 }
0901 else
0902 {
0903 if (Verbosity() >= VERBOSITY_MAX)
0904 {
0905 std::cout << "This is a final state track" << std::endl;
0906 }
0907 const HepMC::FourVector& myFourVector = particle->momentum();
0908
0909 HepMC::GenVertex* thisVtx = particle->production_vertex();
0910 double vtxPos[3] = {thisVtx->point3d().x(), thisVtx->point3d().y(), thisVtx->point3d().z()};
0911 if (m_recalcualteEtaRange)
0912 {
0913 recalculateEta(myFourVector.py(), vtxPos);
0914 }
0915
0916 if (myFourVector.perp() < m_pt_req)
0917 {
0918 trackFailedPT = true;
0919 }
0920 if (!isInRange(m_eta_low_req, myFourVector.eta(), m_eta_high_req))
0921 {
0922 trackFailedETA = true;
0923 }
0924
0925 if (Verbosity() >= VERBOSITY_MAX)
0926 {
0927 std::cout << "pT = " << myFourVector.perp() << ", eta = " << myFourVector.eta() << std::endl;
0928 std::cout << "The track " << passOrFail(myFourVector.perp() >= m_pt_req) << " the pT requirement, ";
0929 std::cout << "the track " << passOrFail(isInRange(m_eta_low_req, myFourVector.eta(), m_eta_high_req)) << " the eta requirement." << std::endl;
0930 }
0931 acceptParticle = true;
0932 }
0933
0934 return acceptParticle;
0935 }
0936
0937
0938
0939
0940 bool DecayFinder::checkIfCorrectGeant4Particle(PHG4Particle* particle, bool& hasPhoton, bool& hasPi0, bool& trackFailedPT, bool& trackFailedETA)
0941 {
0942 bool acceptParticle = false;
0943
0944 std::vector<int> positive_intermediates_ID;
0945 for (int i : m_intermediates_ID)
0946 {
0947 positive_intermediates_ID.push_back(abs(i));
0948 }
0949
0950 if (std::find(positive_intermediates_ID.begin(), positive_intermediates_ID.end(),
0951 abs(particle->get_pid())) != positive_intermediates_ID.end())
0952 {
0953 std::vector<int> requiredIntermediateDecayProducts;
0954 std::vector<int> actualIntermediateDecayProducts;
0955
0956 auto it = std::find(positive_intermediates_ID.begin(), positive_intermediates_ID.end(), abs(particle->get_pid()));
0957 int index = std::distance(positive_intermediates_ID.begin(), it);
0958
0959 unsigned int trackStart = 0, trackStop = 0;
0960 if (index == 0)
0961 {
0962 trackStop = m_nTracksFromIntermediates[0];
0963 }
0964 else
0965 {
0966 for (int i = 0; i < index; ++i)
0967 {
0968 trackStart += m_nTracksFromIntermediates[i];
0969 }
0970 trackStop = trackStart + m_nTracksFromIntermediates[index];
0971 }
0972
0973 std::vector<int> positive_intermediateDecayProducts;
0974 for (unsigned int i = trackStart; i < trackStop; ++i)
0975 {
0976 positive_intermediateDecayProducts.push_back(abs(m_daughters_ID[i]));
0977 requiredIntermediateDecayProducts.push_back(m_daughters_ID[i]);
0978 }
0979
0980 bool fakeBreak = false;
0981 searchGeant4Record(particle->get_barcode(), particle->get_pid(), positive_intermediateDecayProducts, fakeBreak,
0982 hasPhoton, hasPi0, trackFailedPT, trackFailedETA, actualIntermediateDecayProducts);
0983
0984 acceptParticle = compareDecays(requiredIntermediateDecayProducts, actualIntermediateDecayProducts);
0985 }
0986 else if ((particle->get_pid() == 22) || (particle->get_pid() == 111))
0987 {
0988 return false;
0989 }
0990 else
0991 {
0992 if (Verbosity() >= VERBOSITY_MAX)
0993 {
0994 std::cout << "This is a final state track" << std::endl;
0995 }
0996 double px = particle->get_px();
0997 double py = particle->get_py();
0998 double pz = particle->get_pz();
0999
1000 double pt = std::sqrt(std::pow(px, 2) + std::pow(py, 2));
1001 double p = std::sqrt(std::pow(px, 2) + std::pow(py, 2) + std::pow(pz, 2));
1002 double eta = 0.5 * std::log((p + pz) / (p - pz));
1003
1004 PHG4VtxPoint* thisVtx = m_truthinfo->GetVtx(particle->get_vtx_id());
1005 double vtxPos[3] = {thisVtx->get_x(), thisVtx->get_y(), thisVtx->get_z()};
1006 if (m_recalcualteEtaRange)
1007 {
1008 recalculateEta(py, vtxPos);
1009 }
1010
1011 if (pt < m_pt_req)
1012 {
1013 trackFailedPT = true;
1014 }
1015 if (!isInRange(m_eta_low_req, eta, m_eta_high_req))
1016 {
1017 trackFailedETA = true;
1018 }
1019
1020 if (Verbosity() >= VERBOSITY_MAX)
1021 {
1022 std::cout << "pT = " << pt << ", eta = " << eta << std::endl;
1023 std::cout << "The track " << passOrFail(pt >= m_pt_req) << " the pT requirement, ";
1024 std::cout << "the track " << passOrFail(isInRange(m_eta_low_req, eta, m_eta_high_req)) << " the eta requirement." << std::endl;
1025 }
1026
1027 acceptParticle = true;
1028 }
1029
1030 return acceptParticle;
1031 }
1032
1033
1034
1035
1036 bool DecayFinder::compareDecays(std::vector<int> required, std::vector<int> actual)
1037 {
1038 bool accept = false;
1039 multiplyVectorByScalarAndSort(required, +1);
1040 multiplyVectorByScalarAndSort(actual, +1);
1041
1042 if (Verbosity() >= VERBOSITY_MAX)
1043 {
1044 std::cout << "Printing required decay products: ";
1045 for (int i : required)
1046 {
1047 std::cout << i << ", ";
1048 }
1049 std::cout << std::endl;
1050 std::cout << "Printing actual decay products: ";
1051 for (int i : actual)
1052 {
1053 std::cout << i << ", ";
1054 }
1055 std::cout << std::endl;
1056 if (required == actual)
1057 {
1058 std::cout << "*\n* These vectors match\n*\n"
1059 << std::endl;
1060 }
1061 else
1062 {
1063 std::cout << "*\n* These vectors DONT match\n*\n"
1064 << std::endl;
1065 }
1066 }
1067
1068 if (required == actual)
1069 {
1070 accept = true;
1071 }
1072
1073 if (m_getChargeConjugate && !accept)
1074 {
1075 multiplyVectorByScalarAndSort(required, -1);
1076 if (required == actual)
1077 {
1078 accept = true;
1079 }
1080
1081 if (Verbosity() >= VERBOSITY_MAX)
1082 {
1083 std::cout << "Checking CC state" << std::endl;
1084 if (required == actual)
1085 {
1086 std::cout << "*\n* These vectors match\n*\n"
1087 << std::endl;
1088 }
1089 else
1090 {
1091 std::cout << "*\n* These vectors DONT match\n*\n"
1092 << std::endl;
1093 }
1094 }
1095 }
1096 return accept;
1097 }
1098
1099
1100
1101
1102 int DecayFinder::deleteElement(int arr[], int n, int x)
1103 {
1104
1105
1106 int i;
1107 for (i = 0; i < n; i++)
1108 {
1109 if (arr[i] == x)
1110 {
1111 break;
1112 }
1113 }
1114
1115
1116 if (i < n)
1117 {
1118
1119
1120 n = n - 1;
1121 for (int j = i; j < n; j++)
1122 {
1123 arr[j] = arr[j + 1];
1124 }
1125 }
1126
1127 return n;
1128 }
1129
1130 bool DecayFinder::findParticle(const std::string& particle)
1131 {
1132 bool particleFound = true;
1133 if (!TDatabasePDG::Instance()->GetParticle(particle.c_str()))
1134 {
1135 if (Verbosity() >= VERBOSITY_SOME)
1136 {
1137 std::cout << "The particle, " << particle << " is not in the TDatabasePDG particle list" << std::endl;
1138 }
1139 particleFound = false;
1140 }
1141
1142 return particleFound;
1143 }
1144
1145 void DecayFinder::multiplyVectorByScalarAndSort(std::vector<int>& v, int k)
1146 {
1147
1148 std::transform(v.begin(), v.end(), v.begin(), [k](const int& c)
1149 {
1150 int particlesWithNoCC[] = {111, 113, 115, 130, 220, 221, 223, 225, 310, 330, 331, 333, 335, 440, 441, 443,
1151 445, 551, 553, 555, 10111, 10113, 10221, 10223, 10331, 10333, 10441, 10443,
1152 10551, 10553, 20113, 20223, 20333, 20443, 20553, 100443, 100553};
1153 if (std::find(std::begin(particlesWithNoCC), std::end(particlesWithNoCC), c) != std::end(particlesWithNoCC))
1154 {
1155 return c;
1156 }
1157 else
1158 {
1159 return c * k;
1160 }
1161 });
1162 std::sort(v.begin(), v.end());
1163 }
1164
1165 int DecayFinder::get_pdgcode(const std::string& name)
1166 {
1167 if (findParticle(name))
1168 {
1169 return TDatabasePDG::Instance()->GetParticle(name.c_str())->PdgCode();
1170 }
1171 else
1172 {
1173 return 0;
1174 }
1175 }
1176
1177 int DecayFinder::get_charge(const std::string& name)
1178 {
1179 if (findParticle(name))
1180 {
1181 return TDatabasePDG::Instance()->GetParticle(name.c_str())->Charge() / 3;
1182 }
1183 else
1184 {
1185 return -99;
1186 }
1187 }
1188
1189 bool DecayFinder::isInRange(float min, float value, float max)
1190 {
1191 return min <= value && value <= max;
1192 }
1193
1194 void DecayFinder::calculateEffectiveTPCradius(double vertex[3], double& effective_top_r, double& effective_bottom_r)
1195 {
1196 double r_sq = std::pow(m_tpc_r, 2);
1197 double x_sq = std::pow(vertex[0], 2);
1198
1199 effective_top_r = std::sqrt(r_sq - x_sq) - vertex[1];
1200 effective_bottom_r = std::sqrt(r_sq - x_sq) + vertex[1];
1201 }
1202
1203 void DecayFinder::recalculateEta(double py, double vertex[3])
1204 {
1205 calculateEffectiveTPCradius(vertex, m_effective_top_tpc_r, m_effective_bottom_tpc_r);
1206
1207 double z_north = m_tpc_z - vertex[2];
1208 double z_south = m_tpc_z + vertex[2];
1209
1210 double top_north_theta = atan(m_effective_top_tpc_r / z_north);
1211 double bottom_north_theta = atan(m_effective_bottom_tpc_r / z_north);
1212 double top_south_theta = atan(m_effective_top_tpc_r / z_south);
1213 double bottom_south_theta = atan(m_effective_bottom_tpc_r / z_south);
1214
1215 double top_north_eta = -1 * std::log(tan(top_north_theta / 2));
1216 double bottom_north_eta = -1 * std::log(tan(bottom_north_theta / 2));
1217 double top_south_eta = std::log(tan(top_south_theta / 2));
1218 double bottom_south_eta = std::log(tan(bottom_south_theta / 2));
1219
1220 m_eta_high_req = py < 0 ? bottom_north_eta : top_north_eta;
1221 m_eta_low_req = py < 0 ? bottom_south_eta : top_south_eta;
1222 }
1223
1224 int DecayFinder::createDecayNode(PHCompositeNode* topNode)
1225 {
1226 PHNodeIterator iter(topNode);
1227
1228 PHCompositeNode* lowerNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
1229 if (!lowerNode)
1230 {
1231 lowerNode = new PHCompositeNode("DST");
1232 topNode->addNode(lowerNode);
1233 std::cout << "DST node added" << std::endl;
1234 }
1235
1236 std::string baseName;
1237
1238 if (m_container_name.empty())
1239 {
1240
1241 baseName = Name();
1242 }
1243 else
1244 {
1245 baseName = m_container_name;
1246 }
1247
1248
1249 size_t pos;
1250 std::string undrscr = "_";
1251 std::string nothing = "";
1252 std::map<std::string, std::string> forbiddenStrings;
1253 forbiddenStrings["/"] = undrscr;
1254 forbiddenStrings["("] = undrscr;
1255 forbiddenStrings[")"] = nothing;
1256 forbiddenStrings["+"] = "plus";
1257 forbiddenStrings["-"] = "minus";
1258 forbiddenStrings["*"] = "star";
1259 forbiddenStrings[" "] = nothing;
1260 for (auto const& [badString, goodString] : forbiddenStrings)
1261 {
1262 while ((pos = baseName.find(badString)) != std::string::npos)
1263 {
1264 baseName.replace(pos, 1, goodString);
1265 }
1266 }
1267
1268 m_nodeName = baseName + "_DecayMap";
1269
1270 m_decayMap = new DecayFinderContainer_v1();
1271 PHIODataNode<PHObject>* decayNode = new PHIODataNode<PHObject>(m_decayMap, m_nodeName.c_str(), "PHObject");
1272 lowerNode->addNode(decayNode);
1273 std::cout << m_nodeName << " node added" << std::endl;
1274
1275 return Fun4AllReturnCodes::EVENT_OK;
1276 }
1277
1278 void DecayFinder::fillDecayNode(PHCompositeNode* topNode, Decay& decay)
1279 {
1280 m_decayMap = findNode::getClass<DecayFinderContainer_v1>(topNode, m_nodeName.c_str());
1281 m_decayMap->insert(decay);
1282 }
1283
1284 void DecayFinder::printInfo()
1285 {
1286 std::cout << "\n---------------DecayFinder information---------------" << std::endl;
1287 std::cout << "Module name: " << Name() << std::endl;
1288 std::cout << "Decay descriptor: " << m_decayDescriptor << std::endl;
1289 std::cout << "Number of generated decays: " << m_counter << std::endl;
1290 std::cout << " Number of decays that failed pT requirement: " << m_nCandFail_pT << std::endl;
1291 std::cout << " Number of decays that failed eta requirement: " << m_nCandFail_eta << std::endl;
1292 std::cout << " Number of decays that failed pT and eta requirements: " << m_nCandFail_pT_and_eta << std::endl;
1293 std::cout << "Number of decays that could be reconstructed: " << m_nCandReconstructable << std::endl;
1294 std::cout << " Number of reconstructable mothers with associated photon: " << m_nCandHas_Photon << std::endl;
1295 std::cout << " Number of reconstructable mothers with associated Pi0: " << m_nCandHas_Pi0 << std::endl;
1296 std::cout << " Number of reconstructable mothers with associated photon and Pi0: " << m_nCandHas_Photon_and_Pi0 << std::endl;
1297 std::cout << " Number of reconstructable mothers with no associated photons or Pi0s: " << m_nCandHas_noPhoton_and_noPi0 << std::endl;
1298 std::cout << "-----------------------------------------------------\n";
1299 std::cout << std::endl;
1300 }
1301
1302 void DecayFinder::printNode(PHCompositeNode* topNode)
1303 {
1304 std::cout << "----------------";
1305 std::cout << " DecayFinderNode: " << m_nodeName << " information ";
1306 std::cout << "----------------" << std::endl;
1307 DecayFinderContainer_v1* map = findNode::getClass<DecayFinderContainer_v1>(topNode, m_nodeName.c_str());
1308 for (auto& iter : *map)
1309 {
1310 Decay decay = iter.second;
1311 for (auto& i : decay)
1312 {
1313 std::cout << "Particle embedding ID: " << i.first.first << ", particle barcode: " << i.first.second << ", particle PDG ID: " << i.second << std::endl;
1314 }
1315 }
1316 std::cout << "--------------------------------------------------------------------------------------------------" << std::endl;
1317 }
1318
1319 std::string DecayFinder::passOrFail(bool condition)
1320 {
1321 return condition ? "passed" : "failed";
1322 }