Package 'argosTrack'

Title: Fit Movement Models to Argos Data for Marine Animals
Description: Fits various state-space animal movement models to Argos data. Models are fitted by maximum likelihood with the Laplace approximation via 'TMB' and 'nlminb'. Functionality to visualize and simulate the models is available.
Authors: Christoffer Moesgaard Albertsen [aut, cre]
Maintainer: Christoffer Moesgaard Albertsen <[email protected]>
License: GPL-2
Version: 1.2.2-1
Built: 2024-11-23 02:46:46 UTC
Source: https://github.com/calbertsen/argosTrack

Help Index


Track data from the adult ringed seal with tag id 106373.

Description

A dataset containing tracking data from the Argos system for an adult ringed seal (Pusa hispida) released 10-30-2011 near Sanikiluaq, Nunavut, Canada.

Usage

adult_ringed_seal

Format

A data frame with 2698 rows and 5 variables:

id

Tag id

date

Date and time of the observation

lc

Argos location class of the observation

lat

Observed latitude (degrees north) from the Argos system

lon

Observed longitude (degrees east) from the Argos system

Source

Fisheries and Oceans Canada

References

Albertsen, C. M., Whoriskey, K., Yurkowski, D., Nielsen, A., and Flemming, J. M. (2015) Fast fitting of non-Gaussian state-space models to animal movement data via Template Model Builder. Ecology, 96(10), 2598-2604. doi: 10.1890/14-2101.1


Construct an Animal object

Description

Construct an Animal object

Usage

Animal(observation, movement, measurement, name = "", ...)

Arguments

observation

an Observation object

movement

a Movement object

measurement

a Measurement object

name

Name of the animal

...

Other arguments not currently used

Value

An Animal class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

Animal-class Movement-class Observation-class Measurement-class

Examples

d <- subadult_ringed_seal
obs <- Observation(lon = d$lon,
                   lat = d$lat,
                   dates = as.POSIXct(d$date,tz="GMT"),
                   locationclass = d$lc
                   )
meas <- Measurement(model="n")
mov <- RW(unique(obs$dates))
anim <- Animal(obs,mov,meas,"Subadult")
plot(anim)

A Reference Class for combining animal movement model parts.

Description

A Reference Class for combining animal movement model parts.

Fields

name

String identifying the animal

observation

Observation reference class

movement

Movement reference class

measurement

Measurement reference class

optim

List to keep results of optimization through fitTrack

Methods

addToLatPlot(obsarglist = list(pch = 16, col = "grey"), movarglist = list(), sdarglist = list(col = "grey", border = NA, lwd = 3, lty = 2), sd = TRUE)

Function to add estimated movement and observation latitudes to a plot. obsarglist is a list of base graphics arguments to be used for observation data, movarglist is a list of base graphics arguments to be used for movement estimates, and sdarglist is a list of base graphics arguments to be used for movement standard errors (if sd = TRUE).

addToLonPlot(obsarglist = list(pch = 16, col = "grey"), movarglist = list(), sdarglist = list(col = "grey", border = NA, lwd = 3, lty = 2), sd = TRUE)

Function to add estimated movement and observation longitudes to a plot. obsarglist is a list of base graphics arguments to be used for observation data, movarglist is a list of base graphics arguments to be used for movement estimates, and sdarglist is a list of base graphics arguments to be used for movement standard errors (if sd = TRUE).

addToMapPlot(obsarglist = list(type = "b", pch = 16, col = "grey"), movarglist = list())

Function to add estimated movement and observation coordinates to a plot. obsarglist is a list of base graphics arguments to be used for observation data, movarglist is a list of base graphics arguments to be used for movement estimates, and sdarglist is a list of base graphics arguments to be used for movement standard errors (if sd = TRUE).

getRange(sd = FALSE, zoomToState = FALSE)

Returns the range of animal coordinates. If zoomToState = TRUE, only movement data is used. If sd = TRUE, the range includes standard errors of the movement estimates.

getTMBdata()

Function to return a data list for TMB::MakeADFun

getTMBmap(fixcorrection = FALSE, ...)

Function to return a map list for TMB::MakeADFun. If fixcorrection = TRUE, the ratio between measurement errors for argos data is fixed.

getTMBparameters()

Function to return a parameter list for TMB::MakeADFun

initialize(observation, movement, measurement, name = "", ...)

Method to initialize the class. 'observation' is an Observation reference class object, 'movement' is a Movement reference class object, 'measurement' is a Measuremenbt reference class object, and 'name' is an identifier string.

simulate(newObject = FALSE)

Function to simulate movement and measurement data for the Animal. A list containing the movement (X) and measurement (Y) data is returned. If newObject = TRUE, an Animal object is added to the list.

updateFromFit(opt)

Function to save the result from nlminb.

Note

The reference class methods are not intended to be used by the user. Please use appropriate functions such as fitTrack, plot, plotLon, plotLat, plotMap, and simTrack.

Author(s)

Christoffer Moesgaard Albertsen

See Also

Animal Movement-class Observation-class Measurement-class


Defunct functions in the argosTrack package

Description

These functions are provided to guide users from old versions of the package to the new functionality.

Usage

argosTrack(...)

Arguments

...

Not used

Value

Nothing

argosTrack function

The argosTrack function is replaced by a combination of functions. To achieve the similar results as the old function call,

argosTrack <- function(lon,lat,dates,locationclass,
                      include = rep(TRUE,length(dates)),
                      equalbetas = TRUE,
                      timevarybeta = 1,
                      fixgammas = TRUE,
                      fixcorrection = FALSE,
                      dfVals = NULL,
                      dfMap = NULL,
                      minDf = 3.0,
                      errordistribution = "t",
                      movementmodel = "ctcrw",
                      verbose = TRUE,
                      timeunit = "mins",
                      nlminb.control = list(eval.max=2000,
                           iter.max=1500,
                            rel.tol=1e-3,
                            x.tol=1.5e-2))

