Anomalous Transport
HallThruster has a few anomalous transport models built in and allows users to define their own. This page describes these models and the process by which algebraic and multi-equation transport models can be added by the user.
The AnomalousTransportModel interface is not yet finalized and subject to revision. Keep this in mind when using this feature.
Built-in Models
NoAnom()
Model for no anomalous transport (anomalous collision frequency = 0).
anom_model = NoAnom()
Bohm(c)
Model where the anomalous collision frequency scales with the electron cyclotron frequency ωce times some scaling factor c
TwoZoneBohm(c1, c2)
HallThruster's default anomalous transport option. This is a standard model of anomalous transport frequently used in Hall thruster simulations. The anomalous collision frequency is defined as
\[\begin{aligned} +
Anomalous Transport
HallThruster has a few anomalous transport models built in and allows users to define their own. This page describes these models and the process by which algebraic and multi-equation transport models can be added by the user.
The AnomalousTransportModel interface is not yet finalized and subject to revision. Keep this in mind when using this feature.
Built-in Models
!!! warning Incomplete documentation Some models (those relating to pressure-dependent effects) have not been finalized and are not represented in this documentation. Please see the source code for a complete listing.
NoAnom()
Model for no anomalous transport (anomalous collision frequency = 0).
anom_model = NoAnom()
Bohm(c)
Model where the anomalous collision frequency scales with the electron cyclotron frequency ωce times some scaling factor c
TwoZoneBohm(c1, c2)
HallThruster's default anomalous transport option. This is a standard model of anomalous transport frequently used in Hall thruster simulations. The anomalous collision frequency is defined as
\[\begin{aligned} \nu_{AN} &= c_1 \omega_{ce} \quad z < L_{ch} \\ &= c_2 \omega_{ce} \quad z > L_{ch} \end{aligned}\]
In the above expression, $c_1$ and $c_2$ are tunable coefficients, $\omega_{ce} = e B / m_e$ is the electron cyclotron frequency, and $L_{ch}$ is the channel length. A TwoZoneBohm
model is initialized as follows
anom_model = TwoZoneBohm(c1, c2)
The transition between the zones is determined by the user-provided transition function. This defaults to a step function.
MultiLogBohm(z, c)
Model similar to that employed in Hall2De, where the mobility is Bohm-like (i.e. νan(z) = c(z) * ωce(z)
) and z is in meters.
The function c(z)
is defined by a sequence of nodes (z, c)
provided by the user. At z = z[1]
, c(z) = c[1]
, and so forth.
At z[i] < z < z[i+1]
, log(c)
is defined by linearly interpolating between log(c[i])
and log(c[i+1])
.
For z < z[1]
, c = c[1]
and for z > z[end]
, c(z) = c[end]
.
The user may also provide a single array of [z[1], z[2], ..., z[end], c[1], c[2], ..., c[end]]
. The number of z
values must be equal to the number of c values.
The AnomalousTransportModel
interface
Currently, HallThruster.jl expects all models to be written to be callable structs, taking arguments U, params, i
, where U
is the system state vector, params
are the simulation parameters (including the cache of all variables), and i
is the index of the cell.
Custom anomalous transport models
Users of HallThruster may define their own models by defining a custom subtype of AnomalousTransportModel
. Suppose we want to implement $nu_{AN} = \beta \omega_{ce}$ (classic Bohm diffusion). This is a fixed anomalous transport model and does not change as the simulation progresses. We would first define our type:
using HallThruster
@@ -23,7 +23,7 @@
\end{aligned}\]In this expression, $K$ is a tunable coefficient, $u_i$ is the ion velocity, $W = n_e k_B T_e$ is the wave energy density, $m_e$ is the electron mass, $n_e$ is the electron number density, $c_s$ is the ion sound speed, and $v_{de}$ is the electron azimuthal drift speed. Let's say we want to save the wave energy density in addition to the anomalous collision frequency. We begin by defining the model:
struct LafleurModel <: AnomalousTransportModel
K::Float64
end
Next, we add a method to the num_anom_variables
function. Since we want to save the wave energy density, we need 1 additional anomalous transport variable.
num_anom_variables(::LafleurModel) = 1
Now, we define the behavior of the model in a function. Since the model is based on assuming the wave energy convects with the ions, we will use upwind differencing for the gradient.
function (model::LafleurModel)(νan, params)
-
+
(;config, cache) = params
mi = config.propellant.m
K = model.K
@@ -53,7 +53,7 @@
# Save W to cache.anom_variables[1]
anom_variables[1][i] = W
-
+
# Return anomalous collision frequency
νan[i] = abs(K * grad_ui_W / (me * cs * vde * ne))
end
@@ -63,7 +63,7 @@
anom_variables[1][end] = anom_variables[1][end-1]
νan[1] = νan[2]
νan[end] = νan[end-1]
-
+
return νan
end
The saved value of the wave energy density can then be recovered as
solution.savevals[frame].anom_variables[1]
where solution
is the Solution
object resulting from a call to run_simulation
Solving arbitrary PDEs using the AnomalousTransportModel interface
We can use this basic interface to solve PDEs within HallThruster.jl. We demonstrate this by solving the scalar advection equation using first-order upwind differencing in space and forward Euler integration in time, with periodic boundary conditions. The scalar advection equation is given by:
\[ \begin{aligned}
\frac{\partial u}{\partial t} + a \frac{\partial u}{\partial x} = 0
@@ -78,17 +78,17 @@
# 1) Advected quantity u
# 2) gradient of advected quantity (du/dz)
HallThruster.num_anom_variables(::ScalarAdvection) = 2
Next, we define the model function:
function (model::ScalarAdvection)(νan, params)
-
+
ncells = length(νan)
-
- # Extract variables from params
+
+ # Extract variables from params
cache = params.cache
z = params.z_cell
u = cache.anom_variables[1]
du_dz = cache.anom_variables[2]
a = model.advection_velocity
dt = params.dt
-
+
if params.iteration[] < 1
# Initialize
model.initializer(u, z)
@@ -151,9 +151,9 @@
# Run simulation
solution = HallThruster.run_simulation(
config;
- ncells,
+ ncells,
dt,
- duration = nsteps * dt,
+ duration = nsteps * dt,
nsave = nsteps
)
@@ -164,7 +164,7 @@
# Time needed to transit the domain
t_transit = L / advection_velocity
-# Number of periods
+# Number of periods
num_periods = floor(Int, nsteps * dt / t_transit)
# Plot results
@@ -189,17 +189,17 @@
HallThruster.num_anom_variables(::ScalarAdvection) = 2
function (model::ScalarAdvection)(νan, params)
-
+
ncells = length(νan)
-
- # Extract variables from params
+
+ # Extract variables from params
cache = params.cache
z = params.z_cell
u = cache.anom_variables[1]
du_dz = cache.anom_variables[2]
a = model.advection_velocity
dt = params.dt
-
+
if params.iteration[] < 1
# Initialize
model.initializer(u, z)
@@ -277,9 +277,9 @@
# Run simulation
solution = HallThruster.run_simulation(
config;
- ncells,
+ ncells,
dt,
- duration = nsteps * dt,
+ duration = nsteps * dt,
nsave = nsteps
)
@@ -292,7 +292,7 @@
# Time needed to transit the domain
t_transit = L / advection_velocity
-# Number of periods
+# Number of periods
num_periods = floor(Int, nsteps * dt / t_transit)
# Plot results
@@ -304,4 +304,4 @@
linecolor = cgrad(:turbo, num_periods + 1, categorical = true)[i+1]
)
end
-display(p)
Settings
This document was generated with Documenter.jl version 0.27.25 on Thursday 14 March 2024. Using Julia version 1.10.2.