LDMX Software
simcore::APrimePhysics Class Reference

Defines basic APrime physics. More...

#include <APrimePhysics.h>

Public Member Functions

 APrimePhysics (const framework::config::Parameters &params)
 Class constructor.
 
virtual ~APrimePhysics ()=default
 Class destructor.
 
void ConstructParticle ()
 Construct particle.
 
void ConstructProcess ()
 Construct the process.
 

Static Public Attributes

static const std::string NAME = "APrime"
 The name of this physics constructor.
 

Private Attributes

G4double ap_mass_
 the mass of the A' for this run
 
bool enable_
 is dark brem enabled for this run?
 
bool fcp_enable_
 is A' -> fcp conversion enabled for this run?
 
G4double fcp_mass_
 mass of the fcp in MeV
 
G4double fcp_charge_
 charge of the fcp in units of e
 
G4double fcp_xsec_factor_
 cross section biasing factor for A' -> fcp conversion
 
framework::config::Parameters parameters_
 Dark brem parameters to pass to the process (if enabled)
 
std::unique_ptr< G4DarkBremsstrahlung > process_
 
APrimeConversionToFCPsfcp_conversion_process_ {nullptr}
 A' -> fcp conversion process (owned by G4 process manager after registration)
 

Detailed Description

Defines basic APrime physics.

It constructs the APrime particle and links the dark brem process to the electron.

See also
G4APrime
G4DarkBremsstrahlung in G4DarkBreM
Note
This class basically does not do anything except register the custom particle (G4APrime) and custom process (G4DarkBremsstrahlung).

Definition at line 44 of file APrimePhysics.h.

Constructor & Destructor Documentation

◆ APrimePhysics()

simcore::APrimePhysics::APrimePhysics ( const framework::config::Parameters & params)

Class constructor.

Parameters
paramsParameters to configure the dark brem process

Definition at line 30 of file APrimePhysics.cxx.

31 : G4VPhysicsConstructor(APrimePhysics::NAME),
32 parameters_{params},
33 process_{nullptr} {
34 ap_mass_ = parameters_.get<double>("ap_mass", 0.) * MeV;
35 enable_ = parameters_.get<bool>("enable", false);
36 fcp_enable_ = parameters_.get<bool>("fcp_enable", false);
37 fcp_mass_ = parameters_.get<double>("fcp_mass", 0.) * MeV;
38 fcp_charge_ = parameters_.get<double>("fcp_charge", 0.1);
39 fcp_xsec_factor_ = parameters_.get<double>("fcp_xsec_factor", 1.0);
40
41 ldmx_log(info) << "Configuration:" << " A' mass: " << ap_mass_ / MeV << " MeV"
42 << ", Dark brem enabled: " << (enable_ ? "YES" : "NO")
43 << ", FCP enabled: " << (fcp_enable_ ? "YES" : "NO")
44 << ", FCP mass: " << fcp_mass_ / MeV << " MeV"
45 << ", FCP charge: " << fcp_charge_ << " e"
46 << ", FCP xsec factor: " << fcp_xsec_factor_;
47}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
bool fcp_enable_
is A' -> fcp conversion enabled for this run?
G4double fcp_xsec_factor_
cross section biasing factor for A' -> fcp conversion
framework::config::Parameters parameters_
Dark brem parameters to pass to the process (if enabled)
G4double ap_mass_
the mass of the A' for this run
static const std::string NAME
The name of this physics constructor.
G4double fcp_mass_
mass of the fcp in MeV
bool enable_
is dark brem enabled for this run?
G4double fcp_charge_
charge of the fcp in units of e

References ap_mass_, enable_, fcp_charge_, fcp_enable_, fcp_mass_, fcp_xsec_factor_, framework::config::Parameters::get(), and parameters_.

◆ ~APrimePhysics()

virtual simcore::APrimePhysics::~APrimePhysics ( )
virtualdefault

Class destructor.

Nothing right now.

Member Function Documentation

◆ ConstructParticle()

void simcore::APrimePhysics::ConstructParticle ( )

Construct particle.

Insert A' into the Geant4 particle table. Geant4 registers all instances derived from G4ParticleDefinition and deletes them at the end of processing.

Uses the A' mass given by the parameter APrimeMass to inform the G4APrime instance what mass to use.

See also
G4APrime

Insert A-prime into the Geant4 particle table.

Geant4 registers all instances derived from G4ParticleDefinition and deletes them at the end of the run. We configure the A' to have the input mass and the PDG ID number of 622.

If fcp is enabled, the A' decays to fcp+fcp- (decay_id=17), otherwise to e+e- (decay_id=11).

If fcp is enabled, initialize the fractionally charged particles. This registers fcp- and fcp+ in the Geant4 particle table.

Definition at line 49 of file APrimePhysics.cxx.