use:

    obs <- Observation(lat = lat,
                       lon = lon,
                       dates = as.POSIXct(dates),
                       locationclass = locationclass,
                       include = include)
    mov <- CTCRW(dates = unique(as.POSIXct(dates)),
                 timeunit = timeunit)
    meas <- Measurement(model = model,minDf = minDf)
    anim <- Animal(obs,mov,meas)
    fitTrack(anim,
             silent = !verbose,
             fixcorrection = fixcorrection,
             nlminb.control = nlminb.control
             equaldecay = equalbetas,
             fixdrift = fixgammas)

where the input of the functions corresponds to the argument names of the argosTrack function.


Create a CSB movement model object

Description

Create a CSB movement model object

Usage

CSB(dates, pars = numeric(2), varPars = numeric(2),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Vector of movement parameters: \eqnlog(\beta), \eqn\gamma

varPars

Vector of movement variance parameters: \eqnlog(\sigma_S), \eqnlog(\sigma_B)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A CSB object

Author(s)

Christoffer Moesgaard Albertsen

See Also

CSB-class

Examples

d <- subadult_ringed_seal
mov <- argosTrack:::CSB(unique(as.POSIXct(d$date,tz="GMT")))

Create a CTCRW movement model object

Description

Create a CTCRW movement model object

Usage

CTCRW(dates, pars = numeric(4), varPars = numeric(2),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Vector of movement parameters: \eqnlog(\beta_lat), \eqnlog(\beta_lon), \eqn\gamma_lat, \eqn\gamma_lon

varPars

Vector of movement variance parameters: \eqnlog(\sigma_lat), \eqnlog(\sigma_lon)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A CTCRW object

Author(s)

Christoffer Moesgaard Albertsen

See Also

CTCRW-class

Examples

d <- subadult_ringed_seal
mov <- CTCRW(unique(as.POSIXct(d$date,tz="GMT")))

A Reference Class for fitting a Continuous Time Correlated Random Walk model.

Description

The reference class implements a Continuous Time Correlated Random Walk (Johnson et al. 2008). The velocity is modelled by an Ornstein Uhlenbeck process in each coordinate (cc),

νc,tγc=eβcΔt[νc,tΔtγc]+ηc,t,\nu_{c,t}-\gamma_c = e^{-\beta_c\Delta_{t}}[\nu_{c,t-\Delta_{t}}-\gamma_c]+\eta_{c,t},

and the locations are modelled by the integrated velocity,

Xc,t=Xc,tΔt+νc,tΔt(1eβcΔt)/βc+ζc,t.X_{c,t} = X_{c,t-\Delta_{t}}+\nu_{c,t-\Delta_{t}}\left(1-e^{-\beta_c\Delta_t}\right)/\beta_c+\zeta_{c,t}.

Here, Δt\Delta_t is the length of the time steps and ηc,t\eta_{c,t}, ζc,t\zeta_{c,t} are random variables. ηc,t\eta_{c,t} follows a Gaussian distribution with mean zero and variance V(ηc,t)=σc2(1e2βcΔt)/(2βc)V\left(\eta_{c,t}\right) = \sigma_c^2(1-e^{-2\beta_c\Delta_t})/(2\beta_c). The η\eta's are independent in time and coordiantes. ζc,t\zeta_{c,t} follows a Gaussian distribution with mean zero and variance

V(ζc,t)=σc2βc2(Δt2(1eβcΔt)/βc+(1e2βcΔt)/(2βc)).V\left(\zeta_{c,t}\right)=\frac{\sigma_c^2}{\beta_c^2}\left(\Delta_t-2\left(1-e^{-\beta_c\Delta_t}\right)/\beta_c+\left(1-e^{-2\beta_c\Delta_t}\right)/(2\beta_c)\right).

The ζ\zeta's are independent in time and coordiantes. The model correlates ηc,t\eta_{c,t} and ζc,t\zeta_{c',t'} by

Cov(ηc,t,ζc,t)=σc22βc2(12eβcΔt+e2βcΔi).Cov\left(\eta_{c,t},\zeta_{c,t}\right)= \frac{\sigma_c^2}{2\beta_c^2}\left(1-2e^{-\beta_c\Delta_t}+e^{-2\beta_c\Delta_i}\right).

For ccc\neq c' and ttt\neq t', Cov(ηc,t,ζc,t)=0Cov\left(\eta_{c,t},\zeta_{c',t'}\right)=0

Methods

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

simulate(x0 = c(0, 0, 0, 0))

Function to simulate from the movement model. The initial states (latitudianl velocity, longitudinal velocity, latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. If nauticalStates==TRUE, the result is returned in nautical miles.

References

Johnson, D. S., J. M. London, M. A. Lea, and J. W. Durban. (2008) Continuous-time correlated random walk model for animal telemetry data. Ecology 89, 1208-1215.

Albertsen, C. M., Whoriskey, K., Yurkowski, D., Nielsen, A., and Flemming, J. M. (2015) Fast fitting of non-Gaussian state-space models to animal movement data via Template Model Builder. Ecology, 96(10), 2598-2604. doi: 10.1890/14-2101.1

See Also

Movement-class CTCRW

Other "Movement models": CSB-class, DCRW-class, DIRAC-class, DSBHN-class, DSBW-class, GDCRW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class, RW-class


Create a DCRW movement model object

Description

Create a DCRW movement model object

Usage

DCRW(dates, pars = numeric(3), varPars = numeric(2),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct, equidistant, and increasing POSIXct dates

pars

Vector of movement parameters: \eqnlogit_(0,1)(\gamma), \eqn\phi, \eqnlogit_(-1,1)(\rho)

varPars

Vector of movement variance parameters: \eqnlog(\sigma_lat), \eqnlog(\sigma_lat)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A DCRW object

Author(s)

Christoffer Moesgaard Albertsen

See Also

DCRW-class

Examples

d <- subadult_ringed_seal
dates <- unique(as.POSIXct(d$date,tz="GMT"))
dseq <- seq(min(dates),max(dates), "day")
mov <- DCRW(dseq)

A Reference Class for fitting a discrete time first Difference Correlated Random Walk model.

Description

The reference class implements a discrete time first Difference Correlated Random Walk (Jonsen et al. 2005). The locations are modelled by

Xt=XtΔ+γT(ϕ)(XtXtΔ)+ϵtX_t = X_{t-\Delta} + \gamma T(\phi) (X_t - X_{t-\Delta}) + \epsilon_t

Here, ϵt\epsilon_t is zero mean Gaussian noise with covariance \pmatrix{ \sigma_1^2 & \rho\sigma_1\sigma_2 \cr \rho\sigma_1\sigma_2 & \sigma_2^2 }. T(ϕ)T(\phi) is the rotation matrix \pmatrix{ \cos(\phi) & -\sin(\phi) \cr \sin(\phi) & \cos(\phi) }. γ\gamma is a scalar.

Methods

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. If nauticalStates==TRUE, the result is returned in nautical miles.

References

Jonsen, I., J. Mills Flemming, and R. Myers. (2005) Robust state-space modeling of animal movement data. Ecology 86, 2874-2880.

See Also

Movement-class DCRW

Other "Movement models": CSB-class, CTCRW-class, DIRAC-class, DSBHN-class, DSBW-class, GDCRW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class, RW-class


Create a DIRAC movement model object

Description

Create a DIRAC movement model object

Usage

DIRAC(dates, pars = numeric(0), varPars = numeric(0),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Zero length vector of movement parameters

varPars

Zero length vector of movement variance parameters

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A DIRAC object

Author(s)

Christoffer Moesgaard Albertsen

See Also

DIRAC-class

Examples

d <- subadult_ringed_seal
mov <- DIRAC(unique(as.POSIXct(d$date,tz="GMT")))

A Reference Class for fitting a Dirac movement model.

Description

The reference class implements a Dirac movement model with no movement

Xtδ(X0))X_t \sim \delta(X_0))

Methods

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. The function only returns the locations. If nauticalStates==TRUE, the result is returned in nautical miles.

Author(s)

Christoffer Moesgaard Albertsen

See Also

Movement-class DIRAC

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DSBHN-class, DSBW-class, GDCRW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class, RW-class

Examples

d <- subadult_ringed_seal
mov <- DIRAC(unique(as.POSIXct(d$date,tz="GMT")))

Create a DSBHN movement model object

Description

Create a DSBHN movement model object

Usage

DSBHN(dates, pars = numeric(1), varPars = numeric(1),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Vector of movement parameters: \eqnlog(\rho)

varPars

Vector of movement variance parameters: \eqnlog(\sigma)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A DSBHN object

Author(s)

Christoffer Moesgaard Albertsen

See Also

DSBHN-class

Examples

d <- subadult_ringed_seal
dates <- unique(as.POSIXct(d$date,tz="GMT"))
dseq <- seq(min(dates),max(dates), "day")
mov <- DSBHN(dseq)

A Reference Class for fitting a Discrete time Step-length/Bearings model with Half Normal distributed steps.

Description

The reference class implements a Discrete time Step-length/Bearings model with Half Normal distributed steps.
Let XtX_t be the (bivariate) location (in nautical miles) at time tt. The step lengths are

St=(XtXt1)2S_t = \sqrt{\sum\left(X_t-X_{t-1}\right)^2}

and the bearings are

Bt=atan2(Xy,tXy,t1,Xx,tXy,t1)B_t = atan2(X_{y,t}-X_{y,t-1},X_{x,t}-X_{y,t-1})

where xx denotes the east/west coordinate and yy denotes the north/south coordinate. The step lengths are modelled by a half normal distribution

StHN(σ2)S_t \sim HN(\sigma^2)

and the bearings are modelled by a wrapped Cauchy distribution with location parameter Bt1B_{t-1} and concentration parameter ρ\rho.

Methods

getTMBparameters()

Function to return a parameter list for TMB::MakeADFun

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. If nauticalStates==TRUE, the result is returned in nautical miles.

Note

Step-length/Bearings models can be difficult for the Laplace approximation. Estimation problems can often be solved by decreasing the number of states relative to the number of observations.

See Also

Movement-class, link{DSBHN}, link{DSBW}.

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DIRAC-class, DSBW-class, GDCRW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class, RW-class


Create a DSBW movement model object

Description

Create a DSBW movement model object

Usage

DSBW(dates, pars = numeric(3), varPars = numeric(0),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Vector of movement parameters: \eqnlog(\rho), \eqnlog(shape), \eqnlog(scale)

varPars

Zeero length vector of movement variance parameters

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A DSBW object

Author(s)

Christoffer Moesgaard Albertsen

See Also

DSBW-class

Examples

d <- subadult_ringed_seal
dates <- unique(as.POSIXct(d$date,tz="GMT"))
dseq <- seq(min(dates),max(dates), "day")
mov <- DSBW(dseq)

A Reference Class for fitting a Discrete time Step-length/Bearings model with Weibull distributed steps.

Description

The reference class implements a Discrete time Step-length/Bearings model with Weibull distributed steps (McClintock et al. 2015, with small modifications).
Let XtX_t be the (bivariate) location (in nautical miles) at time tt. The step lengths are

St=(XtXt1)2S_t = \sqrt{\sum\left(X_t-X_{t-1}\right)^2}

and the bearings are

Bt=atan2(Xy,tXy,t1,Xx,tXy,t1)B_t = atan2(X_{y,t}-X_{y,t-1},X_{x,t}-X_{y,t-1})

where xx denotes the east/west coordinate and yy denotes the north/south coordinate. The step lengths are modelled by a half normal distribution

StWeibull(scale,shape)S_t \sim Weibull(scale,shape)

and the bearings are modelled by a wrapped Cauchy distribution with location parameter Bt1B_{t-1} and concentration parameter rhorho.

Methods

getTMBparameters()

Function to return a parameter list for TMB::MakeADFun

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. If nauticalStates==TRUE, the result is returned in nautical miles.

Note

Step-length/Bearings models can be difficult for the Laplace approximation. Estimation problems can often be solved by using DSBHN instead or by decreasing the number of states relative to the number of observations.

References

McClintock, B. T., London, J. M., Cameron, M. F. and Boveng, P. L. (2015) Modelling animal movement using the Argos satellite telemetry location error ellipse. Methods Ecol Evol, 6: 266-277. doi:10.1111/2041-210X.12311

See Also

Movement-class, DSBW, DSBHN.

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DIRAC-class, DSBHN-class, GDCRW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class, RW-class


Generic function to fit movement models.

Description

Generic function to fit movement models.

Usage

fitTrack(object, ...)

Arguments

object

The object to be fitted

...

Additional arguments passed to the appropriate method

Value

The optimization result

See Also

fitTrack,Animal-method


Fits a state-space movement model to animal data

Description

Fits the state-space movement determined by object using maximum likelihood estimation with the Laplace approximation.

Usage

## S4 method for signature 'Animal'
fitTrack(object, method = c("laplace", "mcmc"),
  silent = FALSE, fixcorrection = FALSE, nlminb.control = list(),
  inner.control = list(maxit = 1000), ...)

Arguments

object

Animal reference class object

method

Estimation method. Either laplace for maximum likelihood estimation with the Laplace approximation or mcmc for Markov Chain Monte Carlo (not yet implemented).

silent

Disable all tracing information?

fixcorrection

Should the ratio between measurement errors be fixed?

nlminb.control

List of control arguments passed to 'nlminb'

inner.control

List of control arguments passed to 'TMB::newton'

...

Additional arguments passed to object$getTMBmap

Value

The result returned from 'nlminb' is returned invisibly.

Note

The function changes 'object'.

Author(s)

Christoffer Moesgaard Albertsen

See Also

MakeADFun newton nlminb

fitTrack


Create a GDCRW movement model object

Description

Create a GDCRW movement model object

Usage

GDCRW(dates, pars = numeric(6), varPars = numeric(2),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Vector of movement parameters: \eqnlogit_(0,1)(\gamma_lat), \eqnlogit_(0,1)(gamma_lon), \eqn\phi, \eqnlogit_(-1,1)(\rho), \eqn\mu_lat, \eqn\mu_lon

varPars

Vector of movement variance parameters: \eqnlog(\sigma_lat), \eqnlog(\sigma_lon)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A GDCRW object

Author(s)

Christoffer Moesgaard Albertsen

See Also

GDCRW-class

Examples

d <- subadult_ringed_seal
mov <- GDCRW(unique(as.POSIXct(d$date,tz="GMT")))

A Reference Class for fitting a Generalized first Difference Correlated Random Walk model.

Description

The reference class implements an Irregularized first Difference Correlated Random Walk (Albertsen 2016). The velocity is modelled by the (bivariate) stochastic differential equation

dV_t = -\pmatrix{ -\log\gamma_1 & \phi \cr -\phi & -\log\gamma_2 } (V_t-\mu) d_t + S dB_t,

or for convenience

dVt=Θ(Vtμ)dt+SdBt.dV_t = - \Theta (V_t-\mu) d_t + S dB_t.

Since the locations are know from the velocity by Xt=i=0tViX_{t} = \sum_{i=0}^{t} V_i, the increments are known as

XtnXtn1=Vtn.X_{t_n} - X_{t_{n-1}} = V_{t_n}.

Hence,

X_{t_n} = X_{t_{n-1}} + \mu \pmatrix{ \cos ((t_n-t_{n-1})\phi) & -\sin ((t_n-t_{n-1})\phi) \cr \sin ((t_n-t_{n-1})\phi) & \cos ((t_n-t_{n-1})\phi) } \pmatrix{ \gamma_1^{t_n-t_{n-1}} & 0 \cr 0 & \gamma_2^{t_n-t_{n-1}} } (X_{t_{n-1}} - X_{t_{n-2}} - \mu)

Methods

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. If nauticalStates==TRUE, the result is returned in nautical miles.

References

Albertsen, C.M. (2016) Personal communication.

See Also

Movement-class, GDCRW. A generalization of DCRW, RW, and CTCRW (except only the location is modelled here and it is assumed to be known given the velocities)

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DIRAC-class, DSBHN-class, DSBW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class, RW-class


Generic function to extract movement track.

Description

Generic function to extract movement track.

Usage

getTrack(object, ...)

Arguments

object

The object from which to extract the track

...

Additional arguments passed to the appropriate method

Value

A data frame of positions

See Also

getTrack,Animal-method getTrack,Movement-method getTrack,Observation-method


Extracts track from animal data

Description

Extracts the estimated and observed track from an Animal object

Usage

## S4 method for signature 'Animal'
getTrack(object, ...)

Arguments

object

Animal reference class object

...

Additional arguments passed to the appropriate method

Value

A data frame of positions

Author(s)

Christoffer Moesgaard Albertsen


Extracts track from animal data

Description

Extracts the estimated track from a Movement object

Usage

## S4 method for signature 'Movement'
getTrack(object, ...)

Arguments

object

Movement-class reference class object

...

Additional arguments passed to the appropriate method

Value

A data frame of positions

Author(s)

Christoffer Moesgaard Albertsen


Extracts track from animal data

Description

Extracts the observed track from an Observation object

Usage

## S4 method for signature 'Observation'
getTrack(object, ...)

Arguments

object

Observation reference class object

...

Additional arguments passed to the appropriate method

Value

A data frame of positions

Author(s)

Christoffer Moesgaard Albertsen


Construct a Measurement object

Description

Construct a Measurement object

Usage

Measurement(model = c("n", "t", "sh"), logSdObs = c(0, 0),
  corObs = 0, logCorrection = matrix(c(0.6507, 0.8231, 2.3432, 2.0532,
  3.8514, 2.9602, 3.828, 3.0463, 4.4417, 3.6951, 5.5714, 5.5149, 8, 8), 2,
  7), splineXlogSd = 0, knotPars = numeric(5), df = numeric(8),
  minDf = 3, nauticalObs = FALSE, ...)

Arguments

model

Measurement model distribution

logSdObs

Variance parameters on log-scale

corObs

Correlation parameter

logCorrection

Correction parameters between location classes

splineXlogSd

Longitudinal variance parameter for spline model

knotPars

Latitudinal variance parameters for spline model

df

Degrees of freedom parameters for t distribution

minDf

Minimal degrees of freedom

nauticalObs

Transform observations from longitude/latitude to nautical miles?

...

Not currently used

Value

a Measurement object

Author(s)

Christoffer Moesgaard Albertsen

See Also

Measurement-class Movement-class Observation Animal

Examples

d <- subadult_ringed_seal
obs <- Observation(lon = d$lon,
                   lat = d$lat,
                   dates = as.POSIXct(d$date,tz="GMT"),
                   locationclass = d$lc
                   )
meas <- Measurement(model="t")
mov <- RW(unique(obs$dates))
anim <- Animal(obs,mov,meas,"Subadult")
plot(anim)

A Reference Class for specifying a measurement model.

Description

A Reference Class for specifying a measurement model.

Details

Following Albertsen et al. (2015), the reference class describes the measurement part of the state-space model for the movement data:

YtXtϵt.Y_t - X_t \sim \epsilon_t.

Here ,YtY_t is the observation, XtX_t is the state process at the corresponding time and epsilontepsilon_t is random noise (either Gaussian, t or symmetric hyperbolic) with scale matrix diag(τ2Ra(t),2,τ2Ra(t),2)diag(\tau_2 R_{a(t),2},\tau_2 R_{a(t),2}), where a(t)a(t) is the argos location class of the observation. If the time steps of the state process does not correspond to the observation time stamps (e.g. for a discrete time movement model), a linear interpolation between the last state before and the first state after the observation is used as XtX_t (as Jonsen et al. 2005).

The symmetric hyperbolic distribution is a special case of the distribution introduced by Barndorff-Nielsen (1977). The un-normalized density is

f(x;θ,μ,Σ)=exp(θ+θ(xμ)Σ1(xμ)T)f(x;\theta,\mu,\Sigma) = exp( - \sqrt{\theta + \theta (x-\mu)\Sigma^{-1} (x-\mu)^T})

where xx is a vector of observations, θ\theta is a scalar parameter, μ\mu is a vector of location parameters, and Σ\Sigma is a positive definite scale matrix.

Fields

model

Either "n" for multivariate Gaussian measurements, "t" for multivariate t-distributed measurements (Lange et al. 1989), or "sh" for multivariate symmetric hyperbolic distributed measurements of Argos data (Experimental; see Details).

logSdObs

Vector of scale parameters τ\tau (See Details)

corObs

Correlation parameter between coordiantes (not implemented yet)

logCorrection

Matrix of scale matrix ratios RR (See Details)

vcov

Covariance matrix of parameter estimates

splineXlogSd

Not used for Argos data

knotPars

Not used for Argos data

df

Degrees of freedom for model == "t"

minDf

Minimum degrees of freedom to allow in the model

nauticalObs

Should the observations be modelled in nautical miles?

reportSpline

Not used for Argos data

Methods

getTMBdata(dayOfYearSpline = c())

Function to return a data list for TMB::MakeADFun. Intended to be called from an Animal reference class object.

getTMBmap(correlated = FALSE, useSpline = FALSE, fixcorrection = FALSE, usedLevels = c("GPS", "3", "2", "1", "0", "A", "B", "Z"), numPerLevel = c(), needExtraVarPar = FALSE)

Function to return a map list for TMB::MakeADFun. Intended to be called from an Animal reference class object.

getTMBparameters()

Function to return a parameter list for TMB::MakeADFun. Intended to be called from an Animal reference class object.

initialize(model = c("n", "t", "sh"), logSdObs = c(0, 0), corObs = 0, logCorrection = matrix(c(0.6507, 0.8231, 2.3432, 2.0532, 3.8514, 2.9602, 3.828, 3.0463, 4.4417, 3.6951, 5.5714, 5.5149, 8, 8), 2, 7), splineXlogSd = 0, knotPars = numeric(5), df = numeric(8), minDf = 3, nauticalObs = FALSE, ...)

Method to initialize the class. Input corresponds to the reference class fields.

simulate(observation)

Function to simulate measurement data. 'observation' is an Observation reference class object, containing the time points and Argos location classes for the simulation. A matrix containing the measurement data is returned.

updateFromFit(parList, hessian, reportSpline)

Function to save the result from nlminb.

Note

The reference class methods are not intended to be used by the user. Please use appropriate functions such as simTrack.

Author(s)

Christoffer Moesgaard Albertsen

References

Albertsen, C. M., Whoriskey, K., Yurkowski, D., Nielsen, A., and Flemming, J. M. (2015) Fast fitting of non-Gaussian state-space models to animal movement data via Template Model Builder. Ecology, 96(10), 2598-2604. doi: 10.1890/14-2101.1

Barndorff-Nielsen, O. (1977) Exponentially decreasing distributions for the logarithm of particle size. Proc. R. Soc. Lond. A. 353, 401-419.

Jonsen, I., J. Mills Flemming, and R. Myers. (2005) Robust state-space modeling of animal movement data. Ecology 86, 2874-2880.

Lange, K., Roderick J. A. Little, & Jeremy M. G. Taylor. (1989). Robust Statistical Modeling Using the t Distribution. Journal of the American Statistical Association, 84(408), 881-896.

See Also

Measurement Movement-class, Observation-class, Animal-class


An abstract Reference Class for fitting a Movement model.

Description

Reference class implementing common functions for movement models.

Methods

addToLatPlot(arglist = list())

Function to add estimated movement latitudes to a plot. arglist is a list of base graphics arguments to be used for movement estimates.

addToLatPlotSd(arglist = list(col = "grey", border = NA, lwd = 3, lty = 2))

Function to add standard errors of estimated movement latitudes to a plot. arglist is a list of base graphics arguments to be used.

addToLonPlot(arglist = list())

Function to add estimated movement longitudes to a plot. arglist is a list of base graphics arguments to be used for movement estimates.

addToLonPlotSd(arglist = list(col = "grey", border = NA, lwd = 3, lty = 2))

Function to add standard errors of estimated movement longitudes to a plot. arglist is a list of base graphics arguments to be used.

addToMapPlot(arglist = list())

Function to add estimated movement coordinates to a plot. arglist is a list of base graphics arguments to be used.

getRange(sd = FALSE)

Returns the range of estimated coordinates. If sd = TRUE, the range includes standard errors of the movement estimates.

getTMBdata()

Function to return a data list for TMB::MakeADFun

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

getTMBparameters()

Function to return a parameter list for TMB::MakeADFun

updateFromFit(parList, hessian, sdParList = NULL)

Function to save the result from nlminb.

Note

This reference class should not be used as it does not implement a model. The class is only intended to be a superclass for actual implementations of movement models. A movement model implementation must include an initialization method, a simulation method and (preferably) a getTMBmap method.

See Also

Used in Animal. Currently the following movement models are implemented: CTCRW, DCRW, DSBHN, DSBW, RW.

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DIRAC-class, DSBHN-class, DSBW-class, GDCRW-class, MPCTCRW-class, OUL-class, OUV-class, RW-class


A Reference Class for argosTrack observations.

Description

A Reference Class for argosTrack observations.

Fields

lat

Latitude of observations

lon

Longitude of observations

dates

POSIXct vector of time stamps for observations

dayOfYear

Day of year corresponding to dates

qual

Factor of Argos location quality classes

varModelCode

Integer vector of variance model codes (0, 1, 2, or 3)

include

Logical vector of whether observations should be included in likelihood

Methods

addToLatPlot(arglist = list(pch = 16, col = "grey"))

Function to add observed latitudes to a plot. arglist is a list of base graphics arguments to be used for plotting.

addToLonPlot(arglist = list(pch = 16, col = "grey"))

Function to add observed longitudes to a plot. arglist is a list of base graphics arguments to be used for plotting.

addToMapPlot(arglist = list(type = "b", pch = 16, col = "grey"))

Function to add observed coordinates to a plot. arglist is a list of base graphics arguments to be used for plotting.

getRange()

Returns the range of observed coordinates.

getTMBdata()

Function to return a data list for TMB::MakeADFun

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

getTMBparameters()

Function to return a parameter list for TMB::MakeADFun

initialize(lon, lat, dates, locationclass, include = rep(TRUE, length(dates)), ...)

Method to initialize the class. 'lon' is a vector of longitudes, 'lat' is a vector of latitudes, 'dates' is a vector of POSIXct dates, 'locationclass' is a factor of Argos locatin classes, and Observation reference class object, 'movement' is a Movement reference class object, 'measurement' is a Measuremenbt reference class object, and 'include' is a logical vector of whether observations should be included in likelihood.

Author(s)

Christoffer Moesgaard Albertsen

See Also

Movement-class, Measurement, Animal

Examples

d <- subadult_ringed_seal
obs <- Observation(lon = d$lon,
                   lat = d$lat,
                   dates = as.POSIXct(d$date,tz="GMT"),
                   locationclass = d$lc
                   )
plot(obs)

Create an OUL movement model object

Description

Create an OUL movement model object

Usage

OUL(dates, pars = c(1, 0, 0, 1, 0, 0), varPars = numeric(2),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Vector of movement parameters: \eqnB_1,1, \eqnB_2,1, \eqnB_1,2, \eqnB_2,2, \eqn\gamma_lat, \eqn\gamma_lon

varPars

Vector of movement variance parameters: \eqnlog(\sigma_lat), \eqnlog(\sigma_lon)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A OUL object

Author(s)

Christoffer Moesgaard Albertsen

See Also

OUL-class

Examples

d <- subadult_ringed_seal
mov <- OUL(unique(as.POSIXct(d$date,tz="GMT")))

A Reference Class for fitting an Ornstein-Uhlenbeck Location model.

Description

The reference class implements an Ornstein-Uhlenbeck Location (Blackwell 2016). The locations follow an Ornstein-Uhlenbeck process

Xtγ=eBΔt[XtΔtγ]+ηtX_t - \gamma = e^{B\Delta_t}[X_{t-\Delta_t} - \gamma] + \eta_{t}

where ηt\eta_t is a bivariate normal with zero mean and covariance \pmatrix{\sigma_1^2 & 0 \cr 0 & \sigma_2^2}.

Methods

getTMBmap(...)

Function to return a map list for TMB::MakeADFun.

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. If nauticalStates==TRUE, the result is returned in nautical miles.

References

Blackwell, P. G., Niu, M., Lambert, M. S. and LaPoint, S. D. (2016), Exact Bayesian inference for animal movement in continuous time. Methods Ecol Evol, 7: 184-195. doi:10.1111/2041-210X.12460

See Also

Movement-class, OUL.

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DIRAC-class, DSBHN-class, DSBW-class, GDCRW-class, MPCTCRW-class, Movement-class, OUV-class, RW-class


Plotting an Animal object

Description

Plots time ×\times latitude, time ×\times longitude and longitude ×\times latitude for an Animal reference class object.

Usage

## S3 method for class 'Animal'
plot(x, ...)

Arguments

x

Animal reference class object

...

other parameters to be passed through to plotting functions.

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotMap, plotLat, and plotLon


Plotting a Movement object

Description

Plots time ×\times latitude, time ×\times longitude and longitude ×\times latitude for a Movement reference class object.

Usage

## S3 method for class 'Movement'
plot(x, ...)

Arguments

x

Movement reference class object

...

other parameters to be passed through to plotting functions.

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotMap, plotLat, and plotLon


Generic latitude plot

Description

Generic latitude plot

Usage

plotLat(object, ...)

Arguments

object

Object to plot

...

other parameters to be passed through to plotting functions.

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLat,Animal-method, plotLat,Movement-method, plotLat,Observation-method


Latitude plot for an Animal reference class object

Description

Latitude plot for an Animal reference class object

Usage

## S4 method for signature 'Animal'
plotLat(object, plotArgs = list(), args = list(lwd =
  3, col = "red"), add = FALSE, obsArgs = list(pch = 16),
  sdArgs = list(col = "grey", border = NA), sd = FALSE, ...)

Arguments

object

Animal reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting latitude movement data.

add

If FALSE a new plot window is created.

obsArgs

Arguments for plotting latitude observation data.

sdArgs

Arguments for plotting standard errors.

sd

Should standard errors be plotted?

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLat, plotLat,Movement-method, plotLat,Observation-method


Latitude plot for a Movement reference class object

Description

Latitude plot for a Movement reference class object

Usage

## S4 method for signature 'Movement'
plotLat(object, plotArgs = list(), args = list(),
  add = FALSE, sdArgs = list(col = "grey", border = NA), sd = FALSE,
  ...)

Arguments

object

Movement reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting latitude movement data.

add

If FALSE a new plot window is created.

sdArgs

Arguments for plotting standard errors.

sd

Should standard errors be plotted?

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLat, plotLat,Animal-method, plotLat,Observation-method


Latitude plot for an Observation reference class object

Description

Latitude plot for an Observation reference class object

Usage

## S4 method for signature 'Observation'
plotLat(object, plotArgs = list(),
  args = list(), add = FALSE, ...)

Arguments

object

Observation reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting latitude observation data.

add

If FALSE a new plot window is created.

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLat, plotLat,Animal-method, plotLat,Movement-method


Generic longitude plot

Description

Generic longitude plot

Usage

plotLon(object, ...)

Arguments

object

Object to plot

...

other parameters to be passed through to plotting functions.

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLon,Animal-method, plotLon,Movement-method, plotLon,Observation-method


Longitude plot for an Animal reference class object

Description

Longitude plot for an Animal reference class object

Usage

## S4 method for signature 'Animal'
plotLon(object, plotArgs = list(), args = list(lwd =
  3, col = "red"), add = FALSE, obsArgs = list(pch = 16),
  sdArgs = list(col = "grey", border = NA), sd = FALSE, ...)

Arguments

object

Animal reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting longitude movement data.

add

If FALSE a new plot window is created.

obsArgs

Arguments for plotting longitude observation data.

sdArgs

Arguments for plotting standard errors.

sd

Should standard errors be plotted?

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLon, plotLon,Movement-method, plotLon,Observation-method


Longitude plot for a Movement reference class object

Description

Longitude plot for a Movement reference class object

Usage

## S4 method for signature 'Movement'
plotLon(object, plotArgs = list(), args = list(),
  add = FALSE, sdArgs = list(col = "grey", border = NA), sd = TRUE,
  ...)

Arguments

object

Movement reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting longitude movement data.

add

If FALSE a new plot window is created.

sdArgs

Arguments for plotting standard errors.

sd

Should standard errors be plotted?

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotLon, plotLon,Animal-method, plotLon,Observation-method


Longitude plot for an Observation reference class object

Description

Longitude plot for an Observation reference class object

Usage

## S4 method for signature 'Observation'
plotLon(object, plotArgs = list(),
  args = list(), add = FALSE, ...)

Arguments

object

Observation reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting Longitude observation data.

add

If FALSE a new plot window is created.

...

additional arguments

Value

Invisibly returns the reference class object @seealso plotLon, plotLon,Animal-method, plotLon,Movement-method

Author(s)

Christoffer Moesgaard Albertsen


Generic map plot

Description

Generic map plot

Usage

plotMap(object, ...)

Arguments

object

Object to plot

...

other parameters to be passed through to plotting functions.

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotMap,Animal-method, plotMap,Movement-method, plotMap,Observation-method


Map plot for an Animal reference class object

Description

Map plot for an Animal reference class object

Usage

## S4 method for signature 'Animal'
plotMap(object, plotArgs = list(), args = list(lwd =
  3, col = "red"), add = FALSE, obsArgs = list(pch = 16),
  sdArgs = list(col = "grey", border = NA), ...)

Arguments

object

Animal reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting movement data.

add

If FALSE a new plot window is created.

obsArgs

Arguments for plotting observation data.

sdArgs

Arguments for plotting standard errors.

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotMap, plotMap,Movement-method, plotMap,Observation-method


Map plot for a Movement reference class object

Description

Map plot for a Movement reference class object

Usage

## S4 method for signature 'Movement'
plotMap(object, plotArgs = list(), args = list(pch
  = 16, lwd = 3, col = "red"), add = FALSE, sdArgs = list(col = "grey",
  border = NA), ...)

Arguments

object

Movement reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting movement data.

add

If FALSE a new plot window is created.

sdArgs

Arguments for plotting standard errors.

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotMap, plotMap,Animal-method, plotMap,Observation-method


Map plot for an Observation reference class object

Description

Map plot for an Observation reference class object

Usage

## S4 method for signature 'Observation'
plotMap(object, plotArgs = list(),
  args = list(pch = 16, lwd = 3, col = "red"), add = FALSE,
  sdArgs = list(col = "grey", border = NA), ...)

Arguments

object

Observation reference class object

plotArgs

Arguments to setup background plot

args

Arguments for plotting latitude observation data.

add

If FALSE a new plot window is created.

sdArgs

Arguments for plotting standard errors.

...

additional arguments

Value

Invisibly returns the reference class object

Author(s)

Christoffer Moesgaard Albertsen

See Also

plotMap, plotMap,Animal-method, plotMap,Movement-method


Generic roseplot

Description

Generic roseplot

Usage

roseplot(object, ...)

Arguments

object

Object to plot

...

other parameters

Author(s)

Christoffer Moesgaard Albertsen


Roseplot for Animal reference class object

Description

Calls roseplot on object$movement

Usage

## S4 method for signature 'Animal'
roseplot(object, ...)

Arguments

object

Animal reference class object

...

parameters passed to the method for the Movement reference class

Value

As roseplot,Movement-method

Author(s)

Christoffer Moesgaard Albertsen

See Also

roseplot,Movement-method


Roseplot for Movement reference class object

Description

Roseplot for Movement reference class object

Usage

## S4 method for signature 'Movement'
roseplot(object, nclass = 35, prob = TRUE,
  type = c("both", "step", "angle"), bearings = FALSE, show = TRUE,
  ...)

Arguments

object

Movement reference class object

nclass

Number of histogram classes

prob

Scale the frequencies to probabilities?

type

'step' for plotting step lengths, 'angle' for plotting angles and 'both' for plotting both.

bearings

Plot bearings (TRUE) or turning angles (FALSE)

show

Should a plot be produced?

...

not used

Value

Invisibly returns a list of bearings, step lengths and step lengths in km/h.

Author(s)

Christoffer Moesgaard Albertsen

Examples

d <- subadult_ringed_seal
obs <- Observation(lon = d$lon,
                   lat = d$lat,
                   dates = as.POSIXct(d$date,tz="GMT"),
                   locationclass = d$lc
                   )
meas <- Measurement(model="n")
mov <- RW(unique(obs$dates))
roseplot(mov)

Create a RW movement model object

Description

Create a RW movement model object

Usage

RW(dates, pars = numeric(0), varPars = numeric(2),
  nauticalStates = FALSE, timeunit = "hours")

Arguments

dates

Vector of distinct and increasing POSIXct dates

pars

Zero length vector of movement parameters

varPars

Vector of movement variance parameters: \eqnlog(\sigma_lat), \eqnlog(\sigma_lon)

nauticalStates

Should latent states be transformed from longitude/latitude to nautical miles?

timeunit

timeunit used for calculating time steps.

Value

A RW object

Author(s)

Christoffer Moesgaard Albertsen

See Also

RW-class

Examples

d <- subadult_ringed_seal
mov <- RW(unique(as.POSIXct(d$date,tz="GMT")))

A Reference Class for fitting a Random Walk model.

Description

The reference class implements a Random Walk (e.g. Nielsen et al. 2006).

XtN(XtΔt,Δtdiag(σ12,σ22))X_t \sim N(X_{t-\Delta_t},\Delta_t diag(\sigma_1^2,\sigma_2^2))

Methods

simulate(x0 = c(0, 0))

Function to simulate from the movement model. The initial states (latitudinal/y-coordinate location and longitudinal/x-coordinate location) must be given. The function only returns the locations. If nauticalStates==TRUE, the result is returned in nautical miles.

References

Nielsen, A., Bigelow, K. A., Musyl, M. K. and Sibert, J. R. (2006) Improving light-based geolocation by including sea surface temperature. Fisheries Oceanography, 15: 314-325. doi: 10.1111/j.1365-2419.2005.00401.x

See Also

Movement-class RW

Other "Movement models": CSB-class, CTCRW-class, DCRW-class, DIRAC-class, DSBHN-class, DSBW-class, GDCRW-class, MPCTCRW-class, Movement-class, OUL-class, OUV-class

Examples

d <- subadult_ringed_seal
mov <- RW(unique(as.POSIXct(d$date,tz="GMT")))

Generic function to simulate tracks

Description

Generic function to simulate tracks

Usage

simTrack(object, n, ...)

Arguments

object

Object to simulate from

n

Number of replications

...

other parameters

Author(s)

Christoffer Moesgaard Albertsen

See Also

simTrack,Animal-method, simTrack,Measurement-method, simTrack,Movement-method


Simulate from an Animal state-space model

Description

Simulate from an Animal state-space model

Usage

## S4 method for signature 'Animal'
simTrack(object, n = 1, newObject = TRUE, ...)

Arguments

object

Animal reference class describing the state-space model to simulate from

n

Number of replications

newObject

Should a new Animal object be added to the returned matrix?

...

Not used

Value

A (2 + newObject) x n matrix where the first row (X) contains lists with a matrix of simulated movement tracks, the second row (Y) contains lists with a matrix of simulated observations, and the third row (Animal - if present) contains lists with a new Animal object based on the simulated values.

Author(s)

Christoffer Moesgaard Albertsen

See Also

simTrack, simTrack,Measurement-method, simTrack,Movement-method


Simulate from a measurement model with observation information

Description

Simulate from a measurement model with observation information

Usage

## S4 method for signature 'Measurement'
simTrack(object, n = 1, observation, ...)

Arguments

object

Measurement reference class object to simulate from

n

Number of replications

observation

Observation reference class object with time points and Argos location class information to use in the simulation

...

Not used

Value

A 2 x number of time steps in observation x n array of simulated values

Author(s)

Christoffer Moesgaard Albertsen

See Also

simTrack, simTrack,Animal-method, simTrack,Movement-method


Simulate from a movement model

Description

Simulate from a movement model

Usage

## S4 method for signature 'Movement'
simTrack(object, n = 1, x0 = object$mu[, 1], ...)

Arguments

object

Movement reference class object implementing the model to simulate from

n

Number of replications

x0

Initial values

...

not used

Value

A 2 x number of time steps in object x n array of simulated values

Author(s)

Christoffer Moesgaard Albertsen

See Also

simTrack, simTrack,Animal-method, simTrack,Measurement-method


Track data from the adult ringed seal with tag id 43847.

Description

Track data from the adult ringed seal with tag id 43847.

Usage

subadult_ringed_seal

Format

A data frame with 3586 rows and 5 variables:

id

Tag id

date

Date and time of the observation

lc

Argos location class of the observation

lat

Observed latitude (degrees north) from the Argos system

lon

Observed longitude (degrees east) from the Argos system

Source

Fisheries and Oceans Canada

References

Albertsen, C. M., Whoriskey, K., Yurkowski, D., Nielsen, A., and Flemming, J. M. (2015) Fast fitting of non-Gaussian state-space models to animal movement data via Template Model Builder. Ecology, 96(10), 2598-2604. doi: 10.1890/14-2101.1