Title: | Mathematical Modelling of (Dynamic) Microbial Inactivation |
---|---|
Description: | Functions for modelling microbial inactivation under isothermal or dynamic conditions. The calculations are based on several mathematical models broadly used by the scientific community and industry. Functions enable to make predictions for cases where the kinetic parameters are known. It also implements functions for parameter estimation for isothermal and dynamic conditions. The model fitting capabilities include an Adaptive Monte Carlo method for a Bayesian approach to parameter estimation. |
Authors: | Alberto Garre [aut, cre], Pablo S. Fernandez [aut], Jose A. Egea [aut] |
Maintainer: | Alberto Garre <[email protected]> |
License: | GPL-3 |
Version: | 1.3.0.9999 |
Built: | 2024-10-27 05:05:25 UTC |
Source: | https://github.com/albgarre/bioinactivation |
Returns the predicted logarithmic reduction in microbial count according to Arrhenius model for the time, temperature and model parameters given.
Arrhenius_iso(time, temp, k_ref, temp_ref, Ea)
Arrhenius_iso(time, temp, k_ref, temp_ref, Ea)
time |
numeric indicating the time at which the prediction is taken. |
temp |
numeric indicating the temperature of the treatment. |
k_ref |
numeric indicating the inactivation rate at the reference temperature. |
temp_ref |
numeric indicating the reference temperature. |
Ea |
numeric indicating the activation energy. |
A numeric with the predicted logarithmic reduction
().
Returns the predicted logarithmic reduction in microbial count according to Bigelow's model for the time, temperature and model parameters given.
Bigelow_iso(time, temp, z, D_R, temp_ref)
Bigelow_iso(time, temp, z, D_R, temp_ref)
time |
numeric indicating the time at which the prediction is taken. |
temp |
numeric indicating the temperature of the treatment. |
z |
numeric defining the z-value. |
D_R |
numeric defining the D-value at the reference temperature. |
temp_ref |
numeric defining the reference temperature. |
A numeric with the predicted logarithmic reduction
().
Builds an interpolator of the temperature at each time and its first
derivative.
First derivatives are approximated using forward finite
differences (approxfun
). It is assumed that temperature is 0
and constant outside the time interval provided.
build_temperature_interpolator(temperature_data)
build_temperature_interpolator(temperature_data)
temperature_data |
data frame with the values of the temperatures at each value of time. It need to have 2 columns, named time and temperature. |
a list with with two elements:
temp, the interpolator of the temperature and
dtemp, the interpolator of its first derivative
Makes 3 checks of compatibility between the input parameters for the adjustment and the DOF of the inactivation model selected.
Check of equal length of model DOF and input DOF. Raises a warning if failed.
Check that every single one of the input DOF is a model DOF. Raises a warning if failed.
Check that every single one of the model DOF are defined as an input DOF.
check_model_params( simulation_model, known_params, starting_points, adjust_vars )
check_model_params( simulation_model, known_params, starting_points, adjust_vars )
simulation_model |
character with a valid model key. |
known_params |
named vector or list with the dof of the model known. |
starting_points |
named vector or list with the dof of the model to be adjusted. |
adjust_vars |
logical specifying whether the model variables need to be included in the check (not used for isothermal fit) |
Calculates the first derivative of the Arrhenius model with log-linear inactivation for dynamic problems at a given time for the model parameters provided and the environmental conditions given.
dArrhenius_model(t, x, parms, temp_profile)
dArrhenius_model(t, x, parms, temp_profile)
t |
numeric vector indicating the time of the experiment. |
x |
list with the value of N at t. |
parms |
parameters for the secondary model. No explicit check of their validity is performed. |
temp_profile |
a function that provides the temperature at a given time. |
This function is compatible with the function
predict_inactivation
.
The value of the first derivative of at time
t
as a
list.
temp_ref: Reference temperature for the calculation,
k_ref: inactivation rate at the ref. temp.
Ea: Activation energy.
Calculates the first derivative of the linearized version of Bigelow's model for dynamic problems at a given time for the model parameters provided and the environmental conditions given.
dBigelow_model(t, x, parms, temp_profile)
dBigelow_model(t, x, parms, temp_profile)
t |
numeric vector indicating the time of the experiment. |
x |
list with the value of N at t. |
parms |
parameters for the secondary model. No explicit check of their validity is performed. |
temp_profile |
a function that provides the temperature at a given time. |
The model is developed from the isothermal Bigelow's model assuming during
the derivation process that is time independenent.
This function is compatible with the function
predict_inactivation
.
The value of the first derivative of at time
t
as a
list.
temp_ref: Reference temperature for the calculation,
D_R: D-value at the reference temperature,
z: z value.
Calculates the first derivative of the Geeraerd's model at a given time for the model parameters provided and the environmental conditions given.
dGeeraerd_model(t, x, parms, temp_profile)
dGeeraerd_model(t, x, parms, temp_profile)
t |
numeric vector indicating the time of the experiment. |
x |
list with the values of the variables at time |
parms |
parameters for the secondary model. No explicit check of their validity is performed (see section Model Parameters). |
temp_profile |
a function that provides the temperature at a given time. |
This function is compatible with the function
predict_inactivation
.
A list with the value of the first derivatives of N and C_c at time
t
.
temp_ref: Reference temperature for the calculation,
D_R: D-value at the reference temperature,
z: z value,
N_min: Minimum value of N (defines the tail).
To define the Geeraerd model without tail, assign N_min = 0
.
For the model without shoulder, assign C_0 = 0
Calculates the first derivative of Weibull-Mafart model at a given time for the model parameters provided and the environmental conditions given.
dMafart_model(t, x, parms, temp_profile)
dMafart_model(t, x, parms, temp_profile)
t |
numeric vector indicating the time of the experiment. |
x |
list with the value of N at t. |
parms |
parameters for the secondary model. No explicit check of their validity is performed (see section Model Parameters). |
temp_profile |
a function that provides the temperature at a given time. |
The model is developed from the isothermal Weibull-Mafart model without
taking into
account in the derivation the time dependence of for
non-isothermal temperature profiles.
This function is compatible with the function
predict_inactivation
.
The value of the first derivative of N at time t
as a list.
temp_ref: Reference temperature for the calculation.
delta_ref: Value of the scale factor at the reference temperature.
z: z-value.
p: shape factor of the Weibull distribution.
For t=0, dN = 0 unless n=1. Hence, a small shift needs to be introduced to t.
Calculates the first derivative of Metselaar model at a given time for the model parameters provided and the environmental conditions given.
dMetselaar_model(t, x, parms, temp_profile)
dMetselaar_model(t, x, parms, temp_profile)
t |
numeric vector indicating the time of the experiment. |
x |
list with the value of N at t. |
parms |
parameters for the secondary model. No explicit check of their validity is performed (see section Model Parameters). |
temp_profile |
a function that provides the temperature at a given time. |
The model is developed from the isothermal Metselaar model without
taking into
account in the derivation the time dependence of for
non-isothermal temperature profiles.
This function is compatible with the function
predict_inactivation
.
The value of the first derivative of N at time t
as a list.
temp_ref: Reference temperature for the calculation.
D_R: D-value at the reference temperature.
z: z-value.
p: shape factor of the Weibull distribution.
Delta: Scaling parameter
For t=0, dN = 0 unless n=1. Hence, a small shift needs to be introduced to t.
Calculates the first derivative of Weibull-Peleg model at a given time for the model parameters provided and the environmental conditions given.
dPeleg_model(t, x, parms, temp_profile)
dPeleg_model(t, x, parms, temp_profile)
t |
numeric vector indicating the time of the experiment. |
x |
list with the value of logS at t. |
parms |
parameters for the secondary model. No explicit check of their validity is performed (see section Model Parameters). |
temp_profile |
a function that provides the temperature at a given time. |
The model is developed from the isothermal Weibull model without
taking into
account in the derivation the time dependence of for
non-isothermal temperature profiles.
This function is compatible with the function
predict_inactivation
.
The value of the first derivative of logS at time t
as a list.
temp_crit: Temperature below which there is inactivation.
k_b: slope of the b ~ temp line for temperatures above the critical one.
n: shape factor of the Weibull distribution.
For logS=0, dlogS = 0 unless n=1. Hence, a small shift needs to be introduced to logS.
Example of experimental data of the dynamic inactivation process of a microorganism.
data(dynamic_inactivation)
data(dynamic_inactivation)
A data frame with 19 rows and 2 variables.
time: Time in minutes of the measurement.
N: Number of microorganism.
temperature: Observed temperature.
Fits the parameters of an inactivation model to experimental data.
The function modFit
of the package FME
is
used for the adjustment.
fit_dynamic_inactivation( experiment_data, simulation_model, temp_profile, starting_points, upper_bounds, lower_bounds, known_params, ..., minimize_log = TRUE, tol0 = 1e-05 )
fit_dynamic_inactivation( experiment_data, simulation_model, temp_profile, starting_points, upper_bounds, lower_bounds, known_params, ..., minimize_log = TRUE, tol0 = 1e-05 )
experiment_data |
data frame with the experimental data to be adjusted. It must have a column named “time” and another one named “N”. |
simulation_model |
character identifying the model to be used. |
temp_profile |
data frame with discrete values of the temperature for
each time. It must have one column named |
starting_points |
starting values of the parameters for the adjustment. |
upper_bounds |
named numerical vector defining the upper bounds of the parameters for the adjustment. |
lower_bounds |
named numerical vector defining the lower bounds of the parameters for the adjustment. |
known_params |
named numerical vector with the fixed (i.e., not adjustable) model parameters. |
... |
further arguments passed to |
minimize_log |
logical. If |
tol0 |
numeric. Observations at time 0 make Weibull-based models singular. The time for observatins taken at time 0 are changed for this value. |
A list of class FitInactivation
with the following items:
fit_results: a list of class modFit
with the info
of the adjustment.
best_prediction: a list of class SimulInactivation
with prediction made by the adjusted model.
data: a data frame with the data used for the fitting.
## EXAMPLE 1 ------ data(dynamic_inactivation) # The example data set is used. get_model_data() # Retrieve the valid model keys. simulation_model <- "Peleg" # Peleg's model will be used model_data <- get_model_data(simulation_model) model_data$parameters # Set the model parameters dummy_temp <- data.frame(time = c(0, 1.25, 2.25, 4.6), temperature = c(70, 105, 105, 70)) # Dummy temp. profile ## Set known parameters and initial points/bounds for unknown ones known_params = c(temp_crit = 100) starting_points <- c(n = 1, k_b = 0.25, N0 = 1e+05) upper_bounds <- c(n = 2, k_b = 1, N0 = Inf) lower_bounds <- c(n = 0, k_b = 0, N0 = 1e4) dynamic_fit <- fit_dynamic_inactivation(dynamic_inactivation, simulation_model, dummy_temp, starting_points, upper_bounds, lower_bounds, known_params) plot(dynamic_fit) goodness_of_fit(dynamic_fit) ## END EXAMPLE 1 -----
## EXAMPLE 1 ------ data(dynamic_inactivation) # The example data set is used. get_model_data() # Retrieve the valid model keys. simulation_model <- "Peleg" # Peleg's model will be used model_data <- get_model_data(simulation_model) model_data$parameters # Set the model parameters dummy_temp <- data.frame(time = c(0, 1.25, 2.25, 4.6), temperature = c(70, 105, 105, 70)) # Dummy temp. profile ## Set known parameters and initial points/bounds for unknown ones known_params = c(temp_crit = 100) starting_points <- c(n = 1, k_b = 0.25, N0 = 1e+05) upper_bounds <- c(n = 2, k_b = 1, N0 = Inf) lower_bounds <- c(n = 0, k_b = 0, N0 = 1e4) dynamic_fit <- fit_dynamic_inactivation(dynamic_inactivation, simulation_model, dummy_temp, starting_points, upper_bounds, lower_bounds, known_params) plot(dynamic_fit) goodness_of_fit(dynamic_fit) ## END EXAMPLE 1 -----
Fits the parameters of an inactivation model to experimental using
the Markov Chain Monte Carlo fitting algorithm implemented in
the function modMCMC
of the package FME
.
fit_inactivation_MCMC( experiment_data, simulation_model, temp_profile, starting_points, upper_bounds, lower_bounds, known_params, ..., minimize_log = TRUE, tol0 = 1e-05 )
fit_inactivation_MCMC( experiment_data, simulation_model, temp_profile, starting_points, upper_bounds, lower_bounds, known_params, ..., minimize_log = TRUE, tol0 = 1e-05 )
experiment_data |
data frame with the experimental data to be adjusted. It must have a column named “time” and another one named “N”. |
simulation_model |
character identifying the model to be used. |
temp_profile |
data frame with discrete values of the temperature for
each time. It must have one column named |
starting_points |
starting values of the parameters for the adjustment. |
upper_bounds |
named numerical vector defining the upper bounds of the parameters for the adjustment. |
lower_bounds |
named numerical vector defining the lower bounds of the parameters for the adjustment. |
known_params |
named numerical vector with the fixed (i.e., not adjustable) model parameters. |
... |
other arguments for |
minimize_log |
logical. If |
tol0 |
numeric. Observations at time 0 make Weibull-based models singular. The time for observatins taken at time 0 are changed for this value. |
A list of class FitInactivationMCMC
with the following items:
modMCMC
: a list of class modMCMC
with the
information of the adjustment process.
best_prediction: a list of class SimulInactivation
with the prediction generated by the best predictor.
data: a data frame with the data used for the fitting.
## EXAMPLE 1 ------ data(dynamic_inactivation) # The example data set is used. get_model_data() # Retrieve the valid model keys. simulation_model <- "Peleg" # Peleg's model will be used model_data <- get_model_data(simulation_model) model_data$parameters # Set the model parameters dummy_temp <- data.frame(time = c(0, 1.25, 2.25, 4.6), temperature = c(70, 105, 105, 70)) # Dummy temp. profile ## Set known parameters and initial points/bounds for unknown ones known_params = c(temp_crit = 100) starting_points <- c(n = 1, k_b = 0.25, N0 = 1e+05) upper_bounds <- c(n = 2, k_b = 1, N0 = 1e6) lower_bounds <- c(n = 0.5, k_b = 0.1, N0 = 1e4) MCMC_fit <- fit_inactivation_MCMC(dynamic_inactivation, simulation_model, dummy_temp, starting_points, upper_bounds, lower_bounds, known_params, niter = 100) # It is recommended to increase niter plot(MCMC_fit) goodness_of_fit(MCMC_fit) ## END EXAMPLE 1 -----
## EXAMPLE 1 ------ data(dynamic_inactivation) # The example data set is used. get_model_data() # Retrieve the valid model keys. simulation_model <- "Peleg" # Peleg's model will be used model_data <- get_model_data(simulation_model) model_data$parameters # Set the model parameters dummy_temp <- data.frame(time = c(0, 1.25, 2.25, 4.6), temperature = c(70, 105, 105, 70)) # Dummy temp. profile ## Set known parameters and initial points/bounds for unknown ones known_params = c(temp_crit = 100) starting_points <- c(n = 1, k_b = 0.25, N0 = 1e+05) upper_bounds <- c(n = 2, k_b = 1, N0 = 1e6) lower_bounds <- c(n = 0.5, k_b = 0.1, N0 = 1e4) MCMC_fit <- fit_inactivation_MCMC(dynamic_inactivation, simulation_model, dummy_temp, starting_points, upper_bounds, lower_bounds, known_params, niter = 100) # It is recommended to increase niter plot(MCMC_fit) goodness_of_fit(MCMC_fit) ## END EXAMPLE 1 -----
Fits the parameters of the model chosen to a set of isothermal experiments
using nonlinear regression through the function nls
.
fit_isothermal_inactivation( model_name, death_data, starting_point, known_params, adjust_log = TRUE )
fit_isothermal_inactivation( model_name, death_data, starting_point, known_params, adjust_log = TRUE )
model_name |
character specyfing the model to adjust. |
death_data |
data frame with the experiment data where each row is one observation. It must have the following columns:
|
starting_point |
List with the initial values of the parameters for the adjustment. |
known_params |
List of the parameters of the model known. |
adjust_log |
logical. If |
An instance of class IsoFitInactivation
with the results.
This list has four entries:
nls: The object of class nls
with the results
of the adjustment.
parameters: a list with the values of the model parameters, both fixed and adjusted.
model: a string with the key identifying the model used.
data: the inactivation data used for the fit.
## EXAMPLE 1 ----------- data("isothermal_inactivation") # data set used for the example. get_isothermal_model_data() # retrieve valid model keys. model_name <- "Bigelow" # Bigelow's model will be used for the adjustment. model_data <- get_isothermal_model_data(model_name) model_data$params # Get the parameters of the model ## Define the input arguments known_params = list(temp_ref = 100) starting_point <- c(z = 10,D_R = 1) ## Call the fitting function iso_fit <- fit_isothermal_inactivation(model_name, isothermal_inactivation, starting_point, known_params) ## Output of the results plot(iso_fit, make_gg = TRUE) goodness_of_fit(iso_fit) ## END EXAMPLE 1 --------
## EXAMPLE 1 ----------- data("isothermal_inactivation") # data set used for the example. get_isothermal_model_data() # retrieve valid model keys. model_name <- "Bigelow" # Bigelow's model will be used for the adjustment. model_data <- get_isothermal_model_data(model_name) model_data$params # Get the parameters of the model ## Define the input arguments known_params = list(temp_ref = 100) starting_point <- c(z = 10,D_R = 1) ## Call the fitting function iso_fit <- fit_isothermal_inactivation(model_name, isothermal_inactivation, starting_point, known_params) ## Output of the results plot(iso_fit, make_gg = TRUE) goodness_of_fit(iso_fit) ## END EXAMPLE 1 --------
Provides information of the models implemented for fitting of isothermal
data.
This models are valid only for isothermal adjustment with the function
fit_isothermal_inactivation
. To make predictions with the
function predict_inactivation
or adjust dynamic experiments
with fit_dynamic_inactivation
, use
get_model_data
.
get_isothermal_model_data(model_name = "valids")
get_isothermal_model_data(model_name = "valids")
model_name |
Optional string with the key of the model to use. |
If model_name
is missing, a list of the valid model keys.
If model_name
is not a valid key, NULL is returned.
Otherwise, a list with the parameters of the model selected and its
formula
for the nonlinear adjustment.
Provides information about the function for dynamic predictions associated
to a valid simulation_model
key.
If simulation_model
is missing or NULL
, a character vector
of valid model keys is provided.
This function is designed as an assistant for using the functions
predict_inactivation
and
fit_dynamic_inactivation
.
For the adjustment of isothermal experiments with the function
fit_isothermal_inactivation
, use the function
get_isothermal_model_data
.
get_model_data(simulation_model = NULL)
get_model_data(simulation_model = NULL)
simulation_model |
(optional) character with a valid model key or
|
If simulation_model is NULL
or missing, a character vector of
possible names. Otherwise, a list including information of the relevant
function:
ode: Pointer to the function defining the model ode.
cost: Pointer to the function calculating the error of the approximation.
dtemp: logical defining whether the function requires the definition of the first derivative of temperature.
variables: a character vector defining which entry variables are needed by the model.
variables_priv: for internal use only.
parameters: character vector with the parameters needed by the model.
predict_inactivation
,
fit_dynamic_inactivation
Calculates the error of the prediction of microbial inactivation for
the chosen inactivation model and the given parameters with respect to
the experimental data provided.
This function is compatible with the function
fit_dynamic_inactivation
.
get_prediction_cost( data_for_fit, temp_profile, simulation_model, P, known_params )
get_prediction_cost( data_for_fit, temp_profile, simulation_model, P, known_params )
data_for_fit |
A data frame with the experimental data to fit. It must contain a column named “time” and another one named “N”. |
temp_profile |
|
simulation_model |
character key defining the inactivation model. |
P |
list with the unknown parameters of the model to be adjusted. |
known_params |
list with the parameters of the model fixed (i.e., not adjusted) |
An instance of modCost
with the error of the
prediction.
modCost
, fit_dynamic_inactivation
Goodness of fit for Dynamic fits
goodness_dyna(object)
goodness_dyna(object)
object |
An instance of FitInactivation |
Goodness of fit for Isothermal fits
goodness_iso(object)
goodness_iso(object)
object |
An object of class IsoFitInactivation. |
Goodness of fit for MCMC fits
goodness_MCMC(object)
goodness_MCMC(object)
object |
An instance of FitInactivationMCMC |
Generates a table with several statistical indexes describing the quality of the model fit:
ME: Mean Error.
RMSE: Root Mean Squared Error.
loglik: log-likelihood.
AIC: Akaike Information Criterion.
AICc: Akaike Information Criterion with correction for finite sample size.
BIC: Bayesian Infromation Criterion.
Af: Accuracy factor.
Bf: Bias factor.
goodness_of_fit(object)
goodness_of_fit(object)
object |
A model fit generated by bioinactivation |
Tests if an object is of class FitInactivation
.
is.FitInactivation(x)
is.FitInactivation(x)
x |
object to be checked. |
A logic specifying whether x
is of class
FitInactivation
Tests if an object is of class FitInactivationMCMC
.
is.FitInactivationMCMC(x)
is.FitInactivationMCMC(x)
x |
object to be checked. |
A logic specifying whether x
is of class
FitInactivationMCMC
Tests if an object is of class IsoFitInactivation
.
is.IsoFitInactivation(x)
is.IsoFitInactivation(x)
x |
object to be checked. |
A logic specifying whether x
is of class
IsoFitInactivation
Tests if an object is of class PredInactivationMCMC
.
is.PredInactivationMCMC(x)
is.PredInactivationMCMC(x)
x |
object to be checked. |
A logic specifying whether x
is of class
PredInactivationMCMC
Tests if an object is of class SimulInactivation
.
is.SimulInactivation(x)
is.SimulInactivation(x)
x |
object to be checked. |
A logic specifying whether x
is of class
SimulInactivation
Example of experimental data for an isothermal process of a microorganism.
data(isothermal_inactivation)
data(isothermal_inactivation)
A data frame with 36 rows and 3 variables.
time: Time in minutes of the measurement.
temp: Temperature at which the experiment was made.
log_diff: Logarithmic difference.
Example of experimental data of the dynamic inactivation process of Laterosporus
data(laterosporus_dyna)
data(laterosporus_dyna)
A data frame with 20 rows and 3 variables.
time: Time in minutes of the measurement.
temp: observed temperature.
logN: recorded number of microorganism.
Example of experimental data for an isothermal process of Laterosporus.
data(laterosporus_iso)
data(laterosporus_iso)
A data frame with 52 rows and 3 variables.
time: Time in minutes of the measurement.
temp: Temperature at which the experiment was made.
log_diff: Logarithmic difference.
Returns the predicted logarithmic reduction in microbial count according to Metselaars's model for the time, temperature and model parameters given.
Metselaar_iso(time, temp, D_R, z, p, Delta, temp_ref)
Metselaar_iso(time, temp, D_R, z, p, Delta, temp_ref)
time |
numeric indicating the time at which the prediction is taken. |
temp |
numeric indicating the temperature of the treatment. |
D_R |
numeric defining the delta-value at the reference temperature. |
z |
numeric defining the z-value. |
p |
numeric defining shape factor of the Weibull distribution. |
Delta |
numeric reparametrization factor |
temp_ref |
numeric indicating the reference temperature. |
A numeric with the predicted logarithmic reduction
().
Plots a comparison between the experimental data provided and the prediction
produced by the model parameters adjusted for an instance of
FitInactivation
.
## S3 method for class 'FitInactivation' plot( x, y = NULL, ..., make_gg = TRUE, plot_temp = FALSE, label_y1 = "logN", label_y2 = "Temperature", ylims = NULL )
## S3 method for class 'FitInactivation' plot( x, y = NULL, ..., make_gg = TRUE, plot_temp = FALSE, label_y1 = "logN", label_y2 = "Temperature", ylims = NULL )
x |
the object of class |
y |
ignored |
... |
additional arguments passed to |
make_gg |
logical. If |
plot_temp |
logical. Whether the temperature profile will
be added to the plot. |
label_y1 |
Label of the principal y-axis. |
label_y2 |
Label of the secondary y-axis. |
ylims |
Numeric vector of length 2 with the Limits of the
y-axis. |
If make_gg = FALSE
, the plot is created. Otherwise, an
an instance of ggplot
is generated, printed and returned.
Plots a comparison between the experimental data provided and the prediction
produced by the model parameters adjusted for an instance of
FitInactivationMCMC
.
## S3 method for class 'FitInactivationMCMC' plot( x, y = NULL, ..., make_gg = TRUE, plot_temp = FALSE, label_y1 = "logN", label_y2 = "Temperature", ylims = NULL )
## S3 method for class 'FitInactivationMCMC' plot( x, y = NULL, ..., make_gg = TRUE, plot_temp = FALSE, label_y1 = "logN", label_y2 = "Temperature", ylims = NULL )
x |
the object of class |
y |
ignored |
... |
additional arguments passed to |
make_gg |
logical. If |
plot_temp |
logical. Whether the temperature profile will
be added to the plot. |
label_y1 |
Label of the principal y-axis. |
label_y2 |
Label of the secondary y-axis. |
ylims |
Numeric vector of length 2 with the Limits of the
y-axis. |
If make_gg = FALSE
, the plot is created. Otherwise, an
an instance of ggplot
is generated, printed and returned.
For each one of the temperatures studied, plots a comparison between the
predicted result and the experimental one for an instance of
IsoFitInactivation
.
## S3 method for class 'IsoFitInactivation' plot(x, y = NULL, ..., make_gg = FALSE)
## S3 method for class 'IsoFitInactivation' plot(x, y = NULL, ..., make_gg = FALSE)
x |
the object of class |
y |
ignored |
... |
additional arguments passed to |
make_gg |
logical. If |
Plots the prediction interval generated by
predict_inactivation_MCMC
.
## S3 method for class 'PredInactivationMCMC' plot(x, y = NULL, ..., make_gg = TRUE)
## S3 method for class 'PredInactivationMCMC' plot(x, y = NULL, ..., make_gg = TRUE)
x |
the object of class |
y |
ignored |
... |
additional arguments passed to |
make_gg |
logical. If |
The plot generated in ggplot (default) generates a dashed line with the mean of the MC simulations. Moreover, a ribbon with the 2 first quantiles (i.e. columns 3 and 4) is generated.
The plot generated with base R (make_gg = FALSE
) generates a solid line
wth the mean of the MC simulations. Each one of the other quantiles included
in the data frame are added with different
If make_gg = FALSE
, the plot is created. Otherwise, an
an instance of ggplot
is generated, printed and returned.
Plots the predicted evolution of the logarithmic count with time for an
instance of SimulInactivation
.
## S3 method for class 'SimulInactivation' plot( x, y = NULL, ..., make_gg = TRUE, plot_temp = FALSE, label_y1 = "logN", label_y2 = "Temperature", ylims = NULL )
## S3 method for class 'SimulInactivation' plot( x, y = NULL, ..., make_gg = TRUE, plot_temp = FALSE, label_y1 = "logN", label_y2 = "Temperature", ylims = NULL )
x |
The object of class |
y |
ignored |
... |
additional arguments passed to |
make_gg |
logical. If |
plot_temp |
logical. Whether the temperature profile will
be added to the plot. |
label_y1 |
Label of the principal y-axis. |
label_y2 |
Label of the secondary y-axis. |
ylims |
Numeric vector of length 2 with the Limits of the
y-axis. |
If make_gg = FALSE
, the plot is created. Otherwise, an
an instance of ggplot
is generated, printed and returned.
Predicts the inactivation of a microorganism under isothermal or non-isothermal temperature conditions. The thermal resistence of the microorganism are defined with the input arguments.
predict_inactivation( simulation_model, times, parms, temp_profile, ..., tol0 = 1e-05 )
predict_inactivation( simulation_model, times, parms, temp_profile, ..., tol0 = 1e-05 )
simulation_model |
character identifying the model to be used. |
times |
numeric vector of output times. |
parms |
list of parameters defining the parameters of the model. |
temp_profile |
data frame with discrete values of the temperature for
each time. It must have one column named |
... |
Additional arguments passed to |
tol0 |
numeric. Observations at time 0 make Weibull-based models singular. The time for observatins taken at time 0 are changed for this value. By default ('tol0 = 1e-5') |
The value of the temperature is calculated at each value of time by
linear interpolation of the values provided by the input argument
temp_profile
.
The function ode
of the package deSolve
is
used for the resolution of the differential equation.
A list of class SimulInactivation
with the results. It has
the following entries:
model: character defining the model use for the prediction.
model_parameters: named numeric vector with the values of the model parameters used.
temp_approximations: function used for the interpolation of the temperature. For a numeric value of time given, returns the value of the temperature and its first derivative.
simulation: A data frame with the results calculated. Its
first column contains the times at which the
solution has been calculated. The following
columns the values of the variables of the
model. The three last columns provide the
values of logN
, S
and
logS
.
## EXAMPLE 1 ----------- ## Retrieve the model keys available for dynamic models. get_model_data() ## Set the input arguments example_model <- "Geeraerd" # Geeraerd's model will be used times <- seq(0, 5, length=100) # values of time for output model_data <- get_model_data(example_model) # Retrive the data of the model used print(model_data$parameters) print(model_data$variables) model_parms <- c(D_R = 1, z = 10, N_min = 100, temp_ref = 100, N0 = 100000, C_c0 = 1000 ) ## Define the temperature profile for the prediction temperature_profile <- data.frame(time = c(0, 1.25, 2.25, 4.6), temperature = c(70, 105, 105, 70)) ## Call the prediction function prediction_results <- predict_inactivation(example_model, times, model_parms, temperature_profile) ## Show the results head(prediction_results$simulation) plot(prediction_results) time_to_logreduction(1.5, prediction_results) ## END EXAMPLE 1 -----------
## EXAMPLE 1 ----------- ## Retrieve the model keys available for dynamic models. get_model_data() ## Set the input arguments example_model <- "Geeraerd" # Geeraerd's model will be used times <- seq(0, 5, length=100) # values of time for output model_data <- get_model_data(example_model) # Retrive the data of the model used print(model_data$parameters) print(model_data$variables) model_parms <- c(D_R = 1, z = 10, N_min = 100, temp_ref = 100, N0 = 100000, C_c0 = 1000 ) ## Define the temperature profile for the prediction temperature_profile <- data.frame(time = c(0, 1.25, 2.25, 4.6), temperature = c(70, 105, 105, 70)) ## Call the prediction function prediction_results <- predict_inactivation(example_model, times, model_parms, temperature_profile) ## Show the results head(prediction_results$simulation) plot(prediction_results) time_to_logreduction(1.5, prediction_results) ## END EXAMPLE 1 -----------
Given a model adjustment of a dynamic microbial inactivation process
performed using any of the functions in bioinactivation
calculates
probability intervals at each time point using a Monte Carlo method.
predict_inactivation_MCMC( fit_object, temp_profile, n_simulations = 100, times = NULL, quantiles = c(2.5, 97.5), additional_pars = NULL )
predict_inactivation_MCMC( fit_object, temp_profile, n_simulations = 100, times = NULL, quantiles = c(2.5, 97.5), additional_pars = NULL )
fit_object |
An object of classes |
temp_profile |
data frame with discrete values of the temperature for
each time. It must have one column named |
n_simulations |
a numeric indicating how many Monte Carlo simulations
to perform. |
times |
numeric vector specifying the time points when results are
desired. If |
quantiles |
numeric vector indicating the quantiles to calculate in
percentage. By default, it is set to c(2.5, 97.5) which generates a
prediction interval with confidence 0.95. If |
additional_pars |
Additional parameters not included in the adjustment (e.g. the initial number of microorganism in an isothermal fit). |
A data frame of class PredInactivationMCMC
. On its first column,
time at which the calculation has been made is indicated.
If quantiles = NULL
, the following columns contain the
results of each simulation. Otherwise, the second and third columns
provide the mean and median of the simulations at the given time
point. Following columns contain the quantiles of the results.
Function to be called by predict_inactivation_MCMC
. Produces
a random sample of the parameters adjusted from dynamic experiments with non
linear regression.
sample_dynaFit(dynamic_fit, times, n_simulations)
sample_dynaFit(dynamic_fit, times, n_simulations)
dynamic_fit |
An object of class |
times |
numeric vector specifying the time points when results are
desired. If |
n_simulations |
a numeric indicating how many Monte Carlo simulations to perform. |
It is assumed that the parameters follow a Multivariate Normal distribution
with the mean the parameters estimated by modFit
. The unscaled
covariance matrix returned by modFit
is used.
The function produces a random sample using the function mvrnorm
.
Returns a list with 4 components.
par_sample: data frame with the parameter sample.
simulation_model: model key for the simulation
known_pars: parameters of the model known
times: points where the calculation is made.
Function to be called by predict_inactivation_MCMC
. Produces
a random sample of the parameters adjusted from isothermal experiments.
sample_IsoFit(iso_fit, times, n_simulations)
sample_IsoFit(iso_fit, times, n_simulations)
iso_fit |
An object of class |
times |
numeric vector specifying the time points when results are
desired. If |
n_simulations |
a numeric indicating how many Monte Carlo simulations to perform. |
It is assumed that the parameters follow a Multivariate Normal distribution
with the mean and covariance matrix estimated by the adjustment. The function
produces a random sample using the function mvrnorm
.
Returns a list with 4 components.
par_sample: data frame with the parameter sample.
simulation_model: model key for the simulation
known_pars: parameters of the model known
times: points where the calculation is made.
Function to be called by predict_inactivation_MCMC
. Produces
a random sample of the parameters calculated on the iterations of the Monte
Carlo simulation.
sample_MCMCfit(MCMC_fit, times, n_simulations)
sample_MCMCfit(MCMC_fit, times, n_simulations)
MCMC_fit |
An object of class |
times |
numeric vector specifying the time points when results are
desired. If |
n_simulations |
a numeric indicating how many Monte Carlo simulations to perform. |
Returns a list with 4 components.
par_sample: data frame with the parameter sample.
simulation_model: model key for the simulation
known_pars: parameters of the model known
times: points where the calculation is made.
Summary of a FitInactivation object
## S3 method for class 'FitInactivation' summary(object, ...)
## S3 method for class 'FitInactivation' summary(object, ...)
object |
Instance of Fit Inactivation |
... |
ignored |
Summary of a FitInactivationMCMC object
## S3 method for class 'FitInactivationMCMC' summary(object, ...)
## S3 method for class 'FitInactivationMCMC' summary(object, ...)
object |
Instance of FitInactivationMCMC |
... |
ignored |
Summary of a IsoFitInactivation object
## S3 method for class 'IsoFitInactivation' summary(object, ...)
## S3 method for class 'IsoFitInactivation' summary(object, ...)
object |
Instance of IsoFitInactivation |
... |
ignored |
Calculates the treatment time required to reach a given number of log reductions.
time_to_logreduction(n_logs, my_prediction)
time_to_logreduction(n_logs, my_prediction)
n_logs |
Numeric of length one indicating the number of log recutions |
my_prediction |
An object of class SimulInactivation |
The treatement time is calculated by linear interpolation betweent the two points of the simulation whose logS is closer to n_logs
Returns the predicted logarithmic reduction in microbial count according to Weibull-Mafarts's model for the time, temperature and model parameters given.
WeibullMafart_iso(time, temp, delta_ref, z, p, temp_ref)
WeibullMafart_iso(time, temp, delta_ref, z, p, temp_ref)
time |
numeric indicating the time at which the prediction is taken. |
temp |
numeric indicating the temperature of the treatment. |
delta_ref |
numeric defining the delta-value at the reference temperature. |
z |
numeric defining the z-value. |
p |
numeric defining shape factor of the Weibull distribution. |
temp_ref |
numeric indicating the reference temperature. |
A numeric with the predicted logarithmic reduction
().
Returns the predicted logarithmic reduction in microbial count according to Weibull-Peleg's model for the time, temperature and model parameters given.
WeibullPeleg_iso(time, temp, n, k_b, temp_crit)
WeibullPeleg_iso(time, temp, n, k_b, temp_crit)
time |
numeric indicating the time at which the prediction is taken. |
temp |
numeric indicating the temperature of the treatment. |
n |
numeric defining shape factor of the Weibull distribution. |
k_b |
numeric indicating the slope of the b~temp line. |
temp_crit |
numeric with the value of the critical temperature. |
A numeric with the predicted logarithmic reduction
().