49 {
50 auto model{parameters_.get<framework::config::Parameters>("model")};
51 static const std::map<std::string, G4APrime::DecayMode> decay_lut = {
52 {"no_decay", G4APrime::DecayMode::NoDecay},
53 {"flat_decay", G4APrime::DecayMode::FlatDecay},
54 {"geant_decay", G4APrime::DecayMode::GeantDecay}};
55 auto decay_it{
56 decay_lut.find(model.get<std::string>("decay_mode", "no_decay"))};
57 if (decay_it == decay_lut.end()) {
58 EXCEPTION_RAISE(
59 "BadConf",
60 "Unrecognized decay mode '" + model.get<std::string>("decay_mode") +
61 "',"
62 " options are 'no_decay', 'flat_decay', or 'geant_decay'.");
63 }
64
65 double ap_tau = model.get<double>("ap_tau", -1.0);
66
77 int decay_id = fcp_enable_ ? 17 : 11;
78 ldmx_log(info) << "Initializing A' with:" << ", Mass: " << ap_mass_ / MeV
79 << " MeV" << ", PDG ID: 622" << ", Lifetime (tau): " << ap_tau
80 << ", Decay mode: "
81 << model.get<std::string>("decay_mode", "no_decay")
82 << ", Decay product ID: " << decay_id
83 << (fcp_enable_ ? " (fcp+fcp-)" : " (e+e-)");
84
85 // Initialize the A' with the specified parameters
86 G4APrime::Initialize(ap_mass_, 622, ap_tau, decay_it->second, decay_id);
87
92 if (fcp_enable_) {
93 G4FractionallyCharged::Initialize(fcp_mass_, 17, fcp_charge_);
94 }
95}
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29

References ap_mass_, fcp_charge_, fcp_enable_, fcp_mass_, framework::config::Parameters::get(), and parameters_.

◆ ConstructProcess()

void simcore::APrimePhysics::ConstructProcess ( )

Construct the process.

Links the dark brem processs to the electron through the process manager only if the dark brem process is enabled ('enable' is True).

Also sets up the A' -> fcp+ fcp- conversion process if fcp is enabled.

G4ProcessManager registers and cleans up any created processes, so we can forget about it after creating it.

See also
G4DarkBremsstrahlung in G4DarkBreM
APrimeConversionToFCPs

Definition at line 97 of file APrimePhysics.cxx.

97 {
98 ldmx_log(trace) << "=== APrimePhysics::ConstructProcess ===";
99 // add process to electron if we are enabled
100 if (enable_) {
101 ldmx_log(debug)
102 << "Dark brem is enabled, setting up G4DarkBremsstrahlung...";
103 auto model{parameters_.get<framework::config::Parameters>("model")};
104 auto model_name{model.get<std::string>("name")};
105 if (model_name == "vertex_library" or model_name == "g4db") {
106 static const std::map<std::string, g4db::G4DarkBreMModel::ScalingMethod>
107 method_lut = {
108 {"forward_only",
109 g4db::G4DarkBreMModel::ScalingMethod::ForwardOnly},
110 {"cm_scaling", g4db::G4DarkBreMModel::ScalingMethod::CMScaling},
111 {"undefined", g4db::G4DarkBreMModel::ScalingMethod::Undefined}};
112 auto scaling_method_it{method_lut.find(model.get<std::string>("method"))};
113 if (scaling_method_it == method_lut.end()) {
114 EXCEPTION_RAISE(
115 "BadConf",
116 "Unrecognized scaling method '" + model.get<std::string>("method") +
117 "',"
118 " options are 'forward_only', 'cm_scaling', or 'undefined'.");
119 }
120 // Note: The process variable isn't used here, but creating the
121 // G4DarkBremsstahlung object has side-effects
122 process_ = std::make_unique<G4DarkBremsstrahlung>(
123 std::make_shared<g4db::G4DarkBreMModel>(
124 model.get<std::string>("library_path"),
125 false /* dark brem off muons instead of electrons - we
126 always DB off electrons here */
127 ,
128 model.get<double>("threshold"), model.get<double>("epsilon"),
129 scaling_method_it->second,
130 g4db::G4DarkBreMModel::XsecMethod::Auto,
131 model.get<double>("max_R_for_full", 50.0),
132 model.get<int>("aprime_lhe_id", 1023),
133 true, // always load the library
134 model.get<bool>("scale_APrime", false),
135 model.get<double>("dist_decay_min", 0.0),
136 model.get<double>("dist_decay_max", 1.0)),
137 parameters_.get<bool>("only_one_per_event"),
138 1., /* global bias - should use bias operator instead */
139 parameters_.get<bool>("cache_xsec"));
140 process_->RegisterStorageMechanism(storeElementZ);
141 } else {
142 EXCEPTION_RAISE("BadConf",
143 "Unrecognized model name '" + model_name + "'.");
144 }
145 ldmx_log(trace) << "Initialization of dark brem complete";
146 }
147
148 // Add A' -> fcp+ fcp- conversion process if enabled
149 if (fcp_enable_) {
150 ldmx_log(debug)
151 << "FCP is enabled, setting up A' -> fcp+ fcp- conversion process...";
152 G4ProcessManager* aprime_proc_man = G4APrime::APrime()->GetProcessManager();
153 if (aprime_proc_man == nullptr) {
154 ldmx_log(error) << "FATAL: Unable to access process manager for A'!";
155 EXCEPTION_RAISE("APrimePhysics",
156 "Was unable to access the process manager for A', "
157 "something is very wrong!");
158 }
159 ldmx_log(debug) << "Creating APrimeConversionToFCPs process...";
160 fcp_conversion_process_ = new APrimeConversionToFCPs();
161
162 // Apply cross section biasing factor
163 if (fcp_xsec_factor_ != 1.0) {
165 ldmx_log(debug) << "Applied cross section biasing factor: "
167 }
168
169 aprime_proc_man->AddDiscreteProcess(fcp_conversion_process_);
170 aprime_proc_man->SetProcessOrderingToFirst(
171 fcp_conversion_process_, G4ProcessVectorDoItIndex::idxAll);
172 ldmx_log(info)
173 << "A' -> fcp+ fcp- conversion process registered successfully!";
174 } else {
175 ldmx_log(info) << "FCP not enabled, A' conversion process not registered";
176 }
177}
void setCrossSecFactor(G4double fac)
Set a factor to artificially scale the cross section.
APrimeConversionToFCPs * fcp_conversion_process_
A' -> fcp conversion process (owned by G4 process manager after registration)
static void storeElementZ(const G4Element &element)
Store the atomic Z for the element in which the dark brem occurred.

