11#ifndef SIMCORE_APRIMECONVERSIONTOFCPS_H_
12#define SIMCORE_APRIMECONVERSIONTOFCPS_H_
16#include "Framework/Logger.h"
17#include "G4ParticleDefinition.hh"
18#include "G4VDiscreteProcess.hh"
46 const G4String& processName =
"APrimeToFCPPair",
47 G4ProcessType type = fElectromagnetic);
59 G4bool
IsApplicable(
const G4ParticleDefinition& particle)
override;
84 G4double
GetMeanFreePath(
const G4Track& aTrack, G4double previousStepSize,
85 G4ForceCondition* condition)
override;
91 const G4Element* anElement);
97 const G4Step& aStep)
override;
108 const G4Material* aMaterial);
119 const G4Material* aMaterial);
Discrete process for A' -> fcp+ fcp- conversion in a nuclear field.
enableLogging("APrimeConversionToFCPs")
Enable logging for this class.
G4double getCrossSectionPerAtom(const G4DynamicParticle *aDynamicAPrime, const G4Element *anElement)
Return the total cross section per atom.
const G4ParticleDefinition * the_fcp_minus_
Pointer to fcp- definition.
G4double cross_sec_factor_
Factor to artificially scale the cross section.
const G4Element * selectRandomAtom(const G4DynamicParticle *aDynamicAPrime, const G4Material *aMaterial)
Select a random atom from the material, weighted by cross section.
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Generate the final state: A' -> fcp+ fcp-.
G4bool IsApplicable(const G4ParticleDefinition &particle) override
Check if this process applies to the given particle.
APrimeConversionToFCPs(const G4String &processName="APrimeToFCPPair", G4ProcessType type=fElectromagnetic)
Constructor.
G4double computeCrossSectionPerAtom(G4double APrimeEnergy, G4int Z)
Compute the cross section per atom for the given A' energy and Z.
G4double rc_
Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp.
static const std::string PROCESS_NAME
The name of this process.
void printInfoDefinition()
Print process information.
G4double lowest_energy_limit_
Absolute low energy limit of the model (2 * Mfcp)
void setCrossSecFactor(G4double fac)
Set a factor to artificially scale the cross section.
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
Return the mean free path of the A' in the current material.
G4double highest_energy_limit_
High energy limit of the model.
G4double computeMeanFreePath(G4double APrimeEnergy, const G4Material *aMaterial)
Compute the mean free path for the given energy and material.
const G4ParticleDefinition * the_fcp_plus_
Pointer to fcp+ definition.
G4double getCrossSecFactor() const
void BuildPhysicsTable(const G4ParticleDefinition &p) override
Build physics table (called during initialization).
const G4ParticleDefinition * the_a_prime_
Pointer to the A' particle definition.
G4double limit_energy_
Energy limit for accurate cross section (5 * Mfcp)
std::vector< G4double > temp_
Temporary vector for element selection.
~APrimeConversionToFCPs()=default
Destructor Just the default one.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...