LDMX Software
BertiniSingleNeutronModel.cxx
1
2#include "SimCore/PhotoNuclearModels/BertiniSingleNeutronModel.h"
3namespace simcore {
4
5bool BertiniSingleNeutronProcess::acceptEvent() const {
6 int Nhard{0};
7 int Nhard_neutron{0};
8 int secondaries{theParticleChange.GetNumberOfSecondaries()};
9 for (int i{0}; i < secondaries; ++i) {
10 const auto secondary{theParticleChange.GetSecondary(i)->GetParticle()};
11 const auto pdgCode{secondary->GetDefinition()->GetPDGEncoding()};
12 if (skipCountingParticle(pdgCode)) {
13 continue;
14 }
15 const auto energy{secondary->GetKineticEnergy()};
16 if (energy > threshold_) {
17 Nhard++;
18 if (pdgCode == 2112) {
19 Nhard_neutron++;
20 }
21 }
22 }
23 return Nhard == 1 && Nhard_neutron == 1;
24}
25
27 G4ProcessManager* processManager) {
28 auto photoNuclearProcess{
29 new G4HadronInelasticProcess("photonNuclear", G4Gamma::Definition())};
30 auto model{new BertiniSingleNeutronProcess{threshold_, Zmin_, Emin_,
31 count_light_ions_}};
32 model->SetMaxEnergy(15 * CLHEP::GeV);
33 addPNCrossSectionData(photoNuclearProcess);
34 photoNuclearProcess->RegisterMe(model);
35 processManager->AddDiscreteProcess(photoNuclearProcess);
36}
37} // namespace simcore
38
39DECLARE_PHOTONUCLEAR_MODEL(simcore::BertiniSingleNeutronModel)
constexpr bool skipCountingParticle(const int pdgcode) const
Whether or not to include a particular particle type in any counting.
void ConstructGammaProcess(G4ProcessManager *processManager) override
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.