charges module
Module contains the abstract base class Charge
and example subclasses.
Charge
is an abstract base class representing a point charge. The class has
abstract methods for the position, velocity, and acceleration of the charge in
the x, y, and z directions.
Subclasses, such as OscillatingCharge
, define their own time-dependent
trajectories. These classes are used by the Simulation
class to perform
electromagnetism calculations and simulations.
All units are in SI.
Charge
Abstract base class used to define a point charge object.
Subclasses implement their unique position, velocity, and acceleration methods in x, y, and z as functions of time. By default, the velocity and acceleration values are determined using finite difference approximations from the position methods; however, they can also be explicitly specified.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
q |
float |
Charge value, can be positive or negative. |
required |
h |
float |
Limit for finite difference calculations for vel and acc. |
required |
solve_time(self, tr, t, x, y, z)
Return equation to solve for the retarded time of the charge.
Source code in pycharge/charges.py
def solve_time(
self, tr: ndarray, t: ndarray, x: ndarray, y: ndarray, z: ndarray
) -> ndarray:
"""Return equation to solve for the retarded time of the charge."""
# Griffiths Eq. 10.55
return ((x-self.xpos(tr))**2 + (y-self.ypos(tr))**2
+ (z-self.zpos(tr))**2)**0.5 - c*(t-tr)
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return x acceleration of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
x acceleration at time values.
"""
return (self.xvel(t+0.5*self.h)-self.xvel(t-0.5*self.h))/self.h
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
x positions at time values. |
Source code in pycharge/charges.py
@abstractmethod
def xpos(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return x position of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
x positions at time values.
"""
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return x velocity of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
x velocity at time values.
"""
return (self.xpos(t+0.5*self.h)-self.xpos(t-0.5*self.h))/self.h
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return y acceleration of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
y acceleration at time values.
"""
return (self.yvel(t+0.5*self.h)-self.yvel(t-0.5*self.h))/self.h
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
y position at time values. |
Source code in pycharge/charges.py
@abstractmethod
def ypos(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return y position of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
y position at time values.
"""
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return y velocity of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
y velocity at time values.
"""
return (self.ypos(t+0.5*self.h)-self.ypos(t-0.5*self.h))/self.h
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return z acceleration of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
z acceleration at time values.
"""
return (self.zvel(t+0.5*self.h)-self.zvel(t-0.5*self.h))/self.h
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
z position at time values. |
Source code in pycharge/charges.py
@abstractmethod
def zpos(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return z position of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
z position at time values.
"""
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
Union[ndarray, float] |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> Union[ndarray, float]:
"""Return z velocity of the charge at specified array of times.
Args:
t: 3D meshgrid array of time values.
Returns:
z velocity at time values.
"""
return (self.zpos(t+0.5*self.h)-self.zpos(t-0.5*self.h))/self.h
LinearAcceleratingCharge
Linearly accelerating charge along the x-axis.
At t<0, the charge is stationary at the origin (x=0, y=0, z=0) and
begins accelerating at t=0. The charge stops accelerating at the time
given by stop_t
and continues moving at a constant velocity.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
acceleration |
float |
Accleration along x-axis. Positive value is in +x direction, negative is in -x direction. |
required |
stop_t |
Optional[float] |
Time when charge stops accelerating. If
|
required |
q |
float |
Charge value, can be positive or negative. Default is |
required |
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
xacc = self.acceleration*np.ones(t.shape)
xacc[t < 0] = 0
xacc[t > self.stop_t] = 0
return xacc
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x positions at time values. |
Source code in pycharge/charges.py
def xpos(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
xpos = 0.5*self.acceleration*t**2
xpos[t < 0] = 0
xpos[t > self.stop_t] = (
0.5*self.acceleration*self.stop_t**2
+ self.acceleration*self.stop_t * (t[t > self.stop_t]-self.stop_t)
)
return xpos
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
xvel = self.acceleration*t
xvel[t < 0] = 0
xvel[t > self.stop_t] = self.acceleration*self.stop_t
return xvel
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> float:
return 0
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y position at time values. |
Source code in pycharge/charges.py
def ypos(self, t: Union[ndarray, float]) -> float:
return 0
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> float:
return 0
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> float:
return 0
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z position at time values. |
Source code in pycharge/charges.py
def zpos(self, t: Union[ndarray, float]) -> float:
return 0
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> float:
return 0
LinearDeceleratingCharge
Linearly decelerating charge along the x-axis.
At t<0, the charge is moving at the value of initial_speed
. At t=0,
the charge is at the origin (x=0, y=0, z=0) and decelerates until the
time given by stop_t
and continues moving at a constant velocity.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
deceleration |
float |
Deceleration along x-axis. Positive value is in +x direction, negative is in -x direction. |
required |
initial_speed |
float |
Initial speed of charge at t<0 along x-axis. |
required |
stop_t |
Optional[float] |
Time when the charge stops decelerating. If
|
required |
q |
float |
Charge value, can be positive or negative. Default is |
required |
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
xacc = np.zeros(t.shape)
xacc[t > 0] = -self.deceleration
xacc[t > self.stop_t] = 0
return xacc
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x positions at time values. |
Source code in pycharge/charges.py
def xpos(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
xpos = self.initial_speed*t
xpos[t > 0] = (self.initial_speed*t[t > 0]
- 0.5*self.deceleration*t[t > 0]**2)
xpos[t > self.stop_t] = (self.initial_speed * self.stop_t
- 0.5*self.deceleration*self.stop_t**2)
return xpos
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
xvel = self.initial_speed*np.ones(t.shape)
xvel[t > 0] = self.initial_speed - self.deceleration*t[t > 0]
xvel[t > self.stop_t] = (self.initial_speed
- self.deceleration*self.stop_t)
return xvel
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> float:
return 0
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y position at time values. |
Source code in pycharge/charges.py
def ypos(self, t: Union[ndarray, float]) -> float:
return 0
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> float:
return 0
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> float:
return 0
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z position at time values. |
Source code in pycharge/charges.py
def zpos(self, t: Union[ndarray, float]) -> float:
return 0
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> float:
return 0
LinearVelocityCharge
Charge moving with constant velocity along the x-axis.
At t=0, the charge is at the position given by init_pos
with the
magnitude of the velocity given by speed
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
speed |
float |
Speed of the charge along the x-axis. Positive value is in +x direction, negative is in -x direction. |
required |
init_pos |
float |
Initial x position at t=0. |
required |
q |
float |
Charge value, can be positive or negative. Default is |
required |
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> float:
return 0
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x positions at time values. |
Source code in pycharge/charges.py
def xpos(self, t: Union[ndarray, float]) -> ndarray:
t = np.array(t)
return self.speed*t + self.init_pos
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> float:
return self.speed
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> float:
return 0
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y position at time values. |
Source code in pycharge/charges.py
def ypos(self, t: Union[ndarray, float]) -> float:
return 0
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> float:
return 0
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> float:
return 0
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z position at time values. |
Source code in pycharge/charges.py
def zpos(self, t: Union[ndarray, float]) -> float:
return 0
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> float:
return 0
OrbittingCharge
Radially orbitting charge in the x-y plane.
At t=0, the charge is at the position (x=radius
, y=0, z=0) and orbits
counter-clockwise.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
radius |
float |
Radius of the orbitting charge trajectory. |
required |
omega |
float |
Angular frequency of the orbit (units: rad/s). |
required |
start_zero |
bool |
Determines if the charge begins orbitting at t=0.
Defaults to |
required |
stop_t |
Optional[float] |
Time when the charge stops orbitting.
If |
required |
q |
float |
Charge value, can be positive or negative. Default is |
required |
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> ndarray:
xacc = -self.radius*self.omega**2*np.cos(self.omega*t)
if self.start_zero:
xacc[t < 0] = 0
if self.stop_t is not None:
xacc[t > self.stop_t] = 0
return xacc
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x positions at time values. |
Source code in pycharge/charges.py
def xpos(self, t: Union[ndarray, float]) -> ndarray:
xpos = self.radius*np.cos(self.omega*t)
if self.start_zero:
xpos[t < 0] = self.radius
if self.stop_t is not None:
xpos[t > self.stop_t] = (self.radius
* np.cos(self.omega*self.stop_t))
return xpos
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> ndarray:
xvel = -self.radius*self.omega*np.sin(self.omega*t)
if self.start_zero:
xvel[t < 0] = 0
if self.stop_t is not None:
xvel[t > self.stop_t] = 0
return xvel
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> ndarray:
yacc = -self.radius*self.omega**2*np.sin(self.omega*t)
if self.start_zero:
yacc[t < 0] = 0
if self.stop_t is not None:
yacc[t > self.stop_t] = 0
return yacc
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
y position at time values. |
Source code in pycharge/charges.py
def ypos(self, t: Union[ndarray, float]) -> ndarray:
ypos = self.radius*np.sin(self.omega*t)
if self.start_zero:
ypos[t < 0] = 0
if self.stop_t is not None:
ypos[t > self.stop_t] = (self.radius
* np.sin(self.omega*self.stop_t))
return ypos
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> ndarray:
yvel = self.radius*self.omega*np.cos(self.omega*t)
if self.start_zero:
yvel[t < 0] = 0
if self.stop_t is not None:
yvel[t > self.stop_t] = 0
return yvel
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> float:
return 0
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z position at time values. |
Source code in pycharge/charges.py
def zpos(self, t: Union[ndarray, float]) -> float:
return 0
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> float:
return 0
OscillatingCharge
Sinusoidally oscillating charge along a specified axis.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
origin |
Tuple[float, float, float] |
List of x, y, and z values for the oscillating charge's origin. |
required |
direction |
Tuple[float, float, float] |
List of x, y, and z values for the charge direction vector. |
required |
amplitude |
float |
Amplitude of the oscillations. |
required |
omega |
float |
Angular frequency of the oscillations (units: rad/s). |
required |
start_zero |
bool |
Determines if the charge begins oscillating at t=0.
Defaults to |
required |
stop_t |
Optional[float] |
Time when the charge stops oscillating.
If |
required |
q |
float |
Charge value, can be positive or negative. Default is |
required |
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> ndarray:
xacc = -(self.direction[0]*self.amplitude
* self.omega**2*np.cos(self.omega*t))
if self.start_zero:
xacc[t < 0] = 0
if self.stop_t is not None:
xacc[t > self.stop_t] = 0
return xacc
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x positions at time values. |
Source code in pycharge/charges.py
def xpos(self, t: Union[ndarray, float]) -> ndarray:
xpos = (self.direction[0]*self.amplitude*np.cos(self.omega*t)
+ self.origin[0])
if self.start_zero:
xpos[t < 0] = self.origin[0]
if self.stop_t is not None:
xpos[t > self.stop_t] = (self.direction[0]*self.amplitude
* np.cos(self.omega*self.stop_t)
+ self.origin[0])
return xpos
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> ndarray:
xvel = -(self.direction[0]*self.amplitude
* self.omega*np.sin(self.omega*t))
if self.start_zero:
xvel[t < 0] = 0
if self.stop_t is not None:
xvel[t > self.stop_t] = 0
return xvel
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> ndarray:
yacc = -(self.direction[1]*self.amplitude
* self.omega**2*np.cos(self.omega*t))
if self.start_zero:
yacc[t < 0] = 0
if self.stop_t is not None:
yacc[t > self.stop_t] = 0
return yacc
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
y position at time values. |
Source code in pycharge/charges.py
def ypos(self, t: Union[ndarray, float]) -> ndarray:
ypos = (self.direction[1]*self.amplitude*np.cos(self.omega*t)
+ self.origin[1])
if self.start_zero:
ypos[t < 0] = self.origin[1]
if self.stop_t is not None:
ypos[t > self.stop_t] = (self.direction[1]*self.amplitude
* np.cos(self.omega*self.stop_t)
+ self.origin[1])
return ypos
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> ndarray:
yvel = -(self.direction[1]*self.amplitude
* self.omega*np.sin(self.omega*t))
if self.start_zero:
yvel[t < 0] = 0
if self.stop_t is not None:
yvel[t > self.stop_t] = 0
return yvel
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> ndarray:
zacc = -(self.direction[2]*self.amplitude
* self.omega**2*np.cos(self.omega*t))
if self.start_zero:
zacc[t < 0] = 0
if self.stop_t is not None:
zacc[t > self.stop_t] = 0
return zacc
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
z position at time values. |
Source code in pycharge/charges.py
def zpos(self, t: Union[ndarray, float]) -> ndarray:
zpos = (self.direction[2]*self.amplitude*np.cos(self.omega*t)
+ self.origin[2])
if self.start_zero:
zpos[t < 0] = self.origin[2]
if self.stop_t is not None:
zpos[t > self.stop_t] = (self.direction[2]*self.amplitude
* np.cos(self.omega*self.stop_t)
+ self.origin[2])
return zpos
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
ndarray |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> ndarray:
zvel = -(self.direction[2]*self.amplitude
* self.omega*np.sin(self.omega*t))
if self.start_zero:
zvel[t < 0] = 0
if self.stop_t is not None:
zvel[t > self.stop_t] = 0
return zvel
StationaryCharge
Stationary point charge.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
position |
Tuple[float, float, float] |
List of x, y, and z values for the stationary position. |
required |
q |
float |
Charge value, can be positive or negative. Default is |
required |
xacc(self, t)
Return x acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
x acceleration at time values. |
Source code in pycharge/charges.py
def xacc(self, t: Union[ndarray, float]) -> float:
return 0
xpos(self, t)
Return x position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
x positions at time values. |
Source code in pycharge/charges.py
def xpos(self, t: Union[ndarray, float]) -> float:
return self.position[0]
xvel(self, t)
Return x velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
x velocity at time values. |
Source code in pycharge/charges.py
def xvel(self, t: Union[ndarray, float]) -> float:
return 0
yacc(self, t)
Return y acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y acceleration at time values. |
Source code in pycharge/charges.py
def yacc(self, t: Union[ndarray, float]) -> float:
return 0
ypos(self, t)
Return y position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y position at time values. |
Source code in pycharge/charges.py
def ypos(self, t: Union[ndarray, float]) -> float:
return self.position[1]
yvel(self, t)
Return y velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
y velocity at time values. |
Source code in pycharge/charges.py
def yvel(self, t: Union[ndarray, float]) -> float:
return 0
zacc(self, t)
Return z acceleration of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z acceleration at time values. |
Source code in pycharge/charges.py
def zacc(self, t: Union[ndarray, float]) -> float:
return 0
zpos(self, t)
Return z position of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z position at time values. |
Source code in pycharge/charges.py
def zpos(self, t: Union[ndarray, float]) -> float:
return self.position[2]
zvel(self, t)
Return z velocity of the charge at specified array of times.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t |
Union[ndarray, float] |
3D meshgrid array of time values. |
required |
Returns:
Type | Description |
---|---|
float |
z velocity at time values. |
Source code in pycharge/charges.py
def zvel(self, t: Union[ndarray, float]) -> float:
return 0