Skip to content

hallmd.models.thruster

Module for Hall thruster models.

Note

Only current implementation is for the 1d fluid Hallthruster.jl code. Other thruster codes can be implemented similarly here.

Includes:

  • run_hallthruster_jl - General wrapper to run HallThruster.jl for a single set of inputs
  • hallthruster_jl() - PEM wrapper to run HallThruster.jl for a set of PEM inputs
  • get_jl_env - Get the path of the julia environment created for HallThruster.jl for a specific git ref
  • PEM_TO_JULIA - Mapping of PEM variable names to a path in the HallThruster.jl input/output structure (defaults)

get_jl_binary()

use the juliaup config file to find out where the default julia binary is stored. if juliaup is not installed, returns 'julia'. this lets us avoid running juliaup's wrapper executable, which has led to concurrency issues. This is designed for use on the cluster, so on windows it just returns julia.

Source code in src/hallmd/models/thruster.py
def get_jl_binary():
    """use the juliaup config file to find out where the default julia binary is stored.
    if juliaup is not installed, returns 'julia'.
    this lets us avoid running juliaup's wrapper executable, which has led to concurrency issues.
    This is designed for use on the cluster, so on windows it just returns `julia`."""

    if platform.system() == "Windows":
        return "julia"

    home = os.environ["HOME"]
    juliaup_dir = Path(f"{home}/.julia/juliaup")
    juliaup_json = juliaup_dir / "juliaup.json"
    if not os.path.exists(juliaup_json):
        return "julia"

    with open(juliaup_json, "r") as fd:
        jl_config = json.load(fd)

    default = jl_config["Default"]
    versions = jl_config["InstalledVersions"]
    channels = jl_config["InstalledChannels"]

    version = versions[channels[default]["Version"]]
    path = juliaup_dir / version["Path"] / "bin" / "julia"
    return path

get_jl_env(git_ref)

Get the path of the julia environment created for HallThruster.jl for a specific git ref.

PARAMETER DESCRIPTION
git_ref

The git ref (i.e. commit hash, version tag, branch, etc.) of HallThruster.jl to use.

TYPE: str

Source code in src/hallmd/models/thruster.py
def get_jl_env(git_ref: str) -> Path:
    """Get the path of the julia environment created for HallThruster.jl for a specific git ref.

    :param git_ref: The git ref (i.e. commit hash, version tag, branch, etc.) of HallThruster.jl to use.
    """
    global_env_dir = Path("~/.julia/environments/").expanduser()
    env_path = global_env_dir / f"hallthruster_{git_ref}"
    return env_path

hallthruster_jl(thruster_inputs=None, thruster='SPT-100', config=None, simulation=None, postprocess=None, model_fidelity=(2, 2), output_path=None, version=HALLTHRUSTER_VERSION_DEFAULT, pem_to_julia='default', fidelity_function='default', julia_script=None, run_kwargs='default', shock_threshold=None)

Run a single HallThruster.jl simulation for a given set of inputs. This function will write a temporary input file to disk, call HallThruster.run_simulation() in Julia, and read the output file back into Python. Will return time-averaged performance metrics and ion velocity for use with the PEM.

Note that the specific inputs and outputs described here can be configured using the pem_to_julia dict.

Required configuration

You must specify a thruster, a domain, a mass flow rate, and a discharge voltage to run the simulation. The thruster must be defined in the hallmd.devices directory or as a dictionary with the required fields. The mass flow rate and discharge voltage are specified in thruster_inputs as mdot_a (kg/s) and V_a (V), respectively. The domain is specified as a list [left_bound, right_bound] in the config dictionary. See the HallThruster.jl docs for more details.

PARAMETER DESCRIPTION
thruster_inputs

named key-value pairs of thruster inputs: P_b, V_a, mdot_a, T_e, u_n, l_t, a_1, a_2, delta_z, z0, p0, and V_cc for background pressure (Torr), anode voltage, anode mass flow rate (kg/s), electron temperature (eV), neutral velocity (m/s), transition length (m), anomalous transport coefficients, and cathode coupling voltage. Will override the corresponding values in config if provided.

