Detailed interface

These classes manage the low-level runtime interface between Celeritas and Geant4.

class SharedParams

Shared (one instance for all threads) Celeritas problem data.

The CeleritasDisabled accessor queries the CELER_DISABLE environment variable as a global option for disabling Celeritas offloading.

This should be instantiated on the main thread during problem setup, preferably as a shared pointer. The shared pointer should be passed to a thread-local LocalTransporter instance. At the beginning of the run, after Geant4 has initialized physics data, the Initialize method must be called on the main thread to populate the Celeritas data structures (geometry, physics). InitializeWorker must subsequently be invoked on all worker threads to set up thread-local data (including CUDA device initialization and objects that require Geant4 allocators).

Some low-level objects, such as the output diagnostics and Geant4 geometry wrapper, can be created independently of Celeritas being enabled.

class LocalTransporter : public celeritas::TrackOffloadInterface

Manage offloading of tracks to Celeritas.

This class must be constructed locally on each worker thread/task/stream, usually as a shared pointer that’s accessible to:

  • a run action (for initialization),

  • an event action (to set the event ID and flush offloaded tracks at the end of the event)

  • a tracking action (to try offloading every track)

Warning

Due to Geant4 thread-local allocators, this class must be finalized or destroyed on the same CPU thread in which is created and used!

Interface utilities

struct AlongStepFactoryInput

Input argument to the AlongStepFactory interface.

When passed to a factory instance, all member data will be set (so the instance will be ‘true’).

Most of these classes have been forward-declared because they simply need to be passed along to another class’s constructor.

Public Functions

inline explicit operator bool() const

True if all data is assigned.

class ExceptionConverter

Translate Celeritas C++ exceptions into Geant4 G4Exception calls.

This should generally be used when wrapping calls to Celeritas in a user application.

For example, the user event action to transport particles on device could be used as:

void EventAction::EndOfEventAction(const G4Event*)
{
    // Transport any tracks left in the buffer
    celeritas::ExceptionConverter call_g4exception{"celer.event.flush"};
    CELER_TRY_HANDLE(transport_->Flush(), call_g4exception);

    // Debug error checking
    CELER_ENSURE(transport_.GetBufferSize() == 0
                 || call_g4exception.forwarded());
}

class AlongStepFactoryInterface

Helper class for emitting an AlongStep action.

Currently Celeritas accepts a single custom along-step action (i.e., the same stepper is used for charged particles across all energies and regions of the problem). The along-step action is a single GPU kernel that combines the field stepper selection, the magnetic field, slowing-down calculation, multiple scattering, and energy loss fluctuations.

The factory will be called from the thread that initializes SharedParams. Instead of a daughter class, you can provide any function-like object that has the same interface.

Celeritas provides a few “default” configurations of along-step actions in celeritas/alongstep.

Subclassed by celeritas::CartMapFieldAlongStepFactory, celeritas::CylMapFieldAlongStepFactory, celeritas::RZMapFieldAlongStepFactory, celeritas::UniformAlongStepFactory

Classes usable by Geant4

These utilities are based on Celeritas data structures and capabilities but are written to be usable both by the celer-g4 app and potential other users.

Fields

template<class P, class F>
class MagneticField : public G4MagneticField

Wrap a Celeritas field as a Geant4 magnetic field.

Template Parameters:
  • P – params for creating field

  • F – field calculator

using celeritas::RZMapMagneticField = celeritas::MagneticField<RZMapFieldParams, RZMapField>

Geant4 magnetic field class.

using celeritas::CartMapMagneticField = celeritas::MagneticField<CartMapFieldParams, CartAdapterField>

Geant4 magnetic field class for XYZ uniform grid field.

using celeritas::CylMapMagneticField = celeritas::MagneticField<CylMapFieldParams, CylMapField>

Geant4 magnetic field adapter for cylindrical field.

CartMapFieldParams::Input celeritas::MakeCartMapFieldInput(CartMapFieldGridParams const &params)

Generates input for CartMapField params with configurable uniform grid dimensions in native Geant4 units.

This must be called after G4RunManager::Initialize as it will retrieve the G4FieldManager’s field to sample it.

CylMapFieldParams::Input celeritas::MakeCylMapFieldInput(std::vector<G4double> const &r_grid, std::vector<G4double> const &phi_values, std::vector<G4double> const &z_grid)

