Importing PEtab Standard Format

PEtab, from which PEtab.jl gets its name, is a flexible, table-based standard format for specifying parameter estimation problems [7]. If a problem is provided in this standard format, PEtab.jl can import it directly. This tutorial covers how to import PEtab problems.

Input - a Valid PEtab Problem

To import a PEtab problem, a valid PEtab problem is required. A tutorial on creating a PEtab problem can be found in the PEtab documentation, and a linting tool is available in Python for checking correctness. Additionally, PEtab.jl performs several format checks when importing the problem.

A collection of valid PEtab problems is also available in the PEtab benchmark repository. In this tutorial, we will use an already published model from the PEtab benchmark repository. Specifically, we will consider the STAT5 signaling model, referred to here as the Boehm model (after the first author) [8]. The PEtab files for this model can be found here.

Importing a PEtabModel

A PEtab problem consists of five files: an SBML model file, a table with simulation conditions, a table with observables, a table with measurements, and a table with parameters to estimate. These are tied together by a YAML file, and to import a problem, you only need to provide the YAML file path:

using PEtab
# path_yaml depends on where the model is saved
path_yaml = joinpath(@__DIR__, "assets", "boehm", "Boehm_JProteomeRes2014.yaml")
model = PEtabModel(path_yaml)

Given a PEtabModel, it is straightforward to create a PEtabODEProblem:

petab_prob = PEtabODEProblem(model)
PEtabODEProblem: boehm with ODE-states 8 and 9 parameters to estimate
---------------- Problem options ---------------
Gradient method: ForwardDiff
Hessian method: ForwardDiff
ODE-solver nllh: Rodas5P
ODE-solver gradient: Rodas5P

As described in the starting tutorial, this PEtabODEProblem can then be used for parameter estimation, or Bayesian inference. For tunable options when importing a PEtab problem, see the API documentation.

What Happens During PEtab Import (Deep Dive)

When importing a PEtab model, several things happen:

  1. The SBML file is converted into a Catalyst.jl ReactionSystem using the SBMLImporter.jl package. This ReactionSystem is then converted into an ODESystem. During this step, the model is symbolically pre-processed, which includes computing the ODE Jacobian symbolically. The latter typically improves simulation performance.
  2. The observable PEtab table is translated into Julia functions that compute observables (h), measurement noise (σ), and initial values (u0).
  3. Any potential model events are translated into Julia callbacks.

All of these steps happen automatically. By setting write_to_file=true when importing the model, the generated model functions can be found in the dir_yaml/Julia_model_files/ directory.

References

[7]
L. Schmiester, Y. Schälte, F. T. Bergmann, T. Camba, E. Dudkin, J. Egert, F. Fröhlich, L. Fuhrmann, A. L. Hauber, S. Kemmer and others. PEtab—Interoperable specification of parameter estimation problems in systems biology. PLoS computational biology 17, e1008646 (2021).
[8]
M. E. Boehm, L. Adlung, M. Schilling, S. Roth, U. Klingmüller and W. D. Lehmann. Identification of isoform-specific dynamics in phosphorylation-dependent STAT5 dimerization by quantitative mass spectrometry and mathematical modeling. Journal of proteome research 13, 5685–5694 (2014).