bootstrap.hankel: GEVP method based on Hankel matrices.

View source: R/hankel.R

bootstrap.hankelR Documentation

GEVP method based on Hankel matrices.

Description

Alternative method to determine energy levels from correlation matrices. A so-called Hankel matrix is generated from an input cf object and a generalised eigenvalue problem is solved then. This is the function to call. It will perform a bootstrap analysis.

Usage

bootstrap.hankel(cf, t0 = 1, n = 2, N = (cf$Time/2 + 1),
  t0fixed = TRUE, deltat = 1, Delta = 1, submatrix.size = 1,
  element.order = 1)

Arguments

cf

object of type cf

t0

Integer. Initial time value of the GEVP, must be in between 0 and Time/2-n. Default is 1. Used when t0fixed=TRUE.

n

Integer. Size of the Hankel matrices to generate

N

Integer. Maximal time index in correlation function to be used in Hankel matrix

t0fixed

Integer. If set to TRUE, keep t0 fixed and vary deltat, otherwise keep deltat fixed and vary t0.

deltat

Integer. value of deltat used in the hankel GEVP. Default is 1. Used t0fixed=FALSE

Delta

integer. Delta is the time shift used in the Hankel matrix.

submatrix.size

Integer. Submatrix size to be used in build of Hankel matrices. Submatrix.size > 1 is experimental.

element.order

Integer vector. specifies how to fit the n linearly ordered single correlators into the correlator matrix for submatrix.size > 1. element.order=c(1,2,3,4) leads to a matrix matrix(cf[element.order], nrow=2). Matrix elements can occur multiple times, such as c(1,2,2,3) for the symmetric case, for example.

Details

See vignette(name="hankel", package="hadron")

Value

List object of class "hankel". The eigenvalues are stored in a numeric vector t0, the corresonding samples in t. The reference input time t0 is stored as reference_time in the returned list.

See Also

Other hankel: bootstrap.hankel_summed(), gevp.hankel_summed(), gevp.hankel(), hankel2cf(), hankel2effectivemass(), plot_hankel_spectrum()

Examples


data(correlatormatrix)
correlatormatrix <- bootstrap.cf(correlatormatrix, boot.R=99, boot.l=1, seed=132435)
t0 <- 4
correlatormatrix.gevp <- bootstrap.gevp(cf=correlatormatrix, t0=t0, element.order=c(1,2,3,4))
pc1 <- gevp2cf(gevp=correlatormatrix.gevp, id=1)
pc1.hankel <- bootstrap.hankel(cf=pc1, t0=1, n=2)
hpc1 <- hankel2cf(hankel=pc1.hankel, id=1)
plot(hpc1, log="y")
heffectivemass1 <- hankel2effectivemass(hankel=pc1.hankel, id=1)

hadron documentation built on Sept. 9, 2022, 5:06 p.m.