API

ReplicateBE.rbe

ReplicateBE.rbeFunction
rbe(df; dvar::Symbol,
    subject::Symbol,
    formulation::Symbol,
    period::Symbol,
    sequence::Symbol,
    g_tol::Float64 = 1e-8, x_tol::Float64 = 0.0, f_tol::Float64 = 0.0,
    iterations::Int = 100,
    store_trace = false, extended_trace = false, show_trace = false,
    memopt = true,
    init = [],
    postopt = false, maxopttry = 100)

Mixed model fitting function for replicate bioequivalence without data preparation (apply categorical! for each factor and sort! to dataframe).

Mixed model in matrix form:

\[y = X\beta + Zu + \epsilon\]

with covariance matrix for each subject:

\[V_{i} = Z_{i}GZ_i'+R_{i}\]

Arguments

  • df - DataFrame with data;

Keywords:

  • dvar - dependent variable;
  • subject - subject factor;
  • formulation - formulation factor;
  • period - period factor;
  • sequence - sequence factor;
  • g_tol = 1e-8
  • x_tol = 0.0
  • f_tol = 0.0
  • iterations = 100 - maximum iteration for optimization
  • store_trace = false
  • extended_trace = false
  • show_trace = false
  • memopt = true - memory optimization (function cache)
  • init = [] - initial variance paremeters
  • postopt = false - post optimization
source

ReplicateBE.rbe!

ReplicateBE.rbe!Function

This function can convert non-categorical to categorical and sort dataframe.

It can takes more time, but can help to avoid some errors like: "ERROR: type ContinuousTerm has no field contrasts".

source

StatsBase.coef

StatsBase.coefFunction
coef(rbe::RBE)

Return model coefficients.

\[\beta = {(\sum_{i=1}^n X_{i}'V_i^{-1}X_{i})}^{-1}(\sum_{i=1}^n X_{i}'V_i^{-1}y_{i})\]
source

StatsBase.confint

StatsBase.confintFunction
confint(obj::RBE; level::Real=0.95, expci::Bool = false, inv::Bool = false, df = :sat)

Compute confidence intervals for coefficients, with confidence level level (by default 95%).

Arguments

  • expci = true: return exponented CI.

  • inv = true: return -estimate ± t(alpha, df)*SE

  • df = :sat: use Satterthwaite DF approximation.

  • df = :df3 or df = :cont: DF (contain) = N - rank(ZX).

\[CI = estimate ± t(alpha, df)*SE\]
source

StatsBase.dof

StatsBase.dofFunction
dof(rbe::RBE)

Return the number of degrees of freedom for the coefficients of the model.

source

StatsBase.stderror

StatsBase.stderrorFunction
StatsBase.stderror(rbe::RBE)

Return the standard errors for the coefficients of the model.

\[se = \sqrt{LCL'}\]

where

\[C = (\sum_{i=1}^{n} X_i'V_i^{-1}X_i)^{-1}\]
source

ReplicateBE.coefnum

ReplicateBE.contrast

ReplicateBE.contrastFunction
contrast(rbe::RBE, L::Matrix; name = "Contrast", memopt = true)::ContrastTable

Return contrast table for L matrix. Table include:

  • F
  • NumDF
  • DF
  • P|f|
\[F = \frac{\beta'L'(LCL')^{-1}L\beta}{rank(LCL')}\]

where

\[C = (\sum_{i=1}^{n} X_i'V_i^{-1}X_i)^{-1}\]

DF for one-dimetion case:

\[df = \frac{2(LCL')^{2}}{g'Ag}\]

where $A = 2H^{-1}$

where $g = \triangledown _{\theta}(LC_{\theta}L')$

DF for multi-dimention case see Schaalje et al 2002.

p value calculated with:

pval    = ccdf(FDist(numdf, df), F)
source

ReplicateBE.design

ReplicateBE.estimate

ReplicateBE.estimateFunction
estimate(rbe::RBE, L::Matrix; df = :sat, name = "Estimate", memopt = true, alpha = 0.05)

Return estimate table for L 1xp matrix. Table include:

  • estimate value
  • SE
  • DF
  • t
  • P|t|
  • CI Upper
  • CI Lower
\[estimate = L\beta\]
\[SE = \sqrt{LCL'}\]
\[t = estimate / SE\]

For df = :sat:

\[df = \frac{2(LCL')^{2}}{g'Ag}\]

where

\[C = (\sum_{i=1}^{n} X_i'V_i^{-1}X_i)^{-1}\]

where $A = 2H^{-1}$

where $g = \triangledown _{\theta}(LC_{\theta}L')$

For df = :cont (contain):

\[df = N - rank(ZX)\]

CI estimate is:

\[CI = estimate ± t(alpha, df)*SE\]

Example of L matrix if length of fixed effect vector is 6, estimate for 4-th value:

L = [0 0 0 1 0 0]
source

ReplicateBE.fixed

ReplicateBE.optstat

ReplicateBE.randrbeds

ReplicateBE.randrbedsFunction
    randrbeds(;n=24, sequence=[1,1],
        design = ["T" "R" "T" "R"; "R" "T" "R" "T"],
        inter=[0.5, 0.4, 0.9], intra=[0.1, 0.2],
        intercept = 0, seqcoef = [0.0, 0.0], periodcoef = [0.0, 0.0, 0.0, 0.0],
        formcoef = [0.0, 0.0],
        dropobs::Int = 0, seed::Int = 0)

Random dataset generation for bioequivalence.

Keywords

  • $n=24$ - number of subjects;
  • $sequence = [1,1]$ - distribution in sequences [1,1] means 1:1, [1,3] - 1:4 etc.;
  • $design = ["T" "R" "T" "R"; "R" "T" "R" "T"]$ - desin matrix, each line is a sequence, each column - periods, cell - formulation id;
  • $inter=[0.5, 0.4, 0.9]$ - inter-subject variance vector for G matrix (length 3): [σ₁, σ₂, ρ], where σ₁, σ₂ - formulation inter-subject variance, ρ - covariance coefficient;
  • $intra=[0.1, 0.2]$ - intra-subject variance vector for R matrix (length 2): [σ₁, σ₂], where σ₁, σ₂ - formulation intra-subject variance for each formulation;
  • $intercept = 0$ - intercept effect value;
  • $seqcoef = [0.0, 0.0]$ - coefficients of sequences, additive (length(sequence) == length(seqcoef) == size(design, 1));
  • $periodcoef = [0.0, 0.0, 0.0, 0.0]$ - coefficients of periods, additive (length(periodcoef) == size(design, 2));
  • $formcoef = [0.0, 0.0]$ - coefficients of formulations, additive ;
  • $dropobs = 0$ number of randomly dropped observations;
  • $seed = 0$ - seed for random number generator (0 - random seed).

Multivariate normal disribution:

\[f(\mathbf{x}; \boldsymbol{\mu}, \boldsymbol{V}) = \frac{1}{(2 \pi)^{d/2} |\boldsymbol{V}|^{1/2}} \exp \left( - \frac{1}{2} (\mathbf{x} - \boldsymbol{\mu})^T V^{-1} (\mathbf{x} - \boldsymbol{\mu}) \right)\]

Where V:

\[V_{i} = Z_{i}GZ_i'+R_{i}\]
source

Another way to use:

    randrbeds(n::Int, sequence::Vector,
        design::Matrix,
        θinter::Vector, θintra::Vector,
        intercept::Real, seqcoef::Vector, periodcoef::Vector, formcoef::Vector,
        dropsubj::Float64, dropobs::Int, seed::Int)
source

Using with RandRBEDS object:

randrbeds(task::RandRBEDS)
source

ReplicateBE.randrbetask

ReplicateBE.randrbetaskFunction
randrbetask(;n=24, sequence=[1,1],
        design = ["T" "R" "T" "R"; "R" "T" "R" "T"],
        inter=[0.5, 0.4, 0.9], intra=[0.1, 0.2],
        intercept = 0, seqcoef = [0.0, 0.0], periodcoef = [0.0, 0.0, 0.0, 0.0], formcoef = [0.0, 0.0],
        dropsubj = 0, dropobs::Int = 0, seed = 0)::RandRBEDS

Make task for random dataset generation.

source

ReplicateBE.reml2

ReplicateBE.reml2Function
reml2(rbe::RBE, θ::Vector{T}) where T <: AbstractFloat

Returm -2logREML for rbe model with θ variance vector.

source
reml2(rbe::RBE)

Returm -2logREML for rbe model

\[logREML(\theta,\beta) = -\frac{N-p}{2} - \frac{1}{2}\sum_{i=1}^nlog|V_{i}|- -\frac{1}{2}log|\sum_{i=1}^nX_i'V_i^{-1}X_i|-\frac{1}{2}\sum_{i=1}^n(y_i - X_{i}\beta)'V_i^{-1}(y_i - X_{i}\beta)\]
source
-2 REML function for ForwardDiff
source

ReplicateBE.simulation

ReplicateBE.simulationFunction
simulation(task::RandRBEDS; io = stdout, verbose = false, num = 100,
    l = log(0.8), u = log(1.25), seed = 0)

Count successful BE outcomes.

Parameters

  • task - RandRBEDS object

Keywords

  • $io = stdout$ - text output
  • $verbose = false$ - print messages to io
  • $num = 100$ - number of simulations
  • $l = log(0.8)$ - lower bound
  • $u = log(1.25)$ - upper bound
  • $seed = 0$ - seed for random number generator (0 - random seed)
source

ReplicateBE.theta

ReplicateBE.thetaFunction
theta(rbe::RBE)

Return theta (θ) vector (vector of variation parameters from optimization procedure).

source

ReplicateBE.typeiii