LDMX Software
APrimePhysics.cxx
Go to the documentation of this file.
1
9
10#include "G4DarkBreM/G4APrime.h"
11#include "G4DarkBreM/G4DarkBreMModel.h"
12#include "SimCore/UserEventInformation.h"
13
14// Geant4
15#include "G4Electron.hh"
16#include "G4EventManager.hh"
17#include "G4ProcessManager.hh"
18
19namespace simcore {
20
21const std::string APrimePhysics::NAME = "APrime";
22
33static void store_element_z(const G4Element& element) {
34 static_cast<UserEventInformation*>(
35 G4EventManager::GetEventManager()->GetUserInformation())
36 ->setDarkBremMaterialZ(element.GetZ());
37}
38
40 : G4VPhysicsConstructor(APrimePhysics::NAME),
41 parameters_{params},
42 process_{nullptr} {
43 ap_mass_ = parameters_.getParameter<double>("ap_mass", 0.) * MeV;
44 enable_ = parameters_.getParameter<bool>("enable", false);
45}
46
49 static const std::map<std::string, G4APrime::DecayMode> decay_lut = {
50 {"no_decay", G4APrime::DecayMode::NoDecay},
51 {"flat_decay", G4APrime::DecayMode::FlatDecay},
52 {"geant_decay", G4APrime::DecayMode::GeantDecay}};
53 auto decay_it{decay_lut.find(
54 model.getParameter<std::string>("decay_mode", "no_decay"))};
55 if (decay_it == decay_lut.end()) {
56 EXCEPTION_RAISE(
57 "BadConf",
58 "Unrecognized decay mode '" +
59 model.getParameter<std::string>("decay_mode") +
60 "',"
61 " options are 'no_decay', 'flat_decay', or 'geant_decay'.");
62 }
63
64 double ap_tau = model.getParameter<double>("ap_tau", -1.0);
65
73 G4APrime::Initialize(ap_mass_, 622, ap_tau, decay_it->second);
74}
75
77 // add process to electron if we are enabled
78 if (enable_) {
79 auto model{
81 auto model_name{model.getParameter<std::string>("name")};
82 if (model_name == "vertex_library" or model_name == "g4db") {
83 static const std::map<std::string, g4db::G4DarkBreMModel::ScalingMethod>
84 method_lut = {
85 {"forward_only",
86 g4db::G4DarkBreMModel::ScalingMethod::ForwardOnly},
87 {"cm_scaling", g4db::G4DarkBreMModel::ScalingMethod::CMScaling},
88 {"undefined", g4db::G4DarkBreMModel::ScalingMethod::Undefined}};
89 auto scaling_method_it{
90 method_lut.find(model.getParameter<std::string>("method"))};
91 if (scaling_method_it == method_lut.end()) {
92 EXCEPTION_RAISE(
93 "BadConf",
94 "Unrecognized scaling method '" +
95 model.getParameter<std::string>("method") +
96 "',"
97 " options are 'forward_only', 'cm_scaling', or 'undefined'.");
98 }
99 // Note: The process variable isn't used here, but creating the
100 // G4DarkBremsstahlung object has side-effects
101 process_ = std::make_unique<G4DarkBremsstrahlung>(
102 std::make_shared<g4db::G4DarkBreMModel>(
103 model.getParameter<std::string>("library_path"),
104 false /* dark brem off muons instead of electrons - we
105 always DB off electrons here */
106 ,
107 model.getParameter<double>("threshold"),
108 model.getParameter<double>("epsilon"), scaling_method_it->second,
109 g4db::G4DarkBreMModel::XsecMethod::Auto,
110 model.getParameter<double>("max_R_for_full", 50.0),
111 model.getParameter<int>("aprime_lhe_id", 622),
112 true, // always load the library
113 model.getParameter<bool>("scale_APrime", false),
114 model.getParameter<double>("dist_decay_min", 0.0),
115 model.getParameter<double>("dist_decay_max", 1.0)),
116 parameters_.getParameter<bool>("only_one_per_event"),
117 1., /* global bias - should use bias operator instead */
118 parameters_.getParameter<bool>("cache_xsec"));
119 process_->RegisterStorageMechanism(store_element_z);
120 } else {
121 EXCEPTION_RAISE("BadConf",
122 "Unrecognized model name '" + model_name + "'.");
123 }
124 G4cout << "[ APrimePhysics ] : Initialization of dark brem complete"
125 << G4endl;
126 }
127}
128
129} // namespace simcore
static void store_element_z(const G4Element &element)
Store the atomic Z for the element in which the dark brem occurred.
Class which defines basic APrime physics.
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:27
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
Defines basic APrime physics.
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.
bool enable_
is dark brem enabled for this run?
void ConstructProcess()
Construct the process.
Encapsulates user defined information associated with a Geant4 event.