FiniteStraightWire

class plasmapy.formulary.magnetostatics.FiniteStraightWire(p1: Quantity, p2: Quantity, current: Quantity)[source]

Bases: Wire

Finite length straight wire class.

The p1 to p2 direction is the positive current direction.

Parameters:
  • p1 (Quantity) – Three-dimensional Cartesian coordinate of one end of the straight wire.

  • p2 (Quantity) – Three-dimensional Cartesian coordinate of another end of the straight wire.

  • current (astropy.units.Quantity) – Electric current.

Methods Summary

magnetic_field(p)

Calculate magnetic field generated by this wire at position p.

to_GeneralWire()

Convert this Wire into a GeneralWire.

Methods Documentation

magnetic_field(p) Quantity[source]

Calculate magnetic field generated by this wire at position p.

Parameters:

p (astropy.units.Quantity) – Three-dimensional position vector

Returns:

B – Magnetic field at the specified position

Return type:

astropy.units.Quantity

Notes

The magnetic field generated by a straight, finite wire with constant electric current can be found at a point in 3D space using the Biot–Savart law.

Let the point where the magnetic field will be calculated be represented by the point \(p_0\) (or p) and the wire’s beginning and end as \(p_1\) and \(p_2\), respectively (with corresponding position vectors \(\vec{p}_0\), \(\vec{p}_1\), and \(\vec{p}_2\), respectively). Further, the vector from points \(p_i\) to \(p_j\) can be written as \(\vec{p}_{ij} = \vec{p}_j - \vec{p}_i\).

Next, consider the triangle with the points \(p_0\), \(p_1\), and \(p_2\) as vertices. The vector from the vertex \(p_0\) to the perpendicular foot opposite the vertex \(p_0\), which will be used to find the unit vector in the direction of the magnetic field, can be expressed as

\[\vec{p}_f = \vec{p}_1 + \vec{p}_{12} \frac{\vec{p}_{10} \cdot \vec{p}_{12}} {|\vec{p}_{12}|^2}.\]

The magnetic field \(\vec{B}\) generated by the wire with current \(I\) can be found at the point \(p_0\) using the Biot–Savart law which in this case simplifies to

\[\vec{B} = \frac{\mu_0 I}{4π} (\cos θ_1 - \cos θ_2) \hat{B}\]

where \(\mu_0\) is the permeability of free space, \(\theta_1\) (\(\theta_2\)) is the angle between \(\vec{p}_{10}\) (\(\vec{p}_{20}\)) and \(\vec{p}_{12}\) with

\[\cos\theta_1 = \frac{\vec{p}_{10} \cdot \vec{p}_{12}} {|\vec{p}_{10}| |\vec{p}_{12}|}, \quad \cos\theta_2 = \frac{\vec{p}_{20} \cdot \vec{p}_{12}} {|\vec{p}_{20}| |\vec{p}_{12}|},\]

and

\[\hat{B} = \frac{\vec{p}_{12} \times \vec{p}_{f0}} {|\vec{p}_{12} \times \vec{p}_{f0}|}\]

is the unit vector in the direction of the magnetic field at the point \(p_0\).

to_GeneralWire()[source]

Convert this Wire into a GeneralWire.