LDMX Software
GammaConversionToFCPs.h
Go to the documentation of this file.
1
11#ifndef SIMCORE_GAMMACONVERSIONTOFCPS_H_
12#define SIMCORE_GAMMACONVERSIONTOFCPS_H_
13
14#include <vector>
15
16#include "Framework/Logger.h"
17#include "G4ParticleDefinition.hh"
18#include "G4VDiscreteProcess.hh"
19
20namespace simcore {
21
31class GammaConversionToFCPs : public G4VDiscreteProcess {
32 public:
37 static const std::string PROCESS_NAME;
38
45 explicit GammaConversionToFCPs(const G4String& processName = "GammaToFCPPair",
46 G4ProcessType type = fElectromagnetic);
47
53
58 G4bool IsApplicable(const G4ParticleDefinition& particle) override;
59
64 void BuildPhysicsTable(const G4ParticleDefinition& p) override;
65
70
75 void setCrossSecFactor(G4double fac);
76
78 inline G4double getCrossSecFactor() const { return cross_sec_factor_; }
79
83 G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
84 G4ForceCondition* condition) override;
85
89 G4double getCrossSectionPerAtom(const G4DynamicParticle* aDynamicGamma,
90 const G4Element* anElement);
91
95 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
96 const G4Step& aStep) override;
97
101 G4double computeCrossSectionPerAtom(G4double GammaEnergy, G4int Z);
102
106 G4double computeMeanFreePath(G4double GammaEnergy,
107 const G4Material* aMaterial);
108
109 // delete copy constructor and assignment
111 GammaConversionToFCPs& operator=(const GammaConversionToFCPs&) = delete;
112
113 private:
117 const G4Element* selectRandomAtom(const G4DynamicParticle* aDynamicGamma,
118 const G4Material* aMaterial);
119
121 G4double mfcp_;
123 G4double rc_;
131 G4double cross_sec_factor_ = 1.0;
132
134 const G4ParticleDefinition* the_fcp_minus_;
136 const G4ParticleDefinition* the_fcp_plus_;
137
139 std::vector<G4double> temp_;
140
144 enableLogging("GammaConversionToFCPs");
145};
146
147} // namespace simcore
148
149#endif // SIMCORE_GAMMACONVERSIONTOFCPS_H_
Discrete process for gamma -> fcp+ fcp- conversion in a nuclear field.
G4double getCrossSectionPerAtom(const G4DynamicParticle *aDynamicGamma, const G4Element *anElement)
Return the total cross section per atom.
const G4Element * selectRandomAtom(const G4DynamicParticle *aDynamicGamma, const G4Material *aMaterial)
Select a random atom from the material, weighted by cross section.
G4double limit_energy_
Energy limit for accurate cross section (5 * Mfcp)
G4double computeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
Compute the mean free path for the given energy and material.
G4double computeCrossSectionPerAtom(G4double GammaEnergy, G4int Z)
Compute the cross section per atom for the given gamma energy and Z.
const G4ParticleDefinition * the_fcp_plus_
Pointer to fcp+ definition.
enableLogging("GammaConversionToFCPs")
Enable logging for this class.
G4bool IsApplicable(const G4ParticleDefinition &particle) override
Check if this process applies to the given particle.
G4double lowest_energy_limit_
Absolute low energy limit of the model (2 * Mfcp)
G4double rc_
Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp.
G4double highest_energy_limit_
High energy limit of the model.
void BuildPhysicsTable(const G4ParticleDefinition &p) override
Build physics table (called during initialization).
std::vector< G4double > temp_
Temporary vector for element selection.
void printInfoDefinition()
Print process information.
GammaConversionToFCPs(const G4String &processName="GammaToFCPPair", G4ProcessType type=fElectromagnetic)
Constructor.
G4double cross_sec_factor_
Factor to artificially scale the cross section.
static const std::string PROCESS_NAME
The name of this process.
void setCrossSecFactor(G4double fac)
Set a factor to artificially scale the cross section.
const G4ParticleDefinition * the_fcp_minus_
Pointer to fcp- definition.
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
Return the mean free path of the gamma in the current material.
~GammaConversionToFCPs()=default
Destructor Just the default one.
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Generate the final state: gamma -> fcp+ fcp-.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...