Copyright | (c) Scott N. Walck 2014-2018 |
---|---|
License | BSD3 (see LICENSE) |
Maintainer | Scott N. Walck <walck@lvc.edu> |
Stability | experimental |
Safe Haskell | Trustworthy |
Language | Haskell98 |
Functions for learning physics.
- type TheTime = Double
- type TimeStep = Double
- type Velocity = Vec
- type SimpleState = (TheTime, Position, Velocity)
- type SimpleAccelerationFunction = SimpleState -> Vec
- simpleStateDeriv :: SimpleAccelerationFunction -> DifferentialEquation SimpleState
- simpleRungeKuttaStep :: SimpleAccelerationFunction -> TimeStep -> SimpleState -> SimpleState
- data St = St {}
- data DSt = DSt Vec Vec
- type OneParticleSystemState = (TheTime, St)
- type OneParticleAccelerationFunction = OneParticleSystemState -> Vec
- oneParticleStateDeriv :: OneParticleAccelerationFunction -> DifferentialEquation OneParticleSystemState
- oneParticleRungeKuttaStep :: OneParticleAccelerationFunction -> TimeStep -> OneParticleSystemState -> OneParticleSystemState
- oneParticleRungeKuttaSolution :: OneParticleAccelerationFunction -> TimeStep -> OneParticleSystemState -> [OneParticleSystemState]
- type TwoParticleSystemState = (TheTime, St, St)
- type TwoParticleAccelerationFunction = TwoParticleSystemState -> (Vec, Vec)
- twoParticleStateDeriv :: TwoParticleAccelerationFunction -> DifferentialEquation TwoParticleSystemState
- twoParticleRungeKuttaStep :: TwoParticleAccelerationFunction -> TimeStep -> TwoParticleSystemState -> TwoParticleSystemState
- type ManyParticleSystemState = (TheTime, [St])
- type ManyParticleAccelerationFunction = ManyParticleSystemState -> [Vec]
- manyParticleStateDeriv :: ManyParticleAccelerationFunction -> DifferentialEquation ManyParticleSystemState
- manyParticleRungeKuttaStep :: ManyParticleAccelerationFunction -> TimeStep -> ManyParticleSystemState -> ManyParticleSystemState
- type Charge = Double
- data ChargeDistribution
- totalCharge :: ChargeDistribution -> Charge
- type Current = Double
- data CurrentDistribution
- eField :: ChargeDistribution -> VectorField
- electricFlux :: Surface -> ChargeDistribution -> Double
- electricPotentialFromField :: Position -> VectorField -> ScalarField
- electricPotentialFromCharge :: ChargeDistribution -> ScalarField
- bField :: CurrentDistribution -> VectorField
- magneticFlux :: Surface -> CurrentDistribution -> Double
- data Vec
- xComp :: Vec -> Double
- yComp :: Vec -> Double
- zComp :: Vec -> Double
- vec :: Double -> Double -> Double -> Vec
- (^+^) :: AdditiveGroup v => v -> v -> v
- (^-^) :: AdditiveGroup v => v -> v -> v
- (*^) :: VectorSpace v => Scalar v -> v -> v
- (^*) :: (VectorSpace v, (~) * s (Scalar v)) => v -> s -> v
- (^/) :: (VectorSpace v, (~) * s (Scalar v), Fractional s) => v -> s -> v
- (<.>) :: InnerSpace v => v -> v -> Scalar v
- (><) :: Vec -> Vec -> Vec
- magnitude :: (InnerSpace v, (~) * s (Scalar v), Floating s) => v -> s
- zeroV :: AdditiveGroup v => v
- negateV :: AdditiveGroup v => v -> v
- sumV :: (Foldable f, AdditiveGroup v) => f v -> v
- iHat :: Vec
- jHat :: Vec
- kHat :: Vec
- data Position
- type Displacement = Vec
- type ScalarField = Position -> Double
- type VectorField = Position -> Vec
- type Field v = Position -> v
- type CoordinateSystem = (Double, Double, Double) -> Position
- cartesian :: CoordinateSystem
- cylindrical :: CoordinateSystem
- spherical :: CoordinateSystem
- cart :: Double -> Double -> Double -> Position
- cyl :: Double -> Double -> Double -> Position
- sph :: Double -> Double -> Double -> Position
- cartesianCoordinates :: Position -> (Double, Double, Double)
- cylindricalCoordinates :: Position -> (Double, Double, Double)
- sphericalCoordinates :: Position -> (Double, Double, Double)
- displacement :: Position -> Position -> Displacement
- shiftPosition :: Displacement -> Position -> Position
- shiftObject :: Displacement -> (a -> Position) -> a -> Position
- shiftField :: Displacement -> (Position -> v) -> Position -> v
- addFields :: AdditiveGroup v => [Field v] -> Field v
- rHat :: VectorField
- thetaHat :: VectorField
- phiHat :: VectorField
- sHat :: VectorField
- xHat :: VectorField
- yHat :: VectorField
- zHat :: VectorField
- data Curve = Curve {}
- normalizeCurve :: Curve -> Curve
- concatCurves :: Curve -> Curve -> Curve
- concatenateCurves :: [Curve] -> Curve
- reverseCurve :: Curve -> Curve
- evalCurve :: Curve -> Double -> Position
- shiftCurve :: Displacement -> Curve -> Curve
- straightLine :: Position -> Position -> Curve
- simpleLineIntegral :: (InnerSpace v, Scalar v ~ Double) => Int -> Field v -> Curve -> v
- dottedLineIntegral :: Int -> VectorField -> Curve -> Double
- crossedLineIntegral :: Int -> VectorField -> Curve -> Vec
- data Surface = Surface {
- surfaceFunc :: (Double, Double) -> Position
- lowerLimit :: Double
- upperLimit :: Double
- lowerCurve :: Double -> Double
- upperCurve :: Double -> Double
- unitSphere :: Surface
- centeredSphere :: Double -> Surface
- sphere :: Double -> Position -> Surface
- northernHemisphere :: Surface
- disk :: Double -> Surface
- shiftSurface :: Displacement -> Surface -> Surface
- surfaceIntegral :: (VectorSpace v, Scalar v ~ Double) => Int -> Int -> Field v -> Surface -> v
- dottedSurfaceIntegral :: Int -> Int -> VectorField -> Surface -> Double
- data Volume = Volume {}
- unitBall :: Volume
- unitBallCartesian :: Volume
- centeredBall :: Double -> Volume
- ball :: Double -> Position -> Volume
- northernHalfBall :: Volume
- centeredCylinder :: Double -> Double -> Volume
- shiftVolume :: Displacement -> Volume -> Volume
- volumeIntegral :: (VectorSpace v, Scalar v ~ Double) => Int -> Int -> Int -> Field v -> Volume -> v
- class (VectorSpace (Diff p), Fractional (Scalar (Diff p))) => StateSpace p where
- type Diff p
- (.-^) :: StateSpace p => p -> Diff p -> p
- type Time p = Scalar (Diff p)
- type DifferentialEquation state = state -> Diff state
- type InitialValueProblem state = (DifferentialEquation state, state)
- type EvolutionMethod state = DifferentialEquation state -> Time state -> state -> state
- type SolutionMethod state = InitialValueProblem state -> [state]
- stepSolution :: EvolutionMethod state -> Time state -> SolutionMethod state
- eulerMethod :: StateSpace state => EvolutionMethod state
- rungeKutta4 :: StateSpace p => (p -> Diff p) -> Time p -> p -> p
- integrateSystem :: StateSpace p => (p -> Diff p) -> Time p -> p -> [p]
- label :: String -> (Double, Double) -> Attribute
- postscript :: Attribute
- psFile :: FilePath -> Attribute
- polarToCart :: (Float, Float) -> (Float, Float)
- cartToPolar :: (Float, Float) -> (Float, Float)
- arrow :: Point -> Point -> Picture
- thickArrow :: Float -> Point -> Point -> Picture
Mechanics
Simple one-particle state
type SimpleState = (TheTime, Position, Velocity) Source #
A simple one-particle state, to get started quickly with mechanics of one particle.
type SimpleAccelerationFunction = SimpleState -> Vec Source #
An acceleration function gives the particle's acceleration as a function of the particle's state. The specification of this function is what makes one single-particle mechanics problem different from another. In order to write this function, add all of the forces that act on the particle, and divide this net force by the particle's mass. (Newton's second law).
:: SimpleAccelerationFunction | acceleration function for the particle |
-> DifferentialEquation SimpleState | differential equation |
Time derivative of state for a single particle with a constant mass.
:: SimpleAccelerationFunction | acceleration function for the particle |
-> TimeStep | time step |
-> SimpleState | initial state |
-> SimpleState | state after one time step |
Single Runge-Kutta step
One-particle state
The state of a single particle is given by the position of the particle and the velocity of the particle.
The associated vector space for the state of a single particle.
type OneParticleSystemState = (TheTime, St) Source #
The state of a system of one particle is given by the current time, the position of the particle, and the velocity of the particle. Including time in the state like this allows us to have time-dependent forces.
type OneParticleAccelerationFunction = OneParticleSystemState -> Vec Source #
An acceleration function gives the particle's acceleration as a function of the particle's state.
oneParticleStateDeriv Source #
:: OneParticleAccelerationFunction | acceleration function for the particle |
-> DifferentialEquation OneParticleSystemState | differential equation |
Time derivative of state for a single particle with a constant mass.
oneParticleRungeKuttaStep Source #
:: OneParticleAccelerationFunction | acceleration function for the particle |
-> TimeStep | time step |
-> OneParticleSystemState | initial state |
-> OneParticleSystemState | state after one time step |
Single Runge-Kutta step
oneParticleRungeKuttaSolution Source #
:: OneParticleAccelerationFunction | acceleration function for the particle |
-> TimeStep | time step |
-> OneParticleSystemState | initial state |
-> [OneParticleSystemState] | state after one time step |
List of system states
Two-particle state
type TwoParticleSystemState = (TheTime, St, St) Source #
The state of a system of two particles is given by the current time, the position and velocity of particle 1, and the position and velocity of particle 2.
type TwoParticleAccelerationFunction = TwoParticleSystemState -> (Vec, Vec) Source #
An acceleration function gives a pair of accelerations (one for particle 1, one for particle 2) as a function of the system's state.
twoParticleStateDeriv Source #
:: TwoParticleAccelerationFunction | acceleration function for two particles |
-> DifferentialEquation TwoParticleSystemState | differential equation |
Time derivative of state for two particles with constant mass.
twoParticleRungeKuttaStep Source #
:: TwoParticleAccelerationFunction | acceleration function |
-> TimeStep | time step |
-> TwoParticleSystemState | initial state |
-> TwoParticleSystemState | state after one time step |
Single Runge-Kutta step for two-particle system
Many-particle state
type ManyParticleSystemState = (TheTime, [St]) Source #
The state of a system of many particles is given by the current time and a list of one-particle states.
type ManyParticleAccelerationFunction = ManyParticleSystemState -> [Vec] Source #
An acceleration function gives a list of accelerations (one for each particle) as a function of the system's state.
manyParticleStateDeriv Source #
:: ManyParticleAccelerationFunction | acceleration function for many particles |
-> DifferentialEquation ManyParticleSystemState | differential equation |
Time derivative of state for many particles with constant mass.
manyParticleRungeKuttaStep Source #
:: ManyParticleAccelerationFunction | acceleration function |
-> TimeStep | time step |
-> ManyParticleSystemState | initial state |
-> ManyParticleSystemState | state after one time step |
Single Runge-Kutta step for many-particle system
E&M
Charge
data ChargeDistribution Source #
A charge distribution is a point charge, a line charge, a surface charge,
a volume charge, or a combination of these.
The ScalarField
describes a linear charge density, a surface charge density,
or a volume charge density.
PointCharge Charge Position | point charge |
LineCharge ScalarField Curve |
|
SurfaceCharge ScalarField Surface |
|
VolumeCharge ScalarField Volume |
|
MultipleCharges [ChargeDistribution] | combination of charge distributions |
totalCharge :: ChargeDistribution -> Charge Source #
Total charge (in C) of a charge distribution.
Current
data CurrentDistribution Source #
A current distribution is a line current (current through a wire), a surface current,
a volume current, or a combination of these.
The VectorField
describes a surface current density
or a volume current density.
LineCurrent Current Curve | current through a wire |
SurfaceCurrent VectorField Surface |
|
VolumeCurrent VectorField Volume |
|
MultipleCurrents [CurrentDistribution] | combination of current distributions |
Electric Field
eField :: ChargeDistribution -> VectorField Source #
The electric field produced by a charge distribution. This is the simplest way to find the electric field, because it works for any charge distribution (point, line, surface, volume, or combination).
Electric Flux
electricFlux :: Surface -> ChargeDistribution -> Double Source #
The electric flux through a surface produced by a charge distribution.
Electric Potential
electricPotentialFromField Source #
:: Position | position where electric potential is zero |
-> VectorField | electric field |
-> ScalarField | electric potential |
Electric potential from electric field, given a position to be the zero of electric potential.
electricPotentialFromCharge :: ChargeDistribution -> ScalarField Source #
Electric potential produced by a charge distribution. The position where the electric potential is zero is taken to be infinity.
Magnetic Field
bField :: CurrentDistribution -> VectorField Source #
The magnetic field produced by a current distribution. This is the simplest way to find the magnetic field, because it works for any current distribution (line, surface, volume, or combination).
Magnetic Flux
magneticFlux :: Surface -> CurrentDistribution -> Double Source #
The magnetic flux through a surface produced by a current distribution.
Geometry
Vectors
A type for vectors.
Form a vector by giving its x, y, and z components.
(^+^) :: AdditiveGroup v => v -> v -> v infixl 6 #
Add vectors
(^-^) :: AdditiveGroup v => v -> v -> v infixl 6 #
Group subtraction
(*^) :: VectorSpace v => Scalar v -> v -> v infixr 7 #
Scale a vector
(^*) :: (VectorSpace v, (~) * s (Scalar v)) => v -> s -> v infixl 7 #
Vector multiplied by scalar
(^/) :: (VectorSpace v, (~) * s (Scalar v), Fractional s) => v -> s -> v infixr 7 #
Vector divided by scalar
(<.>) :: InnerSpace v => v -> v -> Scalar v infixr 7 #
Inner/dot product
magnitude :: (InnerSpace v, (~) * s (Scalar v), Floating s) => v -> s #
Length of a vector. See also magnitudeSq
.
zeroV :: AdditiveGroup v => v #
The zero element: identity for '(^+^)'
negateV :: AdditiveGroup v => v -> v #
Additive inverse
sumV :: (Foldable f, AdditiveGroup v) => f v -> v #
Sum over several vectors
Position
A type for position. Position is not a vector because it makes no sense to add positions.
type Displacement = Vec Source #
A displacement is a vector.
type ScalarField = Position -> Double Source #
A scalar field associates a number with each position in space.
type VectorField = Position -> Vec Source #
A vector field associates a vector with each position in space.
type Field v = Position -> v Source #
Sometimes we want to be able to talk about a field without saying whether it is a scalar field or a vector field.
type CoordinateSystem = (Double, Double, Double) -> Position Source #
A coordinate system is a function from three parameters to space.
cartesian :: CoordinateSystem Source #
The Cartesian coordinate system. Coordinates are (x,y,z).
cylindrical :: CoordinateSystem Source #
The cylindrical coordinate system. Coordinates are (s,phi,z), where s is the distance from the z axis and phi is the angle with the x axis.
spherical :: CoordinateSystem Source #
The spherical coordinate system. Coordinates are (r,theta,phi), where r is the distance from the origin, theta is the angle with the z axis, and phi is the azimuthal angle.
A helping function to take three numbers x, y, and z and form the appropriate position using Cartesian coordinates.
A helping function to take three numbers s, phi, and z and form the appropriate position using cylindrical coordinates.
A helping function to take three numbers r, theta, and phi and form the appropriate position using spherical coordinates.
cartesianCoordinates :: Position -> (Double, Double, Double) Source #
Returns the three Cartesian coordinates as a triple from a position.
cylindricalCoordinates :: Position -> (Double, Double, Double) Source #
Returns the three cylindrical coordinates as a triple from a position.
sphericalCoordinates :: Position -> (Double, Double, Double) Source #
Returns the three spherical coordinates as a triple from a position.
:: Position | source position |
-> Position | target position |
-> Displacement |
Displacement from source position to target position.
shiftPosition :: Displacement -> Position -> Position Source #
Shift a position by a displacement.
shiftObject :: Displacement -> (a -> Position) -> a -> Position Source #
An object is a map into Position
.
shiftField :: Displacement -> (Position -> v) -> Position -> v Source #
A field is a map from Position
.
addFields :: AdditiveGroup v => [Field v] -> Field v Source #
Add two scalar fields or two vector fields.
rHat :: VectorField Source #
The vector field in which each point in space is associated
with a unit vector in the direction of increasing spherical coordinate
r, while spherical coordinates theta and phi
are held constant.
Defined everywhere except at the origin.
The unit vector rHat
points in different directions at different points
in space. It is therefore better interpreted as a vector field, rather
than a vector.
thetaHat :: VectorField Source #
The vector field in which each point in space is associated with a unit vector in the direction of increasing spherical coordinate theta, while spherical coordinates r and phi are held constant. Defined everywhere except on the z axis.
phiHat :: VectorField Source #
The vector field in which each point in space is associated with a unit vector in the direction of increasing (cylindrical or spherical) coordinate phi, while cylindrical coordinates s and z (or spherical coordinates r and theta) are held constant. Defined everywhere except on the z axis.
sHat :: VectorField Source #
The vector field in which each point in space is associated with a unit vector in the direction of increasing cylindrical coordinate s, while cylindrical coordinates phi and z are held constant. Defined everywhere except on the z axis.
xHat :: VectorField Source #
The vector field in which each point in space is associated with a unit vector in the direction of increasing Cartesian coordinate x, while Cartesian coordinates y and z are held constant. Defined everywhere.
yHat :: VectorField Source #
The vector field in which each point in space is associated with a unit vector in the direction of increasing Cartesian coordinate y, while Cartesian coordinates x and z are held constant. Defined everywhere.
zHat :: VectorField Source #
The vector field in which each point in space is associated with a unit vector in the direction of increasing Cartesian coordinate z, while Cartesian coordinates x and y are held constant. Defined everywhere.
Curves
Curve
is a parametrized function into three-space, an initial limit, and a final limit.
Curve | |
|
normalizeCurve :: Curve -> Curve Source #
Reparametrize a curve from 0 to 1.
Concatenate two curves.
concatenateCurves :: [Curve] -> Curve Source #
Concatenate a list of curves. Parametrizes curves equally.
reverseCurve :: Curve -> Curve Source #
Reverse a curve.
:: Curve | the curve |
-> Double | the parameter |
-> Position | position of the point on the curve at that parameter |
Evaluate the position of a curve at a parameter.
:: Displacement | amount to shift |
-> Curve | original curve |
-> Curve | shifted curve |
Shift a curve by a displacement.
The straight-line curve from one position to another.
Line Integrals
:: (InnerSpace v, Scalar v ~ Double) | |
=> Int | number of intervals |
-> Field v | scalar or vector field |
-> Curve | curve to integrate over |
-> v | scalar or vector result |
Calculates integral f dl over curve, where dl is a scalar line element.
:: Int | number of half-intervals (one less than the number of function evaluations) |
-> VectorField | vector field |
-> Curve | curve to integrate over |
-> Double | scalar result |
A dotted line integral.
Convenience function for compositeSimpsonDottedLineIntegral
.
:: Int | number of half-intervals (one less than the number of function evaluations) |
-> VectorField | vector field |
-> Curve | curve to integrate over |
-> Vec | vector result |
Calculates integral vf x dl over curve.
Convenience function for compositeSimpsonCrossedLineIntegral
.
Surfaces
Surface is a parametrized function from two parameters to space, lower and upper limits on the first parameter, and lower and upper limits for the second parameter (expressed as functions of the first parameter).
Surface | |
|
unitSphere :: Surface Source #
A unit sphere, centered at the origin.
centeredSphere :: Double -> Surface Source #
A sphere with given radius centered at the origin.
northernHemisphere :: Surface Source #
The upper half of a unit sphere, centered at the origin.
shiftSurface :: Displacement -> Surface -> Surface Source #
Shift a surface by a displacement.
Surface Integrals
:: (VectorSpace v, Scalar v ~ Double) | |
=> Int | number of intervals for first parameter, s |
-> Int | number of intervals for second parameter, t |
-> Field v | the scalar or vector field to integrate |
-> Surface | the surface over which to integrate |
-> v | the resulting scalar or vector |
A plane surface integral, in which area element is a scalar.
dottedSurfaceIntegral Source #
:: Int | number of intervals for first parameter, s |
-> Int | number of intervals for second parameter, t |
-> VectorField | the vector field to integrate |
-> Surface | the surface over which to integrate |
-> Double | the resulting scalar |
A dotted surface integral, in which area element is a vector.
Volumes
Volume is a parametrized function from three parameters to space, lower and upper limits on the first parameter, lower and upper limits for the second parameter (expressed as functions of the first parameter), and lower and upper limits for the third parameter (expressed as functions of the first and second parameters).
unitBallCartesian :: Volume Source #
A unit ball, centered at the origin. Specified in Cartesian coordinates.
centeredBall :: Double -> Volume Source #
A ball with given radius, centered at the origin.
Ball with given radius and center.
northernHalfBall :: Volume Source #
Upper half ball, unit radius, centered at origin.
centeredCylinder :: Double -> Double -> Volume Source #
Cylinder with given radius and height. Circular base of the cylinder is centered at the origin. Circular top of the cylinder lies in plane z = h.
shiftVolume :: Displacement -> Volume -> Volume Source #
Shift a volume by a displacement.
Volume Integral
:: (VectorSpace v, Scalar v ~ Double) | |
=> Int | number of intervals for first parameter (s) |
-> Int | number of intervals for second parameter (t) |
-> Int | number of intervals for third parameter (u) |
-> Field v | scalar or vector field |
-> Volume | the volume |
-> v | scalar or vector result |
A volume integral
Differential Equations
class (VectorSpace (Diff p), Fractional (Scalar (Diff p))) => StateSpace p where Source #
An instance of StateSpace
is a data type that can serve as the state
of some system. Alternatively, a StateSpace
is a collection of dependent
variables for a differential equation.
A StateSpace
has an associated vector space for the (time) derivatives
of the state. The associated vector space is a linearized version of
the StateSpace
.
(.-.) :: p -> p -> Diff p infix 6 Source #
Subtract points
(.+^) :: p -> Diff p -> p infixl 6 Source #
Point plus vector
StateSpace Double Source # | |
StateSpace Vec Source # | |
StateSpace Position Source # | Position is not a vector, but displacement (difference in position) is a vector. |
StateSpace St Source # | |
StateSpace p => StateSpace [p] Source # | |
(StateSpace p, StateSpace q, (~) * (Time p) (Time q)) => StateSpace (p, q) Source # | |
(StateSpace p, StateSpace q, StateSpace r, (~) * (Time p) (Time q), (~) * (Time q) (Time r)) => StateSpace (p, q, r) Source # | |
(.-^) :: StateSpace p => p -> Diff p -> p infixl 6 Source #
Point minus vector
type Time p = Scalar (Diff p) Source #
The scalars of the associated vector space can be thought of as time intervals.
type DifferentialEquation state = state -> Diff state Source #
A differential equation expresses how the dependent variables (state) change with the independent variable (time). A differential equation is specified by giving the (time) derivative of the state as a function of the state. The (time) derivative of a state is an element of the associated vector space.
type InitialValueProblem state = (DifferentialEquation state, state) Source #
An initial value problem is a differential equation along with an initial state.
type EvolutionMethod state Source #
= DifferentialEquation state | differential equation |
-> Time state | time interval |
-> state | initial state |
-> state | evolved state |
An evolution method is a way of approximating the state after advancing a finite interval in the independent variable (time) from a given state.
type SolutionMethod state = InitialValueProblem state -> [state] Source #
A (numerical) solution method is a way of converting an initial value problem into a list of states (a solution). The list of states need not be equally spaced in time.
stepSolution :: EvolutionMethod state -> Time state -> SolutionMethod state Source #
Given an evolution method and a time step, return the solution method which applies the evolution method repeatedly with with given time step. The solution method returned will produce an infinite list of states.
eulerMethod :: StateSpace state => EvolutionMethod state Source #
The Euler method is the simplest evolution method. It increments the state by the derivative times the time step.
rungeKutta4 :: StateSpace p => (p -> Diff p) -> Time p -> p -> p Source #
Take a single 4th-order Runge-Kutta step
integrateSystem :: StateSpace p => (p -> Diff p) -> Time p -> p -> [p] Source #
Solve a first-order system of differential equations with 4th-order Runge-Kutta
Visualization
Plotting
label :: String -> (Double, Double) -> Attribute Source #
An Attribute
with a given label at a given position.
postscript :: Attribute Source #
An Attribute
that requests postscript output.
Gloss library
cartToPolar :: (Float, Float) -> (Float, Float) Source #
theta=0 is positive x axis, output angle in radians
A think arrow