Getting Started in AMICI

This notebook is a brief tutorial for new users that explains the first steps necessary for model simulation in AMICI, including pointers to documentation and more advanced notebooks.

Model Compilation

Before simulations can be run, the model must be imported and compiled. In this process, AMICI performs all symbolic manipulations that later enable scalable simulations and efficient sensitivity computation. The first towards model compilation is the creation of an SbmlImporter instance, which requires an SBML Document that specifies the model using the Systems Biology Markup Language (SBML).

For the purpose of this tutorial, we will use model_steadystate_scaled.xml, which is contained in the same directory as this notebook.

[1]:
import amici
sbml_importer = amici.SbmlImporter('model_steadystate_scaled.xml')

Next, we will compile the model as python extension using the amici.SBMLImporter.sbml2amici method. The first two arguments of this method are the name of the model, which will also be the name of the generated python module, and the model directory, which defines the directory in which the model module will be placed. Compilation will take a couple of seconds.

[2]:
model_name = 'model_steadystate'
model_dir = 'model_dir'
sbml_importer.sbml2amici(model_name, model_dir)

Loading the model module

To run simulations, we need to instantiate amici.Model and amici.Solver instances. As simulations requires instances matching the imported model, they have to be imported from the generated model module.

[3]:
# load the model module
model_module = amici.import_model_module(model_name, model_dir)
# instantiate model
model = model_module.getModel()
# instantiate solver
solver = model.getSolver()

The model allows the user to manipulate model related properties of simulations. This includes the values of model parameters that can be set by using amici.Model.setParameterByName. Here, we set the model parameter p1 to a value of 1e-3.

[4]:
model.setParameterByName('p1',1e-3)

In contrast, the solver instance allows the specification of simulation related properties. This includes setting options for the SUNDIALS solver such as absolute tolerances via amici.Solver.setAbsoluteTolerance. Here we set the absolute integration tolerances to 1e-10.

[5]:
solver.setAbsoluteTolerance(1e-10)

Running Model Simulations

Model simulations can be executed using the amici.runAmiciSimulations routine. By default the model does not not contain any timepoints for which the model is to be simulated. Here we define a simulation timecourse with two timepoints at 0 and 1 and then run the simulation.

[6]:
# set timepoints
model.setTimepoints([0,1])
rdata = amici.runAmiciSimulation(model, solver)

Simulation results are returned as ReturnData instance. The simulated SBML species are stored as x attribute, where rows correspond to the different timepoints and columns correspond to different species.

[7]:
rdata.x
[7]:
array([[0.1       , 0.4       , 0.7       ],
       [0.98208413, 0.51167992, 0.10633388]])

All results attributes are always ordered according to the model. For species, this means that the columns of rdata.x match the ordering of species in the model, which can be accessed as amici.Model.getStateNames

[8]:
model.getStateNames()
[8]:
('x1', 'x2', 'x3')

This notebook only explains the basics of AMICI simulations. In general, AMICI simulations are highly customizable and can also be used to simulate sensitivities. The ExampleSteadystate notebook in this folder gives more detail about the model employed here and goes into the basics of sensitivity analysis. The ExampleEquilibrationLogic notebook, builds on this by using a modified version of this model to give detailed insights into the methods and options to compute steady states before and after simulations, as well as respective sensitivities. The ExampleExperimentalConditions example notebook, goes into the details of how even more complex experimental setups, such as addition of drugs at predefined timepoints, can be simulated in AMICI. Finally, the petab notebook explains how standardized definitions of experimental data and conditions in the PEtab format can be imported in AMICI.