LDMX Software
APrimePhysics.cxx
Go to the documentation of this file.
1
9
10namespace simcore {
11
12const std::string APrimePhysics::NAME = "APrime";
13
24static void storeElementZ(const G4Element& element) {
25 static_cast<UserEventInformation*>(
26 G4EventManager::GetEventManager()->GetUserInformation())
27 ->setDarkBremMaterialZ(element.GetZ());
28}
29
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}
48
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}
96
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...";
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}
178
179} // namespace simcore
Class which defines basic APrime physics.
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
Discrete process for A' -> fcp+ fcp- conversion in a nuclear field.
void setCrossSecFactor(G4double fac)
Set a factor to artificially scale the cross section.
Defines basic APrime physics.
bool fcp_enable_
is A' -> fcp conversion enabled for this run?
APrimeConversionToFCPs * fcp_conversion_process_
A' -> fcp conversion process (owned by G4 process manager after registration)
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)
void ConstructParticle()
Construct particle.
G4double ap_mass_
the mass of the A' for this run
static const std::string NAME
The name of this physics constructor.
APrimePhysics(const framework::config::Parameters &params)
Class 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
void ConstructProcess()
Construct the process.
Encapsulates user defined information associated with a Geant4 event.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...
static void storeElementZ(const G4Element &element)
Store the atomic Z for the element in which the dark brem occurred.