magnetic_energy_density¶
-
plasmapy.formulary.parameters.
magnetic_energy_density
(B: Unit("T")) -> Unit("J / m3")¶ Calculate the magnetic energy density.
Aliases:
ub_
Parameters: B (Quantity) – The magnetic field in units convertible to tesla.
Returns: E_B – The magnetic energy density in units of joules per cubic meter.
Return type: Raises: TypeError
– If the input is not a Quantity.UnitConversionError
– If the input is not in units convertible to tesla.ValueError
– If the magnetic field strength does not have an appropriate. value.
Warns: ~astropy.units.UnitsWarning – If units are not provided, SI units are assumed
Notes
The magnetic energy density is given by:
\[E_B = \frac{B^2}{2 \mu_0}\]The motivation behind having two separate functions for magnetic pressure and magnetic energy density is that it allows greater insight into the physics that are being considered by the user and thus more readable code.
See also
magnetic_pressure()
- Returns an equivalent Quantity, except in units of pascals.
Example
>>> from astropy import units as u >>> magnetic_energy_density(0.1*u.T) <Quantity 3978.87... J / m3>