|
LDMX Software
|
Discrete process for A' -> fcp+ fcp- conversion in a nuclear field. More...
#include <APrimeConversionToFCPs.h>
Public Member Functions | |
| APrimeConversionToFCPs (const G4String &processName="APrimeToFCPPair", G4ProcessType type=fElectromagnetic) | |
| Constructor. | |
| ~APrimeConversionToFCPs ()=default | |
| Destructor Just the default one. | |
| G4bool | IsApplicable (const G4ParticleDefinition &particle) override |
| Check if this process applies to the given particle. | |
| void | BuildPhysicsTable (const G4ParticleDefinition &p) override |
| Build physics table (called during initialization). | |
| void | printInfoDefinition () |
| Print process information. | |
| void | setCrossSecFactor (G4double fac) |
| Set a factor to artificially scale the cross section. | |
| G4double | getCrossSecFactor () const |
| G4double | GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override |
| Return the mean free path of the A' in the current material. | |
| G4double | getCrossSectionPerAtom (const G4DynamicParticle *aDynamicAPrime, const G4Element *anElement) |
| Return the total cross section per atom. | |
| G4VParticleChange * | PostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override |
| Generate the final state: A' -> fcp+ fcp-. | |
| G4double | computeCrossSectionPerAtom (G4double APrimeEnergy, G4int Z) |
| Compute the cross section per atom for the given A' energy and Z. | |
| G4double | computeMeanFreePath (G4double APrimeEnergy, const G4Material *aMaterial) |
| Compute the mean free path for the given energy and material. | |
| APrimeConversionToFCPs (const APrimeConversionToFCPs &)=delete | |
| APrimeConversionToFCPs & | operator= (const APrimeConversionToFCPs &)=delete |
Static Public Attributes | |
| static const std::string | PROCESS_NAME = "APrimeToFCPPair" |
| The name of this process. | |
Private Member Functions | |
| const G4Element * | selectRandomAtom (const G4DynamicParticle *aDynamicAPrime, const G4Material *aMaterial) |
| Select a random atom from the material, weighted by cross section. | |
| enableLogging ("APrimeConversionToFCPs") | |
| Enable logging for this class. | |
Private Attributes | |
| G4double | mfcp_ |
| fcp mass | |
| G4double | rc_ |
| Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp. | |
| G4double | limit_energy_ |
| Energy limit for accurate cross section (5 * Mfcp) | |
| G4double | lowest_energy_limit_ |
| Absolute low energy limit of the model (2 * Mfcp) | |
| G4double | highest_energy_limit_ |
| High energy limit of the model. | |
| G4double | cross_sec_factor_ = 1.0 |
| Factor to artificially scale the cross section. | |
| const G4ParticleDefinition * | the_a_prime_ |
| Pointer to the A' particle definition. | |
| const G4ParticleDefinition * | the_fcp_minus_ |
| Pointer to fcp- definition. | |
| const G4ParticleDefinition * | the_fcp_plus_ |
| Pointer to fcp+ definition. | |
| std::vector< G4double > | temp_ |
| Temporary vector for element selection. | |
Discrete process for A' -> fcp+ fcp- conversion in a nuclear field.
This process is the dark-sector analogue of G4GammaConversionToMuons (gamma -> mu+ mu-). A practically massless dark photon (A') converts into a pair of massive fractionally charged particles (fcp+ fcp-) in the electromagnetic field of a nucleus.
Definition at line 31 of file APrimeConversionToFCPs.h.
|
explicit |
Constructor.
| [in] | processName | name of this process (default: "APrimeToFCPPair") |
| [in] | type | Geant4 process type |
Definition at line 24 of file APrimeConversionToFCPs.cxx.
References highest_energy_limit_, limit_energy_, lowest_energy_limit_, mfcp_, rc_, the_a_prime_, the_fcp_minus_, and the_fcp_plus_.
|
override |
Build physics table (called during initialization).
Sets up the temp vector for element selection.
Definition at line 76 of file APrimeConversionToFCPs.cxx.
References printInfoDefinition(), and temp_.
| G4double simcore::APrimeConversionToFCPs::computeCrossSectionPerAtom | ( | G4double | APrimeEnergy, |
| G4int | Z ) |
Compute the cross section per atom for the given A' energy and Z.
Definition at line 141 of file APrimeConversionToFCPs.cxx.
References cross_sec_factor_, lowest_energy_limit_, mfcp_, and rc_.
Referenced by computeMeanFreePath(), getCrossSectionPerAtom(), and selectRandomAtom().
| G4double simcore::APrimeConversionToFCPs::computeMeanFreePath | ( | G4double | APrimeEnergy, |
| const G4Material * | aMaterial ) |
Compute the mean free path for the given energy and material.
Definition at line 105 of file APrimeConversionToFCPs.cxx.
References computeCrossSectionPerAtom(), limit_energy_, and lowest_energy_limit_.
Referenced by GetMeanFreePath().
|
inline |
Definition at line 79 of file APrimeConversionToFCPs.h.
References cross_sec_factor_.
| G4double simcore::APrimeConversionToFCPs::getCrossSectionPerAtom | ( | const G4DynamicParticle * | aDynamicAPrime, |
| const G4Element * | anElement ) |
Return the total cross section per atom.
Definition at line 135 of file APrimeConversionToFCPs.cxx.
References computeCrossSectionPerAtom().
|
override |
Return the mean free path of the A' in the current material.
Definition at line 92 of file APrimeConversionToFCPs.cxx.
References computeMeanFreePath().
|
override |
Check if this process applies to the given particle.
Definition at line 67 of file APrimeConversionToFCPs.cxx.
References the_a_prime_.
|
override |
Generate the final state: A' -> fcp+ fcp-.
Definition at line 188 of file APrimeConversionToFCPs.cxx.
References limit_energy_, lowest_energy_limit_, mfcp_, selectRandomAtom(), the_fcp_minus_, and the_fcp_plus_.
| void simcore::APrimeConversionToFCPs::printInfoDefinition | ( | ) |
Print process information.
Definition at line 455 of file APrimeConversionToFCPs.cxx.
References cross_sec_factor_, highest_energy_limit_, and lowest_energy_limit_.
Referenced by BuildPhysicsTable().
|
private |
Select a random atom from the material, weighted by cross section.
Definition at line 428 of file APrimeConversionToFCPs.cxx.
References computeCrossSectionPerAtom(), limit_energy_, and temp_.
Referenced by PostStepDoIt().
| void simcore::APrimeConversionToFCPs::setCrossSecFactor | ( | G4double | fac | ) |
Set a factor to artificially scale the cross section.
| [in] | fac | multiplicative factor (default 1.0) |
Definition at line 181 of file APrimeConversionToFCPs.cxx.
References cross_sec_factor_.
Referenced by simcore::APrimePhysics::ConstructProcess().
|
private |
Factor to artificially scale the cross section.
Definition at line 132 of file APrimeConversionToFCPs.h.
Referenced by computeCrossSectionPerAtom(), getCrossSecFactor(), printInfoDefinition(), and setCrossSecFactor().
|
private |
High energy limit of the model.
Definition at line 130 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), and printInfoDefinition().
|
private |
Energy limit for accurate cross section (5 * Mfcp)
Definition at line 126 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), computeMeanFreePath(), PostStepDoIt(), and selectRandomAtom().
|
private |
Absolute low energy limit of the model (2 * Mfcp)
Definition at line 128 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), computeCrossSectionPerAtom(), computeMeanFreePath(), PostStepDoIt(), and printInfoDefinition().
|
private |
fcp mass
Definition at line 122 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), computeCrossSectionPerAtom(), and PostStepDoIt().
|
static |
The name of this process.
Used to attach biasing operators to this process.
Definition at line 37 of file APrimeConversionToFCPs.h.
Referenced by simcore::biasoperators::APrimeToFCPPair::getProcessToBias().
|
private |
Coupling-scaled reduced Compton wavelength: q^2 * elm_coupling / Mfcp.
Definition at line 124 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), and computeCrossSectionPerAtom().
|
private |
Temporary vector for element selection.
Definition at line 142 of file APrimeConversionToFCPs.h.
Referenced by BuildPhysicsTable(), and selectRandomAtom().
|
private |
Pointer to the A' particle definition.
Definition at line 135 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), and IsApplicable().
|
private |
Pointer to fcp- definition.
Definition at line 137 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), and PostStepDoIt().
|
private |
Pointer to fcp+ definition.
Definition at line 139 of file APrimeConversionToFCPs.h.
Referenced by APrimeConversionToFCPs(), and PostStepDoIt().