kin-class {TIMP} | R Documentation |
Class "kin" for kinetic model storage.
Description
kin
is the class for kinetic models;
an object
of class "kin" is initialized if
mod_type = "kin"
is an
argument of initModel
.
All objects of class kin
are sub-classes of
class dat
; see documentation for dat
for a description of
these slots.
Details
See dat-class
for an
example of the initialization of a
kin
object via the initModel
function.
Objects from the Class
Objects can be created by calls of the form new("kin", ...)
or
kin(...)
. Slots whose
description are marked with *** may
be specified in the ...
argument of the initModel
function.
Slots
- anipar
- anispec
- autoclp0
- C2
- chinde
- clinde
- clp0
- clpCon
- clpdep
- clpequ
- clpequspecBD
- clpType
- cohcol
- cohirf
- datafile
- datCall
- drel
- dscal
- dscalspec
dummy
:Object of class
"list"
of dummy parameters which can be used in complex relations- E2
- fixed
- fixedkmat
- free
- fvecind
- getX
- getXsuper
- highcon
- inten
- kin2scal
- kinpar2
- kinscalspecial
- kinscalspecialspec
- lclp0
- lclpequ
- title
- parnames
- prel
- prelspec
- psi.df
- psi.weight
- pvecind
- satMat
- scalx
- usecompnames0
- usecompnamesequ
- usekin2
- weight
- weightList
- weightM
- weightpar
- weightsmooth
- x
- x2
- clpequspec
- compnames
- constrained
- iter
- lightregimespec
- lowcon
- makeps
- mhist
- mod_type
- mvecind
- ncomp
- nl
- nt
- nvecind
- outMat
- positivepar
- sigma
- simdata
- speckin2
- kinpar
*** vector of rate constants to be used as starting values for the exponential decay of components; the length of this vector determines the number of components of the kinetic model.
specpar
:*** Object of class
"list"
parameters for spectral constraintsseqmod
:*** Object of class
"logical"
that isTRUE
if a sequential model is to be applied andFALSE
otherwiseirf
:Object of class
"logical"
that isTRUE
is an IRF is modeled andFALSE
otherwisemirf
:Object of class
"logical"
that isTRUE
if a measured IRF is modeled andFALSE
otherwisemeasured_irf
:*** Object of class
"vector"
containing a measured IRFconvalg
:*** Object of class
"numeric"
1-3 determining the numerical convolution algorithm used in the case of modeling a measured IRF; if3
then supply a reference lifetime in the slotreftau
.reftau
:*** Object of class
"numeric"
containing a reference lifetime to be used whenconvalg=3
irffun
:*** Object of class
"character"
describing the function to use to describe the IRF, by default "gaus"irfpar
:*** Object of class
"vector"
of IRF parameters; for the common Gaussian IRF this vector is orderedc(location, width)
dispmu
:Object of class
"logical"
that isTRUE
if dispersion of the parameter for IRF location is to be modeled andFALSE
otherwisedispmufun
:***Object of class
"character"
describing the functional form of the dispersion of the IRF location parameter; if equal to "discrete" then the IRF location is shifted per element ofx2
andparmu
should have the same length asx2
. defaults to a polynomial descriptionparmu
:*** Object of class
"list"
of starting values for the dispersion model for the IRF locationdisptau
:Object of class
"logical"
that isTRUE
if dispersion of the parameter for IRF width is to be modeled andFALSE
otherwisedisptaufun
:*** Object of class
"character"
describing the functional form of the dispersion of the IRF width parameter; if equal to"discrete"
then the IRF width is parameterized per element ofx2
andpartau
should have the same length asx2
. defaults to a polynomial descriptionpartau
:*** Object of class
"vector"
of starting values for the dispersion model for the IRF FWHMfullk
:Object of class
"logical"
that isTRUE
if the data are to be modeled using a compartmental model defined in a K matrix andFALSE
otherwisekmat
:*** Object of class
"array"
containing the K matrix descriptive of a compartmental modeljvec
:*** Object of class
"vector"
containing the J vector descriptive of the inputs to a compartmental modelncolc
:Object of class
"vector"
describing the number of columns of the C matrix for each clp inx2
kinscal
:*** Object of class
"vector"
of starting values for branching parameters in a compartmental modelkmatfit
:Object of class
"array"
of fitted values for a compartmental modelcohspec
:*** Object of class
"list"
describing the model for coherent artifact/scatter component(s) containing the elementtype
and optionally the elementnumdatasets
. The elementtype
can be set as follows:"irf"
:if
type="irf"
, the coherent artifact/scatter has the time profile of the IRF."freeirfdisp"
:if
type="freeirfdisp"
, the coherent artifact/scatter has a Gaussian time profile whose location and width are parameterized in the vectorcoh
."irfmulti"
:if
type="irfmulti"
the time profile of the IRF is used for the coherent artifact/scatter model, but the IRF parameters are taken per dataset (for the multidataset case), and the integer argumentnumdatasets
must be equal to the number of datasets modeled."seq"
:if
type="seq"
a sequential exponential decay model is applied, whose starting value are contained in an additional list elementstart
. This often models oscillating behavior well, where the number of oscillations is the number of parameter starting values given instart
. The starting values after optimization will be found in the slotcoh
of the object of classtheta
corresponding to each dataset modeled."mix"
:if
type="mix"
iftype="mix"
a sequential exponential decay model is applied along with a model that follows the time profile of the IRF; the coherent artifact/scatter is then a linear superposition of these two models; see the above description ofseq
for how to supply the starting values.
coh
:*** Object of class
"vector"
of starting values for the parameterization of a coherent artifactoscspec
:*** Object of class
"list"
describing the model for additional oscillation component(s) containing the elementtype
and optionally the elementstart
. The elementstart
can be used to specify the starting values for the oscillation function. The elementtype
can be set as follows:"harmonic"
:if
type="harmonic"
, the oscillation function is a damped harmonic oscillator.
oscpar
:*** Object of class
"vector"
of starting values for the oscillation parameterswavedep
:Object of class
"logical"
describing whether the kinetic model is dependent onx2
index (i.e., whether there is clp-dependence)lambdac
:*** Object of class
"numeric"
for the center wavelength to be used in a polynomial description ofx2
-dependenceamplitudes
:*** Object of class
"vector"
that may be used to multiply the concentrations by a square diagonal matrix with the number of columns that the concentration matrix has; the diagonal is given inamplitudes
and these values will be treated as parameters to be optimized.streak
:*** Object of class
"logical"
that defaults toFALSE
; ifstreak=TRUE
then the period of the laser is expected viastreakT
.streakT
:*** Object of class
"numeric"
the period of the laser; this will be used to add a backsweep term to the concentration matrix and should be set in conjunctionstreak=TRUE
.doublegaus
:*** Object of class
"logical"
that defaults toFALSE
and determines whether a double Gaussian should be used to model the IRF. Ifdoublegaus=TRUE
thenirfpar
should contain four numeric values corresponding to the location (mean) of the IRF, the FWHM of the first Gaussian, the FWHM of the second Gaussian, and the relative amplitude of the second Gaussian, respectively.multiplegaus
:*** Object of class
"logical"
that defaults toFALSE
and determines whether multiple Gaussians should be used to model the IRF. Ifmultiplegaus=TRUE
thenirfpar
should contain: two numeric values corresponding to the location (mean) and the FWHM of the first Gaussian of the IRF, and three numeric values for each additional gaussian modeled, corresponding to the relative scaling to the first gaussian, the shift (in time) relative to the first gaussian and the FWHM of the additional Gaussian, respectively.numericalintegration
:*** Object of class
"logical"
that defaults toFALSE
and determines whether a kinetic theory model of a reaction mechanism should be numerically integrated (using deSolve) to find the concentrations. Ifnumericalintegration=TRUE
theninitialvals
should specify the initial concentrations andreactantstoichiometrymatrix
andstoichiometrymatrix
should specify the reaction mechanism, as per Puxty et. al. (2006).initialvals
:*** Object of class
"vector"
giving the concentrations at the initial time step.reactantstoichiometrymatrix
:*** Object of class
"vector"
giving the (integer) stoichiometric coefficients for the reactants; this is the matrix Xr of Puxty et. al. (2006) withdim=NULL
.stoichiometrymatrix
:*** Object of class
"vector"
giving the (integer) stoichiometric coefficients for the reactions; this is the matrix X of Puxty et. al. (2006) withdim=NULL
.
Extends
Class dat-class
, directly.
Author(s)
Katharine M. Mullen, David Nicolaides, Ivo H. M. van Stokkum
References
Puxty, G., Maeder, M., and Hungerbuhler, K. (2006) Tutorial on the fitting of kinetics models to multivariate spectroscopic measurements with non-linear least-squares regression, Chemometrics and Intelligent Laboratory Systems 81, 149-164.
See Also
Examples
## Example in modeling second order kinetics, by
## David Nicolaides.
## On simulated data.
##############################
## load TIMP
##############################
library("TIMP")
##############################
## SIMULATE DATA
##############################
## set up the Example problem, a la in-situ UV-Vis spectroscopy of a simple
## reaction.
## A + 2B -> C + D, 2C -> E
cstart <- c(A = 1.0, B = 0.8, C = 0.0, D = 0.0, E = 0.0)
times <- c(seq(0,2, length=21), seq(3,10, length=8))
k <- c(kA = 0.5, k2C = 1)
## stoichiometry matrices
rsmatrix <- c(1,2,0,0,0,0,0,2,0,0)
smatrix <- c(-1,-2,1,1,0,0,0,-2,0,1)
concentrations <- calcD(k, times, cstart, rsmatrix, smatrix)
wavelengths <- seq(500, 700, by=2)
spectra <- matrix(nrow = length(wavelengths), ncol = length(cstart))
location <- c(550, 575, 625, 650, 675)
delta <- c(10, 10, 10, 10, 10)
spectra[, 1] <- exp( - log(2) *
(2 * (wavelengths - location[1])/delta[1])^2)
spectra[, 2] <- exp( - log(2) *
(2 * (wavelengths - location[2])/delta[2])^2)
spectra[, 3] <- exp( - log(2) *
(2 * (wavelengths - location[3])/delta[3])^2)
spectra[, 4] <- exp( - log(2) *
(2 * (wavelengths - location[4])/delta[4])^2)
spectra[, 5] <- exp( - log(2) *
(2 * (wavelengths - location[5])/delta[5])^2)
sigma <- .001
Psi_q <- concentrations %*% t(spectra) + sigma *
rnorm(dim(concentrations)[1] * dim(spectra)[1])
## store the simulated data in an object of class "dat"
kinetic_data <- dat(psi.df=Psi_q , x = times, nt = length(times),
x2 = wavelengths, nl = length(wavelengths))
##############################
## DEFINE MODEL
##############################
## starting values
kstart <- c(kA = 1, k2C = 0.5)
## model definition for 2nd order kinetics
kinetic_model <- initModel(mod_type = "kin", seqmod = FALSE,
kinpar = kstart,
numericalintegration = TRUE,
initialvals = cstart,
reactantstoichiometrymatrix = rsmatrix,
stoichiometrymatrix = smatrix )
##############################
## FIT INITIAL MODEL
## adding constraints to non-negativity of the
## spectra via the opt option nnls=TRUE
##############################
kinetic_fit <- fitModel(data=list(kinetic_data),
modspec = list(kinetic_model),
opt = kinopt(nnls = TRUE, iter=80,
selectedtraces = seq(1,kinetic_data@nl,by=2)))
## look at estimated parameters
parEst(kinetic_fit)
## various results
## concentrations
conRes <- getX(kinetic_fit)
matplot(times, conRes, type="b", col=1,pch=21, bg=1:5, xlab="time (sec)",
ylab="concentrations", main="Concentrations (2nd order kinetics)")
## spectra
specRes <- getCLP(kinetic_fit)
matplot(wavelengths, specRes, type="b", col=1,pch=21, bg=1:5,
xlab="wavelength (nm)",
ylab="amplitude", main="Spectra")
## see help(getResults) for how to get more results information from
## kinetic_fit
##############################
## CLEANUP GENERATED FILES
##############################
# This removes the files that were generated in this example
# (do not run this code if you wish to inspect the output)
file_list_cleanup = c(Sys.glob("*paramEst.txt"), Sys.glob("*.ps"), Sys.glob("Rplots*.pdf"))
# Iterate over the files and delete them if they exist
for (f in file_list_cleanup) {
if (file.exists(f)) {
unlink(f)
}
}