expokit_dgexpv: EXPOKIT DGEXPV matrix exponentiation on a square matrix

Description Usage Arguments Details Author(s) Examples

View source: R/kexpmv.R

Description

This function converts a matrix to CRS format and exponentiates it via the EXPOKIT dgexpv function with the wrapper functions wrapalldgexpv_ and dgexpv_ around DGEXPV. This function can be used to calculate both the matrix exponential in isolation or the product of the matrix exponential with a vector. This can be achieved by modifying the vector variable as shown below.

Usage

1
2
3
expokit_dgexpv(mat = NULL, t = 15, vector = NULL,
  transpose_needed = TRUE, transform_to_crs = TRUE, crs_n = NULL,
  anorm = NULL, mxstep = 10000, tol = as.numeric(1e-10))

Arguments

mat

an input square matrix.

t

a time value to exponentiate by.

vector

If NULL (default), the full matrix exponential is returned. However, in order to fully exploit the efficient speed of EXPOKIT on sparse matrices, this vector argument should be equal to a vector, v. This vector is an n dimensional vector, which in the Markovian case, can represent the starting probabilities.

transpose_needed

If TRUE (default), matrix will be transposed before the exponential operator is performed.

transform_to_crs

If the input matrix is in square format then the matrix will need transformed into CRS format. This is required for EXPOKIT's sparse-matrix functions DMEXPV and DGEXPV. Default TRUE; if FALSE, then the mat argument must be a CRS-formatted matrix.

crs_n

If a CRS matrix is input, crs_n specifies the order (# of rows or # of columns) of the matrix. Default is NULL.

anorm

The expokit_dgexpv routine requires an initial guess at the norm of the matrix. Using the R function norm might get slow with large matrices. Default is NULL.

mxstep

The EXPOKIT code performs integration steps and mxstep is the maximum number of integration steps performed. Default is 10000 steps; May need to increase this value if function outputs a warning message.

tol

the tolerance defined for the calculations.

Details

NOTE: When looking at DGEXPV vs. DMEXPV. According to the EXPOKIT documentation, the DGEXPV routine should be faster than DMEXPV. This is a result of the additional check the DMEXPV routine runs to check whether the output values satisfy the conditions needed for a Markovian model, which is not done in DGEXPV.

From EXPOKIT:

* The method used is based on Krylov subspace projection
* techniques and the matrix under consideration interacts only
* via the external routine 'matvec' performing the matrix-vector
* product (matrix-free method).

It should be noted that both the DMEXPV and DGEXPV functions within EXPOKIT require the matrix-vector product y = A*x = Q'*x i.e, where A is the transpose of Q. Failure to remember this leads to wrong results.

CRS (Compressed Row Storage) format is a compressed format that is required for EXPOKIT's sparse-matrix functions such as DGEXPV and DMEXPV. However this format is not necessary in EXPOKIT's padm-related functions.

This function is recommended for large sparse matrices, however the infinity norm of the matrix proves to be crucial when selecting the most efficient routine. If the infinity norm of the large sparse matrix is approximately >100 may be of benefit to use the expokit_dgpadm function for faster computations.

Author(s)

Meabh G. McCurdy mmccurdy01@qub.ac.uk

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
# Make a square (n x n) matrix A
# Use expokit_dgexpv to calculate both exp(At) and exp(At)v, where t is a 
# time value and v is an n dimensional column vector.
mat=matrix(c(-0.071207, 0.065573, 0.005634, 0, -0.041206, 0.041206, 0, 0, 0), 
nrow=3, byrow=TRUE)

# Set the time variable t
 t=15

# Exponentiate with EXPOKIT's dgexpv to obtain the full matrix exponential
OutputMat = expokit_dgexpv(mat=mat, t=t, transpose_needed=TRUE, vector=NULL)

print(OutputMat$output_mat)
print(OutputMat$message)

# Can also calculate the matrix exponential with the product of a vector.
# Create the n dimensional vector
v = matrix(0,3,1)
v[1] = 1

# Exponentiate with EXPOKIT's dgexpv
OutputMat = expokit_dgexpv(mat=mat, t=t, transpose_needed=TRUE, vector=v)

print(OutputMat$output_probs)
print(OutputMat$message)

# If message is 'NULL' then no error has occured and the number of 
# mxsteps defined in the function is acceptable.

kexpmv documentation built on May 1, 2019, 7:56 p.m.