Generates input for CylMapField params with configurable nonuniform grid dimensions in native Geant4 units, and \(\phi\) should be in the range [0; \(2\times\pi\)].

This must be called after G4RunManager::Initialize as it will retrieve the G4FieldManager’s field to sample it.

Primary generators

class HepMC3PrimaryGenerator : public G4VPrimaryGenerator

HepMC3 reader class for sharing across threads.

This class should be shared among threads so that events can be correctly split up between them. It should be called from a user’s primary generator action:

void MyAction::GeneratePrimaries(G4Event* event)
{
    CELER_TRY_HANDLE(generator_->GeneratePrimaryVertex(event),
                     ExceptionConverter("celer.event.generate"));
}

Note

This class assumes that all threads will be reading all events sequentially and that events in the HepMC3 file are numbered sequentially from zero.

class PGPrimaryGeneratorAction : public G4VUserPrimaryGeneratorAction

Generate events from a particle gun.

This generates primary particles with energy, position, and direction sampled from distributions specified by the user in PrimaryGeneratorOptions.

See also

PrimaryGenerator

Physics lists

Two physics constructors build exclusively processes supported by Celeritas for Geant4:

class SupportedEmStandardPhysics : public G4VPhysicsConstructor

Construct G4EmStandardPhysics processes that are implemented in Celeritas.

This physics list is targeted at HEP experiments and reproduces most of the Geant4 G4EmStandardPhysics.

Limitations:

  • No support for generic ions

  • No hadronic EM interactions

  • Wentzel VI MSC is not supported

  • No polarized gamma processes

Electron/positron processes:

Processes

Model classes

Pair annihilation

G4eeToTwoGammaModel

Ionization

G4MollerBhabhaModel

Bremsstrahlung (low E)

G4SeltzerBergerModel

Bremsstrahlung (high E)

G4eBremsstrahlungRelModel

Coulomb scattering

G4eCoulombScatteringModel

Multiple scattering (low E)

G4UrbanMscModel

Multiple scattering (high E)

G4WentzelVIModel

Gamma processes:

Processes

Model classes

Compton scattering

G4KleinNishinaCompton

Photoelectric effect

G4LivermorePhotoElectricModel

Rayleigh scattering

G4LivermoreRayleighModel

Gamma conversion

G4PairProductionRelModel

If the gamma_general option is enabled, we create a single unified G4GammaGeneralProcess process, which embeds these other processes and calculates a combined total cross section. It’s faster in Geant4 but shouldn’t result in statistically different answers.

Muon processes are disabled by default:

Processes

Model classes

Pair production

G4MuPairProductionModel

Ionization (low E, mu-)

G4ICRU73QOModel

Ionization (low E, mu+)

G4BraggModel

Ionization (high E)

G4MuBetheBlochModel

Bremsstrahlung

G4MuBremsstrahlungModel

Coulomb scattering

G4eCoulombScatteringModel

Multiple scattering

G4WentzelVIModel

Note

Prior to version 11.1.0, Geant4 used the G4BetheBlochModel for muon ionization between 200 keV and 1 GeV and the G4MuBetheBlochModel above 1 GeV. Since version 11.1.0, the G4MuBetheBlochModel is used for all energies above 200 keV.

Warning

For muon-catalyzed fusion physics, the G4MuonMinusAtomicCapture , which is a G4VRestProcess with G4ProcessType::fHadronic , is added for mu- and is the exception to the EM-only rule. It does not require importing any cross section tables, but its inclusion enables muCF physics. This design decision may need revisiting in the future.

class SupportedOpticalPhysics : public G4VPhysicsConstructor

Construct Celeritas-supported optical physics processes.

Two “modular” physics lists (one using Geant4 hadronics, the other using pure Celeritas) are stand-ins for physics factories suitable for sending to G4RunManager::SetUserInitialization.

class EmPhysicsList : public G4VModularPhysicsList

Construct highly configurable EM and/or optical physics.

class FtfpBertPhysicsList : public G4VModularPhysicsList

Construct the FTFP_BERT physics list with configurable EM standard physics.

Physics setup

The input options incorporate process and model selection as well as default EM parameters to send to Geant4.

struct GeantPhysicsOptions

Construction options for Geant physics.

These options attempt to default to our closest match to G4StandardEmPhysics. They are passed to the EmPhysicsList and FtfpBert physics lists to provide an easy way to set up physics options.

Todo:

This will be moved to celeritas::inp

Rename default_cutoff to be consistent (it’s a production cut)

Gamma physics

bool compton_scattering = {true}

Enable Compton scattering.

bool photoelectric = {true}

Enable the photoelectric effect.

bool rayleigh_scattering = {true}

Enable Rayleigh scattering.

bool gamma_conversion = {true}

Enable electron pair production.

bool gamma_general = {false}

Use G4GammaGeneral instead of individual gamma processes.

Electron and positron physics

bool coulomb_scattering = {false}

Enable discrete Coulomb.

bool ionization = {true}

Enable e- and e+ ionization.

bool annihilation = {true}

Enable positron annihilation.

BremsModelSelection brems = {BremsModelSelection::all}

Enable bremsstrahlung and select a model.

MevEnergy seltzer_berger_limit = {1e3}

Upper limit for the Seltzer-Berger bremsstrahlung model.

MscModelSelection msc = {MscModelSelection::urban}

Enable multiple coulomb scattering and select a model.

Electron/positron MSC requires ionization

RelaxationSelection relaxation = {RelaxationSelection::none}

Enable atomic relaxation and select a model.

Physics options

int em_bins_per_decade = {7}

Number of log-spaced bins per factor of 10 in energy.

bool eloss_fluctuation = {true}

Enable universal energy fluctuations.

bool lpm = {true}

Apply relativistic corrections for select models.

bool integral_approach = {true}

See PhysicsParamsOptions::disable_integral_xs.

Cutoff options

MevEnergy min_energy = {0.1 * 1e-3}

Lowest energy of any EM physics process.

MevEnergy max_energy = {100 * 1e6}

Highest energy of any EM physics process.

double linear_loss_limit = {0.01}

See PhysicsParamsOptions::linear_loss_limit.

MevEnergy lowest_electron_energy = {0.001}

Tracking cutoff kinetic energy for e-/e+.

MevEnergy lowest_muhad_energy = {0.001}

Tracking cutoff kinetic energy for muons/hadrons.

bool apply_cuts = {false}

Kill secondaries below the production cut.

double default_cutoff = {0.1 * units::centimeter}

Set the default production cut for all particle types [len].

Multiple scattering configuration

double msc_range_factor = {0.04}

e-/e+ range factor for MSC models

double msc_muhad_range_factor = {0.2}

Muon/hadron range factor for MSC models.

double msc_safety_factor = {0.6}

Safety factor for MSC models.

double msc_lambda_limit = {0.1 * units::centimeter}

Lambda limit for MSC models [len].

double msc_theta_limit = {constants::pi}

Polar angle limit between single and multiple Coulomb scattering.

double angle_limit_factor = {1}

Factor for dynamic computation of angular limit between SS and MSC.

bool msc_displaced = {true}

Whether lateral displacement is enabled for e-/e+ MSC.

bool msc_muhad_displaced = {false}

Whether lateral displacement is enabled for muon/hadron MSC.

MscStepLimitAlgorithm msc_step_algorithm = {MscStepLimitAlgorithm::safety}

Step limit algorithm for e-/e+ MSC models.

MscStepLimitAlgorithm msc_muhad_step_algorithm{MscStepLimitAlgorithm::minimal}

Step limit algorithm for muon/hadron MSC models.

NuclearFormFactorType form_factor = {NuclearFormFactorType::exponential}

Nuclear form factor model for Coulomm scattering.

Public Members

GeantMuonPhysicsOptions muon = {GeantMuonPhysicsOptions::deactivated()}

Muon EM physics.

bool mucf_physics = {false}

Muon-catalyzed fusion physics.

bool verbose = {false}

Print detailed Geant4 output.

GeantOpticalPhysicsOptions optical = {GeantOpticalPhysicsOptions::deactivated()}

Optical physics options.

Public Static Functions

static inline GeantPhysicsOptions deactivated()

Initialize with no physics.

Detector construction

class DetectorConstruction : public G4VUserDetectorConstruction

Load a GDML file and construct sensitive detectors.

  • In Construct on the main thread, load the GDML file (including detectors)

  • In ConstructSDandField on each worker thread, call the SDBuilder for each distinct SD name, for all LVs that share the SD name