dce_bayes: Bayesian Methods for Pharmacokinetic Modeling of Dynamic...

Description Usage Arguments Details Value Author(s) References See Also Examples

Description

Bayesian analysis of contrast agent concentration time curves from DCE-MRI.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
dcemri.bayes(conc, time, img.mask, model="extended", aif=NULL,
             user=NULL, nriters=3500, thin=3, burnin=500,
             tune=267, tau.ktrans=1, tau.kep=tau.ktrans,
             ab.vp=c(1,19), ab.tauepsilon=c(1,1/1000),
             samples=FALSE, multicore=FALSE, verbose=FALSE, ...)
dcemri.bayes.single(conc, time, nriters=3500, thin=3,
                    burnin=500, tune=267, tau.gamma=1,
                    tau.theta=1, ab.vp=c(1,19),
                    ab.tauepsilon=c(1,1/1000), aif.model=0,
                    aif.parameter=c(2.4,0.62,3,0.016), vp=1)

Arguments

conc

Matrix or array of concentration time series (last dimension must be time).

time

Time in minutes.

img.mask

Mask matrix or array. Voxels with mask=0 will be excluded.

model

is a character string that identifies the type of compartmental model to be used. Acceptable models include:

  • “weinmann”Tofts & Kermode AIF convolved with single compartment model

  • “extended”Weinmann model extended with additional vascular compartment (default)

  • “orton.exp”Extended model using Orton's exponential AIF

aif

is a character string that identifies the parameters of the type of arterial input function (AIF) used with the above model. Acceptable values are: tofts.kermode (default) or fritz.hansen for the weinmann and extended models; orton.exp (default) or user for the orton.exp model.

user

Vector of AIF parameters. For Tofts and Kermode: a_1, m_1, a_2, m_2; for Orton et al.: A_b, μ_b, A_g, μ_g.

nriters

Total number of iterations.

thin

Thining factor.

burnin

Number of iterations for burn-in.

tune

Number for iterations for tuning. The algorithm will be tuned to an acceptance rate between 0.3 and 0.6.

tau.ktrans

Variance parameter for prior on \log(K^{trans}).

tau.kep

Variance parameter for prior on \log(k_{ep}).

ab.vp

Hyper-prior parameters for the Beta prior on vp.

ab.tauepsilon

Hyper-prior parameters for observation error Gamma prior.

samples

If TRUE output includes samples drawn from the posterior distribution for all parameters.

multicore

If TRUE algorithm is parallelized using multicore.

verbose
tau.gamma
tau.theta
aif.model
aif.parameter
vp
...

Details

See Schmid et al. (2006) for details.

Value

Parameter estimates and their standard errors are provided for the masked region of the multidimensional array. They include

ktrans

Transfer rate from plasma to the extracellular, extravascular space (EES).

ktranserror

Error on ktrans.

kep

Rate parameter for transport from the EES to plasma.

keperror

Error on kep.

ve

Fractional occupancy by EES (the ratio between ktrans and kep).

vperror

Error on ve.

vp

Fractional occupancy by plasma.

sigma2

The residual sum-of-squares from the model fit.

time

Acquisition times (for plotting purposes).

Note, not all parameters are available under all models choices.

Author(s)

Volker Schmid

References

Schmid, V., Whitcher, B., Padhani, A.R., Taylor, N.J. and Yang, G.-Z. (2006) Bayesian methods for pharmacokinetic models in dynamic contrast-enhanced magnetic resonance imaging, IEEE Transactions on Medical Imaging, 25 (12), 1627-1636.

See Also

dcemri.lm, dcemri.map, dcemri.spline

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
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
data("buckley")
xi <- seq(5, 300, by=5)
img <- array(t(breast$data)[,xi], c(13,1,1,60))
mask <- array(TRUE, dim(img)[1:3])
time <- buckley$time.min[xi]

## Bayesian estimation with Fritz-Hansen default AIF
fit.bayes <- dcemri.bayes(img, time, mask, aif="fritz.hansen", nriters=1500)

## Bayesian estimation with "orton.exp" function fit to Buckley's AIF
aif <- buckley$input[xi]
aifparams <- orton.exp.lm(time, aif)
aifparams$D <- 1
fit.bayes.aif <- dcemri.bayes(img, time, mask, model="orton.exp",
                              aif="user", user=aifparams, nriters=1500)

plot(breast$ktrans, fit.bayes$ktrans, xlim=c(0,1), ylim=c(0,1),
     xlab=expression(paste("True ", K^{trans})),
     ylab=expression(paste("Estimated ", K^{trans}, " (Bayesian)")))
points(breast$ktrans, fit.bayes.aif$ktrans, pch=2)
abline(0, 1, lwd=1.5, col=2)
legend("right", c("orton.exp", "fritz.hansen"), pch=1:2)

cbind(breast$ktrans, fit.bayes$ktrans[,,1], fit.bayes.aif$ktrans[,,1])

## Not run: 
fit.lm <- dcemri.lm(img, time, mask, aif="fritz.hansen")
fit.lm.aif <- dcemri.lm(img, time, mask, model="orton.exp", aif="user",
                        user=aifparams)

plot(breast$ktrans, fit.bayes$ktrans, xlim=c(0,1), ylim=c(0,1),
     xlab=expression(paste("True ", K^{trans})),
     ylab=expression(paste("Estimated ", K^{trans})))
points(breast$ktrans, fit.bayes.aif$ktrans, pch=2)
points(breast$ktrans, fit.lm$ktrans, pch=3)
points(breast$ktrans, fit.lm.aif$ktrans, pch=4)
abline(0, 1, lwd=1.5, col="red")
legend("bottomright", c("Bayesian Estimation (fritz-hansen)",
                        "Bayesian Estimation (orton.exp)",
                        "Levenburg-Marquardt (fritz.hansen)",
                        "Levenburg-Marquardt (orton.exp)"), pch=1:4)

## End(Not run)

dcemri documentation built on May 2, 2019, 5:27 p.m.