TYPE: Dataset DEFAULT: None

thruster

the name of the thruster to simulate (must be importable from hallmd.devices, see load_device), or a dictionary that provides geometry and magnetic field information of the thruster to simulate; see the Hallthruster.jl docs. Will override thruster in config if provided. If None, will defer to config. Defaults to the SPT-100.

TYPE: Literal['SPT-100'] | str | dict DEFAULT: 'SPT-100'

config

dictionary of configs for HallThruster.jl, see the Hallthruster.jl docs for options and formatting.

TYPE: Optional[dict] DEFAULT: None

simulation

dictionary of simulation parameters for HallThruster.jl

TYPE: Optional[dict] DEFAULT: None

postprocess

dictionary of post-processing parameters for Hallthruster.jl

TYPE: Optional[dict] DEFAULT: None

model_fidelity

tuple of integers that determine the number of cells and the number of charge states to use via ncells = model_fidelity[0] * 50 + 100 and ncharge = model_fidelity[1] + 1. Will override ncells and ncharge in simulation and config if provided.

TYPE: tuple DEFAULT: (2, 2)

output_path

base path to save output files, will write to current directory if not specified

TYPE: Optional[str | Path] DEFAULT: None

version

version of HallThruster.jl to use; will search for a global hallthruster_{version} environment in the ~/.julia/environments/ directory. Can also specify a specific git ref (i.e. branch, commit hash, etc.) to use from GitHub. If the hallthruster_{version} environment does not exist, an error will be raised -- you should create this environment first before using it.

TYPE: str DEFAULT: HALLTHRUSTER_VERSION_DEFAULT

pem_to_julia

a dict mapping of PEM shorthand variable names to a list of keys that maps into the HallThruster.jl input/output data structure. Defaults to the provided PEM_TO_JULIA dict defined in hallmd.models.thruster. For example, {'P_b': ['config', 'background_pressure']} will set config['background_pressure'] = P_b. If specified, will override and extend the default mapping.

TYPE: dict | str DEFAULT: 'default'

fidelity_function

a callable that takes a tuple of integers and returns a dictionary of simulation parameters. Defaults to _default_model_fidelity which sets ncells and ncharge based on the input tuple. The returned simulation parameters must be convertable to Julia via the pem_to_julia mapping. The callable should also take in the current json config dict.

TYPE: Callable[[tuple[int, ...]], dict] DEFAULT: 'default'

julia_script

path to a custom Julia script to run. The script should accept the input json file path as a command line argument. Defaults to just calling HallThruster.run_simulation(input_file).

TYPE: Optional[str | Path] DEFAULT: None

run_kwargs

additional keyword arguments to pass to subprocess.run when calling the Julia script. Defaults to check=True.

TYPE: dict DEFAULT: 'default'

shock_threshold

if provided, an error will be raised if the ion velocity reaches a maximum before this threshold axial location (in m) - used to detect and filter unwanted "shock-like" behavior, for example by providing a threshold of half the domain length. If not provided, then no filtering is performed (default).

TYPE: float DEFAULT: None

RETURNS DESCRIPTION
Dataset

dict of Hallthruster.jl outputs: I_B0, I_d, T, eta_c, eta_m, eta_v, and u_ion for ion beam current (A), discharge current (A), thrust (N), current efficiency, mass efficiency, voltage efficiency, and singly-charged ion velocity profile (m/s), all time-averaged.