References enable_, fcp_conversion_process_, fcp_enable_, fcp_xsec_factor_, framework::config::Parameters::get(), parameters_, simcore::APrimeConversionToFCPs::setCrossSecFactor(), and simcore::storeElementZ().

Member Data Documentation

◆ ap_mass_

G4double simcore::APrimePhysics::ap_mass_
private

the mass of the A' for this run

Definition at line 101 of file APrimePhysics.h.

Referenced by APrimePhysics(), and ConstructParticle().

◆ enable_

bool simcore::APrimePhysics::enable_
private

is dark brem enabled for this run?

Definition at line 104 of file APrimePhysics.h.

Referenced by APrimePhysics(), and ConstructProcess().

◆ fcp_charge_

G4double simcore::APrimePhysics::fcp_charge_
private

charge of the fcp in units of e

Definition at line 113 of file APrimePhysics.h.

Referenced by APrimePhysics(), and ConstructParticle().

◆ fcp_conversion_process_

APrimeConversionToFCPs* simcore::APrimePhysics::fcp_conversion_process_ {nullptr}
private

A' -> fcp conversion process (owned by G4 process manager after registration)

Definition at line 131 of file APrimePhysics.h.

131{nullptr};

Referenced by ConstructProcess().

◆ fcp_enable_

bool simcore::APrimePhysics::fcp_enable_
private

is A' -> fcp conversion enabled for this run?

Definition at line 107 of file APrimePhysics.h.

Referenced by APrimePhysics(), ConstructParticle(), and ConstructProcess().

◆ fcp_mass_

G4double simcore::APrimePhysics::fcp_mass_
private

mass of the fcp in MeV

Definition at line 110 of file APrimePhysics.h.

Referenced by APrimePhysics(), and ConstructParticle().

◆ fcp_xsec_factor_

G4double simcore::APrimePhysics::fcp_xsec_factor_
private

cross section biasing factor for A' -> fcp conversion

Definition at line 116 of file APrimePhysics.h.

Referenced by APrimePhysics(), and ConstructProcess().

◆ NAME

const std::string simcore::APrimePhysics::NAME = "APrime"
static

The name of this physics constructor.

Passed into Geant4 to register this physics, should not conflict with any other Geant4 physics.

Definition at line 53 of file APrimePhysics.h.

◆ parameters_

framework::config::Parameters simcore::APrimePhysics::parameters_
private

Dark brem parameters to pass to the process (if enabled)

Note
This can't be a reference because we pass it to the process after the configuration step from our POV is done. Thus we need our own copy that won't be destroyed.

Definition at line 125 of file APrimePhysics.h.

Referenced by APrimePhysics(), ConstructParticle(), and ConstructProcess().

◆ process_

std::unique_ptr<G4DarkBremsstrahlung> simcore::APrimePhysics::process_
private

Definition at line 127 of file APrimePhysics.h.


The documentation for this class was generated from the following files: