|
FairShip
|
Functions | |
| dict[str, str] | _merge_env (Mapping[str, str|None]|None env_vars) |
| subprocess.CompletedProcess | _run (Sequence[str] args, Mapping[str, str|None]|None env_vars=None, *bool check=True) |
| str | get_1D_flux_name (int nupdg) |
| str | get_2D_flux_name (int nupdg) |
| subprocess.CompletedProcess | make_splines (Sequence[int] nupdglist, str targetcode, int|float emax, int nknots, PathLike outputfile, Mapping[str, str|None]|None env_vars=None) |
| subprocess.CompletedProcess | generate_genie_events (int nevents, int nupdg, int|float emin, int|float emax, str targetcode, PathLike inputflux, PathLike spline, PathLike outputfile, *str|None process=None, int|None seed=None, int|None irun=None, Mapping[str, str|None]|None env_vars=None) |
| subprocess.CompletedProcess | make_ntuples (PathLike inputfile, PathLike outputfile, Mapping[str, str|None]|None env_vars=None) |
| None | add_hists (PathLike inputflux, PathLike simfile, int nupdg) |
| dict[str, str|None] | _parse_env_kv (Sequence[str] pairs) |
| int | main (Sequence[str]|None argv=None) |
Variables | |
| list | __all__ |
| str | __version__ = "0.1.0" |
| str | PathLike = str | os.PathLike |
| logging | logger = logging.getLogger(__name__) |
Utilities to run GENIE tools with per-call environment overrides.
This module wraps common GENIE command-line tools (``gmkspl``, ``gevgen``,
``gntpc``) and a ROOT helper to copy flux histograms into simulation files.
It provides:
- Safe subprocess execution (no ``shell=True``)
- Per-call environment overrides (set/unset variables like ``GXMLPATH``)
- Helpful printing/logging of the exact command line
Notes:
-----
* A Python process cannot modify its parent shell environment. Use the
``env_vars=...`` argument to set/unset variables **for the child process**
only. If you need to persist environment variables, put them in your shell
rc (e.g. ``~/.bashrc``) or wrap execution in ``source/env`` scripts.
* The path to ``Messenger_laconic.xml`` is auto-detected from ``$GENIE`` if
defined and readable; otherwise the ``--message-thresholds`` option is
omitted.
Example:
-------
>>> from genie_interface import generate_genie_events
>>> generate_genie_events(
... nevents=1000, nupdg=14, emin=1, emax=100, targetcode="1000080160",
... inputflux="/path/flux.root", spline="/path/xsec.xml", outputfile="out",
... env_vars={"GXMLPATH": "/opt/genie/config", "GENIE": "/cvmfs/.../genie"}
... )
|
protected |
Merge ``env_vars`` with the current environment.
Parameters
----------
env_vars:
Mapping of KEY -> VALUE to set for the child process. Use ``None`` to
**unset** a variable (remove it) for the child environment.
Returns
-------
dict
The environment dict to pass into ``subprocess.run``.
Definition at line 64 of file genie_interface.py.
|
protected |
Parse ``KEY=VAL`` pairs for ``--env`` CLI flag.
Definition at line 374 of file genie_interface.py.
|
protected |
Run a command with optional environment overrides.
Parameters
----------
args:
Command and arguments, each as a separate list element. (No shell.)
env_vars:
Optional mapping of KEY -> VALUE to inject into the child environment.
Use ``None`` value to **unset** a variable.
check:
If True (default), raise ``CalledProcessError`` on non-zero exit.
Returns
-------
subprocess.CompletedProcess
The completed process object.
Raises
------
subprocess.CalledProcessError
If the command exits with non-zero status and ``check=True``.
Definition at line 89 of file genie_interface.py.
Copy the 2D p–pT flux histogram from ``inputflux`` into ``simfile``.
Parameters
----------
inputflux:
Path to ROOT file containing TH2D flux histograms.
simfile:
Path to the simulation ROOT file to **update** (opened in 'UPDATE' mode).
nupdg:
Neutrino PDG code used to select the histogram (see ``get_2D_flux_name``).
Raises
------
FileNotFoundError
If the requested histogram is not found in ``inputflux``.
RuntimeError
If ROOT fails to open either file.
Definition at line 330 of file genie_interface.py.
| subprocess.CompletedProcess genie_interface.generate_genie_events | ( | int | nevents, |
| int | nupdg, | ||
| int | float | emin, | ||
| int | float | emax, | ||
| str | targetcode, | ||
| PathLike | inputflux, | ||
| PathLike | spline, | ||
| PathLike | outputfile, | ||
| *str | None | process = None, |
||
| int | None | seed = None, |
||
| int | None | irun = None, |
||
| Mapping[str, str | None] | None | env_vars = None |
||
| ) |
Run GENIE ``gevgen`` to generate events.
Parameters
----------
nevents:
Number of events to generate.
nupdg:
Neutrino PDG code (e.g. 12, 14, 16 or negatives for anti-neutrinos).
emin, emax:
Energy range in GeV (inclusive, inclusive as accepted by GENIE).
targetcode:
GENIE nuclear target code.
inputflux:
Path to ROOT file with flux histograms.
spline:
Path to GENIE cross-section spline XML file.
outputfile:
Output prefix/file for GENIE (as accepted by ``-o``).
process:
Optional GENIE generator list (e.g. ``'CCDIS'``, ``'NC'``).
seed:
Optional RNG seed for reproducibility.
irun:
Optional integer run number (``--run``).
env_vars:
Per-call environment overrides (e.g. set/unset ``GXMLPATH``, ``GENIE``).
Returns
-------
subprocess.CompletedProcess
Notes
-----
If ``GENIE`` is defined in the child environment, this function will try
to pass ``--message-thresholds $GENIE/config/Messenger_laconic.xml`` if
the file exists; otherwise that option is omitted.
Definition at line 210 of file genie_interface.py.
| str genie_interface.get_1D_flux_name | ( | int | nupdg | ) |
Return the TH1D flux histogram name for a given neutrino PDG. Naming convention (as used by input flux ROOT files): - Neutrino (PDG > 0): ``'10' + |pdg|`` - Anti-neutrino (PDG < 0): ``'20' + |pdg|`` Examples -------- >>> get_1D_flux_name(12) # nu_e '1012' >>> get_1D_flux_name(-12) # anti-nu_e '2012'
Definition at line 125 of file genie_interface.py.
| str genie_interface.get_2D_flux_name | ( | int | nupdg | ) |
Return the TH2D p–pT flux histogram name for a given neutrino PDG. Naming convention (as used by input flux ROOT files): - Neutrino (PDG > 0): ``'12' + |pdg|`` - Anti-neutrino (PDG < 0): ``'22' + |pdg|`` Examples -------- >>> get_2D_flux_name(12) # nu_e '1212' >>> get_2D_flux_name(-12) # anti-nu_e '2212'
Definition at line 144 of file genie_interface.py.
| int genie_interface.main | ( | Sequence[str] | None | argv = None | ) |
Entry point for the command-line interface.
Definition at line 385 of file genie_interface.py.
| subprocess.CompletedProcess genie_interface.make_ntuples | ( | PathLike | inputfile, |
| PathLike | outputfile, | ||
| Mapping[str, str | None] | None | env_vars = None |
||
| ) |
Convert a GENIE ``.ghep.root`` file to GST via ``gntpc``.
Parameters
----------
inputfile:
Input GHEP file (e.g. ``gntp.0.ghep.root``).
outputfile:
Output GST file path.
env_vars:
Per-call environment overrides.
Returns
-------
subprocess.CompletedProcess
Definition at line 305 of file genie_interface.py.
| subprocess.CompletedProcess genie_interface.make_splines | ( | Sequence[int] | nupdglist, |
| str | targetcode, | ||
| int | float | emax, | ||
| int | nknots, | ||
| PathLike | outputfile, | ||
| Mapping[str, str | None] | None | env_vars = None |
||
| ) |
Generate GENIE cross-section splines via ``gmkspl``.
Parameters
----------
nupdglist:
Iterable of neutrino PDG codes (e.g. ``[12, 14, 16, -12, -14, -16]``).
targetcode:
GENIE nuclear target code (e.g. ``'1000080160'`` for O-16).
emax:
Maximum neutrino energy (GeV).
nknots:
Number of knots for the spline generation.
outputfile:
Path to output XML spline file.
env_vars:
Per-call environment overrides (e.g. set ``GXMLPATH`` or ``GENIE``).
Returns
-------
subprocess.CompletedProcess
Definition at line 163 of file genie_interface.py.
|
private |
Definition at line 46 of file genie_interface.py.
|
private |
Definition at line 57 of file genie_interface.py.
| logging genie_interface.logger = logging.getLogger(__name__) |
Definition at line 61 of file genie_interface.py.
| str genie_interface.PathLike = str | os.PathLike |
Definition at line 59 of file genie_interface.py.