nullclines: Nullclines

View source: R/nullclines.R

nullclinesR Documentation

Nullclines

Description

Plots nullclines for two-dimensional autonomous ODE systems. Can also be used to plot horizontal lines at equilibrium points for one-dimensional autonomous ODE systems.

Usage

nullclines(
  deriv,
  xlim,
  ylim,
  parameters = NULL,
  system = "two.dim",
  points = 101,
  col = c("blue", "cyan"),
  add = TRUE,
  add.legend = TRUE,
  state.names = if (system == "two.dim") c("x", "y") else "y",
  ...
)

Arguments

deriv

A function computing the derivative at a point for the ODE system to be analysed. Discussion of the required structure of these functions can be found in the package vignette, or in the help file for the function ode.

xlim

In the case of a two-dimensional system, this sets the limits of the first dependent variable in which gradient reflecting line segments should be plotted. In the case of a one-dimensional system, this sets the limits of the independent variable in which these line segments should be plotted. Should be a numeric vector of length two.

ylim

In the case of a two-dimensional system this sets the limits of the second dependent variable in which gradient reflecting line segments should be plotted. In the case of a one-dimensional system, this sets the limits of the dependent variable in which these line segments should be plotted. Should be a numeric vector of length two.

parameters

Parameters of the ODE system, to be passed to deriv. Supplied as a numeric vector; the order of the parameters can be found from the deriv file. Defaults to NULL.

system

Set to either "one.dim" or "two.dim" to indicate the type of system being analysed. Defaults to "two.dim".

points

Sets the density at which derivatives are computed; points x points derivatives will be computed. Levels of zero gradient are identified using these computations and the function contour. Increasing the value of points improves identification of nullclines, but increases computation time. Defaults to 101.

col

In the case of a two-dimensional system, sets the colours used for the x- and y-nullclines. In the case of a one-dimensional system, sets the colour of the lines plotted horizontally along the equilibria. Should be a character vector of length two. Will be reset accordingly if it is of the wrong length. Defaults to c("blue", "cyan").

add

Logical. If TRUE, the nullclines are added to an existing plot. If FALSE, a new plot is created. Defaults to TRUE.

add.legend

Logical. If TRUE, a legend is added to the plots. Defaults to TRUE.

state.names

The state names for ode functions that do not use positional states.

...

Additional arguments to be passed to either plot or contour.

Value

Returns a list with the following components (the exact make up is dependent on the value of system):

add

As per input.

add.legend

As per input.

col

As per input, but with possible editing if a character vector of the wrong length was supplied.

deriv

As per input.

dx

A numeric matrix. In the case of a two-dimensional system, the values of the derivative of the first dependent derivative at all evaluated points.

dy

A numeric matrix. In the case of a two-dimensional system, the values of the derivative of the second dependent variable at all evaluated points. In the case of a one-dimensional system, the values of the derivative of the dependent variable at all evaluated points.

parameters

As per input.

points

As per input.

system

As per input.

x

A numeric vector. In the case of a two-dimensional system, the values of the first dependent variable at which the derivatives were computed. In the case of a one-dimensional system, the values of the independent variable at which the derivatives were computed.

xlim

As per input.

y

A numeric vector. In the case of a two-dimensional system, the of values of the second dependent variable at which the derivatives were computed. In the case of a one-dimensional system, the values of the dependent variable at which the derivatives were computed.

ylim

As per input.

Note

In order to ensure a nullcline is plotted, set xlim and ylim strictly enclosing its location. For example, to ensure a nullcline is plotted along x = 0, set ylim to, e.g., begin at -1.

Author(s)

Michael J Grayling

See Also

contour, plot

Examples

# Plot the flow field, nullclines and several trajectories for the
# one-dimensional autonomous ODE system logistic.
logistic_flowField  <- flowField(logistic,
                                 xlim       = c(0, 5),
                                 ylim       = c(-1, 3),
                                 parameters = c(1, 2),
                                 points     = 21,
                                 system     = "one.dim",
                                 add        = FALSE)
logistic_nullclines <- nullclines(logistic,
                                  xlim       = c(0, 5),
                                  ylim       = c(-1, 3),
                                  parameters = c(1, 2),
                                  system     = "one.dim")
logistic_trajectory <- trajectory(logistic,
                                  y0         = c(-0.5, 0.5, 1.5, 2.5),
                                  tlim       = c(0, 5),
                                  parameters = c(1, 2),
                                  system     = "one.dim")

# Plot the velocity field, nullclines and several trajectories for the
# two-dimensional autonomous ODE system simplePendulum.
simplePendulum_flowField  <- flowField(simplePendulum,
                                       xlim       = c(-7, 7),
                                       ylim       = c(-7, 7),
                                       parameters = 5,
                                       points     = 19,
                                       add        = FALSE)
y0                        <- matrix(c(0, 1, 0, 4, -6, 1, 5, 0.5, 0, -3),
                                    5, 2, byrow = TRUE)

simplePendulum_nullclines <- nullclines(simplePendulum,
                                        xlim       = c(-7, 7),
                                        ylim       = c(-7, 7),
                                        parameters = 5,
                                        points     = 500)

simplePendulum_trajectory <- trajectory(simplePendulum,
                                        y0         = y0,
                                        tlim       = c(0, 10),
                                        parameters = 5)

phaseR documentation built on Sept. 2, 2022, 5:07 p.m.