LDMX Software
BertiniAtLeastNProductsModel.cxx
1
2#include "SimCore/PhotoNuclearModels/BertiniAtLeastNProductsModel.h"
3namespace simcore {
4
5bool BertiniAtLeastNProductsProcess::acceptEvent() const {
6 int secondaries{theParticleChange.GetNumberOfSecondaries()};
7 int matching_count{0};
8 for (int i{0}; i < secondaries; ++i) {
9 const auto secondary{theParticleChange.GetSecondary(i)->GetParticle()};
10 const auto pdg_code{secondary->GetDefinition()->GetPDGEncoding()};
11 const auto energy{secondary->GetKineticEnergy()};
12 if (std::find(std::begin(pdg_ids_), std::end(pdg_ids_), pdg_code) !=
13 std::end(pdg_ids_)) {
14 if (energy > threshold_) {
15 ++matching_count;
16 }
17 }
18 if (matching_count >= min_products_) {
19 return true;
20 }
21 }
22 return false;
23}
24
26 G4ProcessManager* processManager) {
27 auto photo_nuclear_process{
28 new G4HadronInelasticProcess("photonNuclear", G4Gamma::Definition())};
29 auto model{new BertiniAtLeastNProductsProcess{threshold_, zmin_, emin_,
30 pdg_ids_, min_products_}};
31 model->SetMaxEnergy(15 * CLHEP::GeV);
32 addPNCrossSectionData(photo_nuclear_process);
33 photo_nuclear_process->RegisterMe(model);
34 processManager->AddDiscreteProcess(photo_nuclear_process);
35}
36} // namespace simcore
37
38DECLARE_PHOTONUCLEAR_MODEL(simcore::BertiniAtLeastNProductsModel)
void constructGammaProcess(G4ProcessManager *processManager)
The primary part of the model interface, responsible for adding the desired G4HadronicInteraction to ...
virtual void addPNCrossSectionData(G4HadronInelasticProcess *process) const
Default implementation for adding XS data for the process.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...