mass_density
- plasmapy.formulary.densities.mass_density(
- density: Unit('1 / m3'), Unit('kg / m3'),
- particle: str | int | integer | Particle | CustomParticle | Quantity,
- z_ratio: float | None = 1,
Calculate the mass density from a number density.
\[ρ = \left| \frac{Z_s}{Z_{particle}} \right| n_s m_{particle} = | Z_{ratio} | n_s m_{particle}\]where \(m_{particle}\) is the particle mass, \(n_s\) is a number density for plasma species \(s\), \(Z_s\) is the charge number of species \(s\), and \(Z_{particle}\) is the charge number of
particle. For example, if the electron density is given for \(n_s\) andparticleis a doubly ionized atom, then \(Z_{ratio} = -1 / 2\).- Parameters:
density (
Quantity) – Either a particle number density (in units of m-3 or equivalent) or a mass density (in units of kg / m3 or equivalent). Ifdensityis a mass density, then it will be passed through and returned without modification.particle (
Particle) – The particle for which the mass density is being calculated. Must be aParticleor a value convertible to aParticle(e.g.,'p+'for protons,'D+'for deuterium, or'He-4 +1'for singly ionized helium-4).z_ratio (
int|float, default:1) – The ratio of the charge numbers corresponding to the plasma species represented bydensityand theparticle. For example, if the givendensityis an electron density andparticleis doubly ionizedHe, thenz_ratio = -0.5.
- Raises:
UnitTypeError – If the
densitydoes not have units equivalent to a number density or mass density.TypeError – If
densityis not of typeQuantity, or convertible.TypeError – If
particleis not of type or convertible toParticle.ValueError – If
densityis negative.
- Returns:
The mass density for the plasma species represented by
particle.- Return type:
Examples
>>> import astropy.units as u >>> mass_density(1 * u.m**-3, "p+") <Quantity 1.67262...e-27 kg / m3> >>> mass_density(4 * u.m**-3, "D+") <Quantity 1.33743...e-26 kg / m3> >>> mass_density(2.0e12 * u.cm**-3, "He") <Quantity 1.32929...e-08 kg / m3> >>> mass_density(2.0e12 * u.cm**-3, "He", z_ratio=0.5) <Quantity 6.64647...e-09 kg / m3> >>> mass_density(1.0 * u.g * u.m**-3, "") <Quantity 0.001 kg / m3>