Source code in src/hallmd/models/thruster.py
def hallthruster_jl(
    thruster_inputs: Dataset = None,
    thruster: Literal["SPT-100"] | str | dict = "SPT-100",
    config: Optional[dict] = None,
    simulation: Optional[dict] = None,
    postprocess: Optional[dict] = None,
    model_fidelity: tuple = (2, 2),
    output_path: Optional[str | Path] = None,
    version: str = HALLTHRUSTER_VERSION_DEFAULT,
    pem_to_julia: dict | str = "default",
    fidelity_function: Callable[[tuple[int, ...]], dict] = "default",
    julia_script: Optional[str | Path] = None,
    run_kwargs: dict = "default",
    shock_threshold: float = None,
) -> Dataset:
    """Run a single `HallThruster.jl` simulation for a given set of inputs. This function will write a temporary
    input file to disk, call `HallThruster.run_simulation()` in Julia, and read the output file back into Python. Will
    return time-averaged performance metrics and ion velocity for use with the PEM.

    Note that the specific inputs and outputs described here can be configured using the `pem_to_julia` dict.

    !!! Warning "Required configuration"
        You must specify a thruster, a domain, a mass flow rate, and a discharge voltage to run the simulation. The
        thruster must be defined in the `hallmd.devices` directory or as a dictionary with the required fields.
        The mass flow rate and discharge voltage are specified in `thruster_inputs` as `mdot_a` (kg/s) and
        `V_a` (V), respectively. The domain is specified as a list `[left_bound, right_bound]` in the
        `config` dictionary. See the
        [HallThruster.jl docs](https://um-pepl.github.io/HallThruster.jl/dev/reference/config/) for more details.

    :param thruster_inputs: named key-value pairs of thruster inputs: `P_b`, `V_a`, `mdot_a`, `T_e`, `u_n`, `l_t`,
                            `a_1`, `a_2`, `delta_z`, `z0`, `p0`, and `V_cc` for background pressure (Torr), anode
                            voltage, anode mass flow rate (kg/s), electron temperature (eV), neutral velocity (m/s),
                            transition length (m), anomalous transport coefficients, and cathode coupling voltage. Will
                            override the corresponding values in `config` if provided.
    :param thruster: the name of the thruster to simulate (must be importable from `hallmd.devices`, see
                     [`load_device`][hallmd.utils.load_device]), or a dictionary that provides geometry and
                     magnetic field information of the thruster to simulate; see the
                     [Hallthruster.jl docs](https://um-pepl.github.io/HallThruster.jl/dev/tutorials/simulation//run/).
                     Will override `thruster` in `config` if provided. If None, will defer to `config`.
                     Defaults to the SPT-100.
    :param config: dictionary of configs for `HallThruster.jl`, see the
                   [Hallthruster.jl docs](https://um-pepl.github.io/HallThruster.jl/dev/reference/config/) for
                   options and formatting.
    :param simulation: dictionary of simulation parameters for `HallThruster.jl`
    :param postprocess: dictionary of post-processing parameters for `Hallthruster.jl`
    :param model_fidelity: tuple of integers that determine the number of cells and the number of charge states to use
                           via `ncells = model_fidelity[0] * 50 + 100` and `ncharge = model_fidelity[1] + 1`.
                           Will override `ncells` and `ncharge` in `simulation` and `config` if provided.
    :param output_path: base path to save output files, will write to current directory if not specified
    :param version: version of HallThruster.jl to use; will
                    search for a global `hallthruster_{version}` environment in the `~/.julia/environments/` directory.
                    Can also specify a specific git ref (i.e. branch, commit hash, etc.) to use from GitHub. If the
                    `hallthruster_{version}` environment does not exist, an error will be raised -- you should create
                    this environment first before using it.
    :param pem_to_julia: a `dict` mapping of PEM shorthand variable names to a list of keys that maps into the
                         `HallThruster.jl` input/output data structure. Defaults to the provided PEM_TO_JULIA dict
                         defined in [`hallmd.models.thruster`][hallmd.models.thruster]. For example,
                         `{'P_b': ['config', 'background_pressure']}` will set `config['background_pressure'] = P_b`.
                         If specified, will override and extend the default mapping.
    :param fidelity_function: a callable that takes a tuple of integers and returns a dictionary of simulation
                              parameters. Defaults to `_default_model_fidelity` which sets `ncells` and `ncharge` based
                              on the input tuple. The returned simulation parameters must be convertable to Julia via
                              the `pem_to_julia` mapping. The callable should also take in the current json config dict.
    :param julia_script: path to a custom Julia script to run. The script should accept the input json file path as
                         a command line argument. Defaults to just calling `HallThruster.run_simulation(input_file)`.
    :param run_kwargs: additional keyword arguments to pass to `subprocess.run` when calling the Julia script.
                       Defaults to `check=True`.
    :param shock_threshold: if provided, an error will be raised if the ion velocity reaches a maximum before this
                            threshold axial location (in m) - used to detect and filter unwanted "shock-like" behavior,
                            for example by providing a threshold of half the domain length. If not provided, then no
                            filtering is performed (default).
    :returns: `dict` of `Hallthruster.jl` outputs: `I_B0`, `I_d`, `T`, `eta_c`, `eta_m`, `eta_v`, and `u_ion` for ion
              beam current (A), discharge current (A), thrust (N), current efficiency, mass efficiency, voltage
              efficiency, and singly-charged ion velocity profile (m/s), all time-averaged.
    """
    if pem_to_julia is None or pem_to_julia == "default":
        pem_to_julia = copy.deepcopy(PEM_TO_JULIA)
    else:
        tmp = copy.deepcopy(PEM_TO_JULIA)
        tmp.update(pem_to_julia)
        pem_to_julia = tmp

    thruster_inputs = thruster_inputs or {}

    # Format PEM inputs for HallThruster.jl
    json_data = _format_hallthruster_jl_input(
        thruster_inputs,
        thruster=thruster,
        config=config,
        simulation=simulation,
        postprocess=postprocess,
        model_fidelity=model_fidelity,
        output_path=output_path,
        pem_to_julia=pem_to_julia,
        fidelity_function=fidelity_function,
    )
    # Get julia environment
    jl_environment = get_jl_env(version) if version is not None else None

    if run_kwargs is None:
        run_kwargs = {}
    elif run_kwargs == "default":
        run_kwargs = {"check": True}

    # Run Julia
    t1 = time.time()
    sim_results = run_hallthruster_jl(json_data, jl_env=jl_environment, jl_script=julia_script, **run_kwargs)
    t2 = time.time()

    # Format QOIs for PEM
    thruster_outputs = _convert_to_pem(sim_results, pem_to_julia)

    # Raise an exception if thrust or beam current are negative (non-physical cases)
    thrust = thruster_outputs.get("T", 0)
    beam_current = thruster_outputs.get("I_B0", 0)
    if thrust < 0 or beam_current < 0:
        raise ValueError(f"Exception due to non-physical case: thrust={thrust} N, beam current={beam_current} A")

    # Raise an exception for ion velocity that exhibits "shock" behavior
    if shock_threshold is not None:
        z_coords = thruster_outputs.get("u_ion_coords")
        ion_velocity = thruster_outputs.get("u_ion")
        if z_coords is not None and ion_velocity is not None:
            if (z_max := z_coords[np.argmax(ion_velocity)]) < shock_threshold:
                raise ValueError(f"Exception due to shock-like behavior: max ion velocity occurs at z={z_max:.3f} m")

    thruster_outputs["model_cost"] = t2 - t1  # seconds

    if output_path is not None:
        output_file = Path(json_data["postprocess"].get("output_file"))
        thruster_outputs["output_path"] = output_file.relative_to(Path(output_path).resolve()).as_posix()

    return thruster_outputs

