HarrisSheet
- class plasmapy.plasma.equilibria1d.HarrisSheet(B0, delta, P0=<Quantity 0. Pa>)[source]
Bases:
object
Define a Harris Sheet Equilibrium.
Magnetic field will be in the \(±x\) direction and the current density will be in the \(±z\) direction in a \(\hat{x} × \hat{y} = \hat{z}\) coordinate system.
- Parameters:
Notes
A current sheet is current limited to a surface.
A Harris sheet is a 1D ideal MHD equilibrium. In resistive MHD if there is any resistivity, it won’t be a true equilibrium since the resistivity will gradually smooth the profile out over time.
A Harris sheet is often used as the initial condition for simulations of magnetic reconnection.
Examples
>>> import astropy.units as u >>> harris_sheet = HarrisSheet(delta = 3 * u.m, B0 = 2 * u.T) >>> harris_sheet.magnetic_field(y = 5 * u.m) <Quantity 1.8622... T>
Methods Summary
Compute the current density.
Compute the magnetic field.
Compute plasma pressure.
Methods Documentation
- current_density(y: Unit('m'))[source]
Compute the current density.
\[J_z(y) = - \frac{B_0}{δ μ_0) \mathrm{sech}^2 \left( \frac{y}{δ} \right)\]- Parameters:
y (
Quantity
) – Orthogonal distance from the current sheet center.
- magnetic_field(y: Unit('m'))[source]
Compute the magnetic field.
In this equation, \(B_0\) is the asymptotic magnitude of the magnetic field for \(y → ±∞\) and \(δ\) is the thickness of the sheet.
\[B_x(y) = B_0 \tanh \left( \frac{y}{δ} \right)\]- Parameters:
y (
Quantity
) – Orthogonal distance from the current sheet center.