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?
 
framework::config::Parameters parameters_
 Dark brem parameters to pass to the process (if enabled)
 
std::unique_ptr< G4DarkBremsstrahlung > process_
 

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 42 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}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
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.
bool enable_
is dark brem enabled for this run?

References ap_mass_, enable_, 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.

Definition at line 38 of file APrimePhysics.cxx.

38 {
39 auto model{parameters_.get<framework::config::Parameters>("model")};
40 static const std::map<std::string, G4APrime::DecayMode> decay_lut = {
41 {"no_decay", G4APrime::DecayMode::NoDecay},
42 {"flat_decay", G4APrime::DecayMode::FlatDecay},
43 {"geant_decay", G4APrime::DecayMode::GeantDecay}};
44 auto decay_it{
45 decay_lut.find(model.get<std::string>("decay_mode", "no_decay"))};
46 if (decay_it == decay_lut.end()) {
47 EXCEPTION_RAISE(
48 "BadConf",
49 "Unrecognized decay mode '" + model.get<std::string>("decay_mode") +
50 "',"
51 " options are 'no_decay', 'flat_decay', or 'geant_decay'.");
52 }
53
54 double ap_tau = model.get<double>("ap_tau", -1.0);
55
63 G4APrime::Initialize(ap_mass_, 622, ap_tau, decay_it->second);
64}
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29

References ap_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).

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

See also
G4DarkBremsstrahlung in G4DarkBreM

Definition at line 66 of file APrimePhysics.cxx.

66 {
67 // add process to electron if we are enabled
68 if (enable_) {
69 auto model{parameters_.get<framework::config::Parameters>("model")};
70 auto model_name{model.get<std::string>("name")};
71 if (model_name == "vertex_library" or model_name == "g4db") {
72 static const std::map<std::string, g4db::G4DarkBreMModel::ScalingMethod>
73 method_lut = {
74 {"forward_only",
75 g4db::G4DarkBreMModel::ScalingMethod::ForwardOnly},
76 {"cm_scaling", g4db::G4DarkBreMModel::ScalingMethod::CMScaling},
77 {"undefined", g4db::G4DarkBreMModel::ScalingMethod::Undefined}};
78 auto scaling_method_it{method_lut.find(model.get<std::string>("method"))};
79 if (scaling_method_it == method_lut.end()) {
80 EXCEPTION_RAISE(
81 "BadConf",
82 "Unrecognized scaling method '" + model.get<std::string>("method") +
83 "',"
84 " options are 'forward_only', 'cm_scaling', or 'undefined'.");
85 }
86 // Note: The process variable isn't used here, but creating the
87 // G4DarkBremsstahlung object has side-effects
88 process_ = std::make_unique<G4DarkBremsstrahlung>(
89 std::make_shared<g4db::G4DarkBreMModel>(
90 model.get<std::string>("library_path"),
91 false /* dark brem off muons instead of electrons - we
92 always DB off electrons here */
93 ,
94 model.get<double>("threshold"), model.get<double>("epsilon"),
95 scaling_method_it->second,
96 g4db::G4DarkBreMModel::XsecMethod::Auto,
97 model.get<double>("max_R_for_full", 50.0),
98 model.get<int>("aprime_lhe_id", 1023),
99 true, // always load the library
100 model.get<bool>("scale_APrime", false),
101 model.get<double>("dist_decay_min", 0.0),
102 model.get<double>("dist_decay_max", 1.0)),
103 parameters_.get<bool>("only_one_per_event"),
104 1., /* global bias - should use bias operator instead */
105 parameters_.get<bool>("cache_xsec"));
106 process_->RegisterStorageMechanism(storeElementZ);
107 } else {
108 EXCEPTION_RAISE("BadConf",
109 "Unrecognized model name '" + model_name + "'.");
110 }
111 ldmx_log(trace) << "Initialization of dark brem complete";
112 }
113}
static void storeElementZ(const G4Element &element)
Store the atomic Z for the element in which the dark brem occurred.

References enable_, framework::config::Parameters::get(), parameters_, 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 96 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 99 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 51 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 108 of file APrimePhysics.h.

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

◆ process_

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

Definition at line 110 of file APrimePhysics.h.


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