Competition Models¶
The competitions
module contains functions for predicting
the results of competitions between different species or strains
and analysing the competitions results to extract the fitness
or selection coefficients of the competing strains.
The model contains:
Competition function (
compete()
) performs a numerical simulation of a competition.Competition models (
dobule_*_ode*
) are Python implementations of ODE systems, consistent with Scipy’sodeint
requirement. These ODE systems define competition models in which two strains/species compete for resources while growing in a well-mixed homogeneous vessel.Fitness and selection coefficient analysis (
fitness_*
,selection_coefs_*
) are functions for infering fitness and selection coefficients from competitions time-series.
How competitions work¶
Consider the following double logistic competition model in the form of a system of ordinary differential equations:
\(N_i\): population size of strain i.
\(r_i\): initial per capita growth rate of strain i
\(K_i\): maximum population size of strain i
To simulate a competition between two strains we need to set the values of
\(N_i(0), r_i, K_i\) for \(i \in {1,2}\)
where 1 and 2 define two strains.
These values can be estimated from growth curve data using curveball.models.fit_model()
,
which fits growth models to data and returns model fit results.
These model fit results can be passed to the compete()
function to simulate a competition:
>>> import pandas as pd
>>> import curveball
>>> plate = pd.read_csv('plate_templates/G-RG-R.csv')
>>> df = curveball.ioutils.read_tecan_xlsx('data/Tecan_280715.xlsx', label='OD', plate=plate)
>>> green = curveball.models.fit_model(df[df.Strain == 'G'], PLOT=False, PRINT=False)[0]
>>> red = curveball.models.fit_model(df[df.Strain == 'R'], PLOT=False, PRINT=False)[0]
>>> t,y,fig,ax = curveball.competitions.compete(green, red, PLOT=True, colors=['green', 'red'])
The result of compete()
, t
and y
are the time array and the frequencies array of the competition:
specifically, y
is a 2-dimensional array, containing the frequency of each strain at each time point.
Members¶
-
curveball.competitions.
baranyi_roberts_gd
(y, t, *args)[source]¶ Fujikawa, Hiroshi, and Mohammad Z Sakha. 2014. Prediction of Microbial Growth in Mixed Culture with a Competition Model. Biocontrol Science 19 (2): 89-92.
-
curveball.competitions.
baranyi_roberts_lv
(y, t, *args)[source]¶ Fujikawa, Hiroshi, and Mohammad Z Sakha. 2014. Prediction of Microbial Growth in Mixed Culture with a Competition Model. Biocontrol Science 19 (2): 89-92.
-
curveball.competitions.
compete
(m1, m2, p0=(0.5, 0.5), y0=None, t=None, hours=24, num_of_points=100, nsamples=1, lag_phase=True, ode=<function double_baranyi_roberts_ode1>, params1=None, params2=None, ci=95, colors=None, ax=None, PLOT=False, sampler='covar', df1=None, df2=None)[source]¶ Simulate competitions between two strains using growth parameters estimated by fitting growth models to growth curves data.
Integrates a 2-dimensional ODE system given by ode with parameters extracted from
lmfit.model.ModelResult
instances m1 and m2. This implementation includes plotting (if required by PLOT); resampling from the distribution of model parameters (when nsamples > 1); changing the ODE system (by providing a different function in ode); and more.The function competes two strains/species; therefore it expects two
lmfit.model.ModelResult
objects, two initial values in y0, etc.- m1, m2lmfit.model.ModelResult
model fitting results of growth models defined in
curveball.models
.- p0, y0tuple, optional
p0 is the initial relative frequencies; y0 is the initial population sizes. if y0 is given than
y0[0]
andy0[0]
are the initial population size for m1 and m2, respectively. if y0 is not given, then it will be set to the average estimated y0 in m1 and m2 multiplied by p0.- tnumpy.ndarray or None, optional
array of time points in which to calculate the population sizes.
- hoursint, optional
if t is not given, determines how many hours should the competition proceed, defaults to 24.
- num_of_pointsint, optional
if t is not given, determines the number of time points to use, defaults to 100.
- nsamplesint, optional
how many replicates of the competition should be simulated; if nsamples = 1, only one competition is simulated with the estimated parameters; otherwise nsamples competitions are simulated with parameters drawn from a parameter distribution determined by sampler. Defaults to 1.
- lag_phasebool, optional
if True, use lag phase as given by m1 and m2. Otherwise, override the lag phase parameters to prevent a lag phase. Defaults to
True
.- odefunc, optional
an ordinary differential systems system defined by a function that accepts
y
,t
, and additional arguments, and returns the derivate ofy
with respect tot
. Defaults todouble_baranyi_roberts_ode0()
.- params1, params2dict, optional
dictionaries of model parameter values; if given, overrides values from m1 and m2.
- cifloat, optional
confidence interval size, in (0, 100), only applicable when PLOT is
True
, defaults to 95%.- colorssequence of str, optional
if PLOT is
True
, this sets the colors of the drawn lines. colors[0] will be used for m1; colors[1] for m2. If not provided, defaults to the current pallete.- axmatplotlib.axes.Axes, optional
if PLOT is
True
, an axes to plot into; if not provided, a new one is created.- PLOTbool, optional
if
True
, the function will plot the curves of y as a function of t. Defaults toFalse
.- samplerstr, optional
if
covar
, the parameters will be sampled using the covariance matrix (curveball.models.sample_params()
); ifbootstrap
the parameters will be sampled by resampling the growth curves (curveball.models.bootstrap_params()
).- df1, df2pandas.DataFrame, optional
the data frames used to fit m1 and m2, only used when sampler is
bootstrap
.
- tnumpy.ndarray
1d (or 2d, if nsamples>1) array of time points, in hours.
- y: numpy.ndarray
2d (or 3d, if nsamples>1) array of strain frequencies. First axis is time, second axis is strain, third axis (if applicable) is sample.
- figmatplotlib.figure.Figure
figure object.
- axnumpy.ndarray
array of
matplotlib.axes.Axes
objects.
- TypeError
if m1 or m2 are not
lmfit.model.ModelResult
.- AssertionError
if an intermediate calculation produced an invalid result (not guaranteed).
>>> import pandas as pd >>> import curveball >>> plate = pd.read_csv('plate_templates/G-RG-R.csv') >>> df = curveball.ioutils.read_tecan_xlsx('data/Tecan_280715.xlsx', label='OD', plate=plate) >>> green = curveball.models.fit_model(df[df.Strain == 'G'], PLOT=False, PRINT=False)[0] >>> red = curveball.models.fit_model(df[df.Strain == 'R'], PLOT=False, PRINT=False)[0] >>> t, y = curveball.competitions.compete(green, red, PLOT=False)
To debug, uncomment lines to return the
infodict
that is returned fromscipy.integrate.odeint()
is run withfull_output=1
.
-
curveball.competitions.
double_baranyi_roberts_gimenez_delgado_ode
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model with competition model inspired by Gimenez and Delgado (2004). The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \prod_j{\Big(1 - \Big(\frac{N_j}{K_j}\Big)^{\nu_j}\Big)}\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode0
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model 1. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \Big(1 - \Big(\frac{\sum_{j}{N_j}}{K_i}\Big)^{\nu_i}\Big)\]\(N_i\): population size of strain i.
\(r_i\): initial per capita growth rate of strain i
\(K_i\): maximum population size of strain i
\(\nu_i\): curvature of the logsitic term of strain i
\(\alpha_i(t)= \frac{q_{0,i}}{q_{0,i} + e^{-v_i t}}\)
\(q_{0,i}\): initial adjustment of strain i to current environment
\(v_i\): adjustment rate of strain i
- yfloat, float
population size
- tfloat
time, usually in hours
- rfloat, float
initial per capita growth rate
- Kfloat, float
maximum population size (\(K>0\))
- nufloat, float
curvature of the logsitic term (\(\nu>0\))
- q0float, float
initial adjustment to current environment (\(0<q_0<1\))
- vfloat, float
adjustment rate (\(v>0\))
- float, float
population growth rate.
- 1
Baranyi, J., Roberts, T. A., 1994. A dynamic approach to predicting bacterial growth in food. Int. J. Food Microbiol.
-
curveball.competitions.
double_baranyi_roberts_ode1
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \Big(1 - \Big(\sum_{j}{\frac{N_j}{K_j}}\Big)^{\nu_i}\Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode2
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \Big(1 - \sum_{j}{\Big(\frac{N_j}{K_j}\Big)^{\nu_j}}\Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode3
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \Big( 1 - \Big(\frac{\sum_j{N_j}}{\bar{K}}\Big)^{\nu_i} \Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode4
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode5
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \Big( 1 - \Big(\frac{N_i}{K_i}\Big)^{\nu_i}\Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode6
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i \alpha_i(t) N_i \Big(1 - \Big(\frac{\sum_j{N_j}}{K_i}\Big)^{\bar{\nu}}\Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode7
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_1}{dt} = 0 \;\;\; \frac{dN_2}{dt} = r_2 \alpha_2(t) N_2 \Big(1 - \Big(\frac{N_2}{K_2}\Big)^{\nu_1}\Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
double_baranyi_roberts_ode8
(y, t, r, K, nu, q0, v)[source]¶ A two species Baranyi-Roberts model. The function calculates the population growth rate at given time points.
\[\frac{dN_i}{dt} = r_i N_i \Big(1 - \Big(\frac{\sum_{j}{N_j}}{K_i}\Big)^{\nu_i}\Big)\]curveball.competitions.double_baranyi_roberts_ode0
-
curveball.competitions.
fitness_LTEE
(y, ref_strain=0, assay_strain=1, t0=0, t1=- 1, ci=0)[source]¶ Calculate relative fitness according to the definition used in the Long Term Evolutionary Experiment (LTEE) 3, where \(A(t), B(t)\) are population densities of assay strain A and reference strain B at time t:
\[\omega = \frac{\log{(A(t) / A(0))}}{\log{(B(t) / B(0))}}\]If ci > 0, treats the third axis of y as replicates (from a parameteric boostrap) to calculate the confidence interval on the fitness.
- ynumpy.ndarray
array of population densities, as produced by
compete()
, where the first axis is time and the second axis is strain. A third axis is also applicable if ci>0.- ref_strainint, optional
the index of the reference strain within y. This strain’s fitness is set to 1 by definition. Defaults to 0 (first).
- assay_strainint, optional
the index of the assay strain within y. The result will be the fitness of this strain relative to the fitness of the reference strain. Defaults to 1 (second).
- t0int, optional
the index of the time point from which to start the calculation of the relative fitness, defaults to 0 (first).
- t1int, optional
the index of the time point at which to end the calculation of the relative fitness, defaults to -1 (last).
- cifloat between 0 and 1, optional
if not zero, a confidence interval will be calculated using the third axis of y as replicates.
- wfloat
the fitness of the assay strain relative to the reference strain.
- lowfloat
if ci > 0 and y.ndim = 3, this is the low limit of the confidence interval of the fitness
- highfloat
if ci > 0 and y.ndim = 3, this is the higher limit of the confidence interval of the fitness
- ValueError
if confidence interval is requested (ci > 0) but there are no replicates (y.ndim != 3).
The result may depend on the choice of t0 and t1 as well as the strain designations (ref_strain and assay_strain).
- 3
Wiser, M. J. & Lenski, R. E. 2015 A Comparison of Methods to Measure Fitness in Escherichia coli. PLoS One.
-
curveball.competitions.
selection_coefs_ts
(t, y, ax=None, PLOT=False)[source]¶ Calculate selection coefficient according to the following formula[2]_, where \(A(t), B(t)\) are population densities of assay strain A and reference strain B at time t:
\[s = \frac{d}{dt} \log{\frac{A(t)}{B(t)}}\]- tnumpy.ndarray
array of time points, as produced by
compete()
- ynumpy.ndarray
array of population densities, as produced by
compete()
, where the first axis is time and the second axis is strain.- axmatplotlib.axes.Axes, optional
if PLOT is
True
, an axes to plot into; if not provided, a new one is created.- PLOTbool, optional
if
True
, the function will plot the curve of s as a function of t.
- svalsnumpy.ndarray
the selection coefficients of the assay strain relative to the reference strain over time.
- figmatplotlib.figure.Figure
figure object.
- axnumpy.ndarray of matplotlib.axes.Axes
array of axes objects.
This formula assumes that the frequencies of the strains follow a logistic curve. Lag phases, interactions, etc. may cause this formula to become irrelevant.
- 12
Chevin, L-M. 2011. On Measuring Selection in Experimental Evolution. Biology Letters.