bootstrap.hankel | R Documentation |
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.
bootstrap.hankel(cf, t0 = 1, n = 2, N = (cf$Time/2 + 1), t0fixed = TRUE, deltat = 1, Delta = 1, submatrix.size = 1, element.order = 1)
cf |
object of type cf |
t0 |
Integer. Initial time value of the GEVP, must be in between 0 and
|
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 |
deltat |
Integer. value of deltat used in the hankel GEVP. Default is 1. Used
|
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 |
See vignette(name="hankel", package="hadron")
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.
Other hankel:
bootstrap.hankel_summed()
,
gevp.hankel_summed()
,
gevp.hankel()
,
hankel2cf()
,
hankel2effectivemass()
,
plot_hankel_spectrum()
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)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.