run_hallthruster_jl(json_input, jl_env=None, jl_script=None, **kwargs)

Python wrapper for HallThruster.run_simulation(json_input) in Julia.

PARAMETER DESCRIPTION
json_input

either a dictionary containing config, simulation, and postprocess options for HallThruster.jl, or a string/Path containing a path to a JSON file with those inputs.

TYPE: dict | str | Path

jl_env

The julia environment containing HallThruster.jl. Defaults to global Julia environment.

TYPE: Optional[str | Path] DEFAULT: None

jl_script

path to a custom Julia script to run. The script should accept the input json file path as a command line argument. Defaults to just calling HallThruster.run_simulation(input_file).

TYPE: Optional[str | Path] DEFAULT: None

kwargs

additional keyword arguments to pass to subprocess.run when calling the Julia script.

DEFAULT: {}

RETURNS DESCRIPTION
dict

dict of Hallthruster.jl outputs. The specific outputs depend on the settings provided in the postprocess dict in the input. If postprocess['output_file'] is present, this function will also write the requested outputs and restart information to that file.

Source code in src/hallmd/models/thruster.py
def run_hallthruster_jl(
    json_input: dict | str | Path, jl_env: Optional[str | Path] = None, jl_script: Optional[str | Path] = None, **kwargs
) -> dict:
    """Python wrapper for `HallThruster.run_simulation(json_input)` in Julia.

    :param json_input: either a dictionary containing `config`, `simulation`, and `postprocess` options for
            HallThruster.jl, or a string/Path containing a path to a JSON file with those inputs.
    :param jl_env: The julia environment containing HallThruster.jl. Defaults to global Julia environment.
    :param jl_script: path to a custom Julia script to run. The script should accept the input json file path as
                      a command line argument. Defaults to just calling `HallThruster.run_simulation(input_file)`.
    :param kwargs: additional keyword arguments to pass to `subprocess.run` when calling the Julia script.

    :returns: `dict` of `Hallthruster.jl` outputs. The specific outputs depend on the settings
              provided in the `postprocess` dict in the input. If `postprocess['output_file']` is present,
              this function will also write the requested outputs and restart information to that file.
    """
    # Read JSON input from file if path provided
    if isinstance(json_input, str | Path):
        with open(json_input, "r") as fp:
            json_input: dict = json.load(fp)

    tempfile_args = dict(suffix=".json", prefix="hallthruster_jl_", mode="w", delete=False, encoding="utf-8")

    # Get output file path. If one not provided, create a temporary
    temp_out = False
    if "output_file" in json_input.get("postprocess", {}):
        output_file = Path(json_input["postprocess"].get("output_file"))
    elif "output_file" in json_input.get("input", {}).get("postprocess", {}):
        output_file = Path(json_input["input"]["postprocess"].get("output_file"))
    else:
        temp_out = True
        fd_out = tempfile.NamedTemporaryFile(**tempfile_args)
        output_file = Path(fd_out.name)
        fd_out.close()

        if json_input.get("input"):
            json_input["input"].setdefault("postprocess", {})
            json_input["input"]["postprocess"]["output_file"] = str(output_file.resolve())
        else:
            json_input.setdefault("postprocess", {})
            json_input["postprocess"]["output_file"] = str(output_file.resolve())

    # Dump input to temporary file
    fd = tempfile.NamedTemporaryFile(**tempfile_args)
    input_file = fd.name
    json.dump(json_input, fd, ensure_ascii=False, indent=4)
    fd.close()

    # Run HallThruster.jl on input file
    if jl_script is None:
        cmd = [
            JL_BINARY,
            "--startup-file=no",
            "-e",
            f'using HallThruster; HallThruster.run_simulation(raw"{input_file}")',
        ]
    else:
        cmd = [JL_BINARY, "--startup-file=no", "--", str(Path(jl_script).resolve()), input_file]

    if jl_env is not None:
        if Path(jl_env).exists():
            cmd.insert(1, f"--project={Path(jl_env).resolve()}")
        else:
            raise ValueError(
                f"Could not find Julia environment {jl_env}. Please create it first. "
                f"See https://github.com/JANUS-Institute/HallThrusterPEM/blob/main/scripts/install_hallthruster.py"
            )

    try:
        subprocess.run(cmd, **kwargs)
    finally:
        # Delete temporary input file
        os.unlink(input_file)

    # Load output data
    with open(output_file, "r") as fp:
        output_data = json.load(fp)

    if temp_out:
        os.unlink(output_file)
        if d := output_data.get("postprocess"):
            if "output_file" in d:
                del d["output_file"]
        if d := output_data.get("input"):
            if d2 := d.get("postprocess"):
                if "output_file" in d2:
                    del d2["output_file"]

    return output_data