{ "cells": [ { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "# pyPESTO: Getting started\n", "\n", "This notebook takes you through the first steps to get started with [pyPESTO](https://github.com/ICB-DCM/pyPESTO).\n", "\n", "\"pyPESTO\n", "\n", "pyPESTO is a python package for parameter inference, offering a unified interface to various optimization and sampling methods.\n", "pyPESTO is highly modular and customizable, e.g., with respect to objective function definition and employed inference algorithms." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "import logging\n", "\n", "import amici.sim.sundials as asd\n", "import matplotlib as mpl\n", "import numpy as np\n", "\n", "import pypesto.optimize as optimize\n", "import pypesto.petab\n", "import pypesto.visualize as visualize\n", "\n", "np.random.seed(1)\n", "\n", "# increase image resolution\n", "mpl.rcParams[\"figure.dpi\"] = 300" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "## 1. Objective Definition\n", "\n", "pyPESTO allows the definition of custom objectives and offers support for objectives defined in the [PEtab](https://github.com/PEtab-dev/PEtab) format.\n", "\n", "### Custom Objective Definition\n", "\n", "You can define an objective via a python function. Also providing an analytical gradient (and potentially also a Hessian) improves the performance of Gradient/Hessian-based optimizers. When accessing parameter uncertainties via profile-likelihoods/sampling, pyPESTO interprets the objective function as the negative-log-likelihood/negative-log-posterior. A more in-depth construction of a custom objective function can be found in [a designated example notebook.](./custom_objective_function.ipynb)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# define objective function\n", "def f(x: np.ndarray) -> np.ndarray:\n", " return x[0] ** 2 + x[1] ** 2\n", "\n", "\n", "# define gradient\n", "def grad(x: np.ndarray) -> np.ndarray:\n", " return 2 * x\n", "\n", "\n", "# define objective\n", "custom_objective = pypesto.Objective(fun=f, grad=grad)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "Define lower and upper parameter bounds and create an optimization problem." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# define optimization bounds\n", "lb = np.array([-10, -10])\n", "ub = np.array([10, 10])\n", "\n", "# create problem\n", "custom_problem = pypesto.Problem(objective=custom_objective, lb=lb, ub=ub)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "Now choose an optimizer to perform the optimization. `minimize` uses multi-start optimization, meaning that the optimization is run `n_start` times from different initial values, in case the problem contains multiple local optima (which of course is not the case for this toy problem)." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# choose optimizer\n", "optimizer = optimize.ScipyOptimizer()\n", "\n", "# do the optimization\n", "result_custom_problem = optimize.minimize(\n", " problem=custom_problem, optimizer=optimizer, n_starts=10\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "`result_custom_problem.optimize_result` now stores a list, that contains the results and metadata of the individual optimizer runs (ordered by function value)." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# E.g., The best model fit was obtained by the following optimization run:\n", "result_custom_problem.optimize_result.list[0]" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# Objective function values of the different optimizer runs:\n", "result_custom_problem.optimize_result.get_for_key(\"fval\")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Problem Definition via PEtab\n", "\n", "#### Background on PEtab\n", "\n", "\"PEtab\n", "\n", "PyPESTO supports the [PEtab](https://github.com/PEtab-dev/PEtab) standard. PEtab is a data format for specifying parameter estimation problems in systems biology.\n", "\n", "A PEtab problem consist of an [SBML](https://sbml.org) file, defining the model topology and a set of `.tsv` files, defining experimental conditions, observables, measurements and parameters (and their optimization bounds, scale, priors...). All files that make up a PEtab problem can be structured in a `.yaml` file. The `pypesto.Objective` coming from a PEtab problem corresponds to the negative-log-likelihood/negative-log-posterior distribution of the parameters.\n", "\n", "For more details on PEtab, the interested reader is referred to [PEtab's format definition](https://petab.readthedocs.io/en/latest/documentation_data_format.html), for examples, the reader is referred to the [PEtab benchmark collection](https://github.com/Benchmarking-Initiative/Benchmark-Models-PEtab). The Model from _[Böhm et al. JProteomRes 2014](https://pubs.acs.org/doi/abs/10.1021/pr5006923)_ is part of the benchmark collection and will be used as the running example throughout this notebook.\n", "\n", "PyPESTO provides an interface to the model simulation tool [AMICI](https://github.com/AMICI-dev/AMICI) for the simulation of Ordinary Differential Equation (ODE) models specified in the SBML format." ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "#### Basic Model Import and Optimization\n", "\n", "The first step is to import a PEtab problem and create a `pypesto.problem` object:" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "%%capture\n", "# directory of the PEtab problem\n", "petab_yaml = \"./conversion_reaction/conversion_reaction.yaml\"\n", "\n", "importer = pypesto.petab.PetabImporter.from_yaml(petab_yaml)\n", "problem = importer.create_problem(verbose=False)" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "petab_yaml" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "Next, we choose an `optimizer` to perform the multi start optimization." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "%%time\n", "%%capture\n", "\n", "# choose optimizer\n", "optimizer = optimize.ScipyOptimizer()\n", "\n", "# do the optimization\n", "result = optimize.minimize(problem=problem, optimizer=optimizer, n_starts=5)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "`result.optimize_result` contains a list with the ordered optimization results." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# E.g., best model fit was obtained by the following optimization run:\n", "result.optimize_result.list[0]" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# Objective function values of the different optimizer runs:\n", "result.optimize_result.get_for_key(\"fval\")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "## 2. Optimizer Choice\n", "\n", "PyPESTO provides a unified interface to a variety of optimizers of different types.\n", "Examples are:\n", "\n", "* All [scipy optimizer](https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html#scipy.optimize.minimize) (`optimize.ScipyOptimizer(method=)`)\n", " * _function-value_ or _least-squares_-based optimizers\n", " * _gradient_ or _hessian_-based optimizers\n", "* [IpOpt](https://pypi.org/project/ipopt/) (`optimize.IpoptOptimizer()`)\n", " * Interior point method\n", "* [Dlib](https://dlib.net) (`optimize.DlibOptimizer(options={'maxiter': })`)\n", " * Global optimizer\n", " * Gradient-free\n", "* [FIDES](https://github.com/fides-dev/fides/) (`optimize.FidesOptimizer()`)\n", " * Interior Trust Region optimizer\n", "* [Particle Swarm](https://github.com/ljvmiranda921/pyswarms) (`optimize.PyswarmsOptimizer()`)\n", " * Particle swarm algorithm\n", " * Gradient-free\n", "* [CMA-ES](https://pypi.org/project/cma-es/) (`optimize.CmaOptimizer()`)\n", " * Covariance Matrix Adaptation Evolution Strategy\n", " * Evolutionary Algorithm" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "optimizer_scipy_lbfgsb = optimize.ScipyOptimizer(method=\"L-BFGS-B\")\n", "optimizer_scipy_powell = optimize.ScipyOptimizer(method=\"Powell\")\n", "\n", "optimizer_fides = optimize.FidesOptimizer(verbose=logging.ERROR)\n", "optimizer_pyswarms = optimize.PyswarmsOptimizer(par_popsize=10)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "The following performs 10 multi-start runs with different optimizers in order to compare their performance. For a faster execution of this notebook, we run only 10 starts. In application, one would use many more optimization starts: around 100-1000 in most cases.\n", "\n", "_Note_: Some of those optimizers need to be installed separately for this section to run (e.g., via `pip install pypesto[fides,pyswarms]`). Furthermore, the computation time is in the order of minutes, so you might want to skip the execution and jump to the section on large scale models." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "%%time\n", "%%capture --no-display\n", "n_starts = 10\n", "\n", "# Due to run time we already use parallelization.\n", "# This will be introduced in more detail later.\n", "engine = pypesto.engine.MultiProcessEngine()\n", "\n", "# Scipy: L-BFGS-B\n", "result_lbfgsb = optimize.minimize(\n", " problem=problem,\n", " optimizer=optimizer_scipy_lbfgsb,\n", " engine=engine,\n", " n_starts=n_starts,\n", ")\n", "\n", "# Scipy: Powell\n", "result_powell = optimize.minimize(\n", " problem=problem,\n", " optimizer=optimizer_scipy_powell,\n", " engine=engine,\n", " n_starts=n_starts,\n", ")\n", "\n", "# Fides\n", "result_fides = optimize.minimize(\n", " problem=problem,\n", " optimizer=optimizer_fides,\n", " engine=engine,\n", " n_starts=n_starts,\n", ")\n", "\n", "# PySwarms\n", "result_pyswarms = optimize.minimize(\n", " problem=problem,\n", " optimizer=optimizer_pyswarms,\n", " engine=engine,\n", " n_starts=1, # Global optimizers are usually run once. The number of particles (par_popsize) is usually the parameter that is adapted.\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Optimizer Convergence\n", "\n", "\n", "A common visualization of optimizer convergence are waterfall plots. Waterfall plots show the (ordered) results of the individual optimization runs. In general, we hope to obtain clearly visible plateaus, as they indicate optimizer convergence to local minima." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "optimizer_results = [\n", " result_lbfgsb,\n", " result_powell,\n", " result_fides,\n", " result_pyswarms,\n", "]\n", "optimizer_names = [\"Scipy: L-BFGS-B\", \"Scipy: Powell\", \"Fides\", \"PySwarms\"]\n", "\n", "pypesto.visualize.waterfall(optimizer_results, legends=optimizer_names);" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Optimizer run time\n", "\n", "Optimizer run time vastly differs among the different optimizers, as can be seen below:" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "print(\"Average Run time per start:\")\n", "print(\"-------------------\")\n", "\n", "for optimizer_name, optimizer_result in zip(\n", " optimizer_names, optimizer_results, strict=True\n", "):\n", " t = np.sum(optimizer_result.optimize_result.get_for_key(\"time\")) / n_starts\n", " print(f\"{optimizer_name}: {t:f} s\")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "## 3. Fitting of large scale models\n", "\n", "When fitting large scale models (i.e. with >100 parameters and accordingly also more data), two important issues are efficient gradient computation and parallelization.\n", "\n", "### Efficient gradient computation\n", "\n", "As seen in the example above and as can be confirmed from own experience: If fast and reliable gradients can be provided, gradient-based optimizers are favourable with respect to optimizer convergence and run time.\n", "\n", "It has been shown that adjoint sensitivity analysis is a fast and reliable method to compute gradients for large scale models, since their run time is (asymptotically) independent of the number of parameters ([Fröhlich et al. PlosCB 2017](https://journals.plos.org/ploscompbiol/article/file?id=10.1371/journal.pcbi.1005331&type=printable)).\n", "\n", "\"pyPESTO\n", "\n", "(Figure from Fröhlich et al. PlosCB 2017) Adjoint sensitivities are implemented in AMICI." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# Set gradient computation method to adjoint\n", "problem.objective.amici_solver.set_sensitivity_method(\n", " asd.SensitivityMethod.adjoint\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Parallelization\n", "\n", "Multi-start optimization can easily be parallelized by using `engines`." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "%%time\n", "%%capture\n", "\n", "# Parallelize\n", "engine = pypesto.engine.MultiProcessEngine()\n", "\n", "# Optimize\n", "result = optimize.minimize(\n", " problem=problem,\n", " optimizer=optimizer_scipy_lbfgsb,\n", " engine=engine,\n", " n_starts=25,\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "## 4. Uncertainty quantification\n", "\n", "PyPESTO focuses on two approaches to assess parameter uncertainties:\n", "\n", "* Profile likelihoods\n", "* Sampling\n", "\n", "### Profile Likelihoods\n", "\n", "[Profile likelihoods](https://academic.oup.com/bioinformatics/article/25/15/1923/213246) compute confidence intervals via a [likelihood ratio test](https://en.wikipedia.org/wiki/Likelihood-ratio_test). Profile likelihoods perform a maximum-projection of the likelihood function on the parameter of interest. The likelihood ratio test then gives a cut-off criterion via the $\\chi^2_1$ distribution.\n", "\n", "In pyPESTO, the maximum projection is solved as a maximization problem and can be obtained via" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "%%time\n", "%%capture\n", "\n", "import scipy as sp\n", "\n", "import pypesto.profile as profile\n", "\n", "options = profile.ProfileOptions(\n", " ratio_min=np.exp(-sp.stats.chi2.ppf(0.99, 1) / 2)\n", ")\n", "result = profile.parameter_profile(\n", " problem=problem,\n", " result=result,\n", " optimizer=optimizer_scipy_lbfgsb,\n", " profile_index=[0, 1],\n", " profile_options=options,\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "The maximum projections can now be inspected via:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "problem.x_scales" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# adapt x_labels.\n", "x_labels = [f\"Log({name})\" for name in problem.x_names]\n", "\n", "ax = visualize.profiles(result, x_labels=x_labels, show_bounds=True)\n", "# visualize optimal parameter values\n", "ax[0].plot(\n", " [result_fides.optimize_result.x[0][0]] * 2,\n", " [0, 1.1],\n", " linestyle=\"--\",\n", " color=\"g\",\n", ")\n", "ax[1].plot(\n", " [result_fides.optimize_result.x[0][1]] * 2,\n", " [0, 1.1],\n", " linestyle=\"--\",\n", " color=\"g\",\n", ")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The 2D profile visualization shows joint parameter paths across all profile pairs. Diagonal plots show 1D profiles (likelihood ratio vs. parameter value). Off-diagonal plots show how each non-profiled parameter moves while another is profiled, with color indicating the likelihood ratio.\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "fig, axes = pypesto.visualize.visualize_2d_profile(result)" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "Furthermore, pyPESTO allows to visualize approximate confidence intervals directly via `profile_cis`. Confidence intervals are computed by finding parameter values for which the log posterior ratio is above the approximate threshold assuming a $\\chi^2$ distribution of the likelihood test statistic. The plot shows that both model parameters are identifiable, since the confidence regions are finite and do not span the whole estimation space (from lower to upper estimation boundary)." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "ax = pypesto.visualize.profile_cis(\n", " result, confidence_level=0.95, show_bounds=True\n", ")\n", "\n", "ax.set_xlabel(\"Log(Parameter value)\");" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Parameter profiles need to be computed until the likelihood ratio is small enough to compute confidence intervals with the selected confidence level (see `ratio_min` attribute of `ProfileOptions`). It is also possible to visualize nested confidence intervals corresponding to different confidence levels." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ax = pypesto.visualize.profile_nested_cis(\n", " result, confidence_levels=[0.99, 0.95, 0.9]\n", ")\n", "\n", "ax.set_xlabel(\"Log(Parameter value)\");" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Sampling\n", "\n", "In pyPESTO, sampling from the posterior distribution can be performed as" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "import pypesto.sample as sample\n", "\n", "n_samples = 1000\n", "\n", "sampler = sample.AdaptiveMetropolisSampler()\n", "\n", "result = sample.sample(\n", " problem, n_samples=n_samples, sampler=sampler, result=result\n", ")" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "Sampling results are stored in `result.sample_result` and can be accessed e.g., via" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "result.sample_result[\"trace_x\"]" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Sampling Diagnostics\n", "\n", "Geweke's test assesses convergence of a sampling run and computes the burn-in of a sampling result. The effective sample size indicates the strength of the correlation between different samples." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "sample.geweke_test(result=result)\n", "result.sample_result[\"burn_in\"]" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "sample.effective_sample_size(result=result)\n", "result.sample_result[\"effective_sample_size\"]" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Visualization of Sampling Results" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# scatter plots\n", "visualize.sampling_scatter(result)\n", "\n", "# marginals\n", "visualize.sampling_1d_marginals(result);" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "#### Sampler Choice:\n", "\n", "Similarly to parameter optimization, pyPESTO provides a unified interface to several sampler/sampling toolboxes, as well as own implementations of samplers:\n", "\n", "* Adaptive Metropolis: `pypesto.sample.AdaptiveMetropolisSampler()`\n", "* Parallel tempering: `pypesto.sample.ParallelTemperingSampler()`\n", "* Adaptive parallel tempering: `pypesto.sample.AdaptiveParallelTemperingSampler()`\n", "* Interface to `pymc` [🔗](https://pypi.org/project/pymc/) via `pypesto.sample.PymcSampler()`" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "## 5. Storage\n", "\n", "You can store and load the results of an analysis via the `pypesto.store` module to a `.hdf5` file.\n", "\n", "### Store result" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "import tempfile\n", "\n", "import pypesto.store as store\n", "\n", "# create a temporary file, for demonstration purpose\n", "f_tmp = tempfile.NamedTemporaryFile(suffix=\".hdf5\", delete=False)\n", "result_file_name = f_tmp.name\n", "\n", "# store the result\n", "store.write_result(result, result_file_name)\n", "f_tmp.close()" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "### Load result file\n", "\n", "You can re-load a result, e.g. for visualizations:" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "pycharm": { "name": "#%%\n" } }, "outputs": [], "source": [ "# read result\n", "result_loaded = store.read_result(result_file_name)\n", "\n", "# e.g. do some visualisation\n", "visualize.waterfall(result_loaded);" ] }, { "cell_type": "markdown", "metadata": { "pycharm": { "name": "#%% md\n" } }, "source": [ "## Software Development Standards:\n", "\n", "PyPESTO is developed with the following standards:\n", "\n", "* **Open source**, code on [GitHub](https://github.com/ICB-DCM/pyPESTO).\n", "* [**Pip installable**](https://pypi.org/project/pypesto/) via: `pip install pypesto`.\n", "* **Documentation** as [RTD](https://pypesto.readthedocs.io/en/stable/) and [example **jupyter notebooks**](https://github.com/ICB-DCM/pyPESTO/tree/master/doc/example) are available.\n", "* Has **continuous integration** & extensive automated **testing**.\n", "* **Code reviews** before merging into the develop/main branch.\n", "\n", "* Currently, **5–10 people are using, extending** and (most importantly) **maintaining** pyPESTO in their \"daily business\"." ] }, { "cell_type": "markdown", "metadata": { "collapsed": false, "jupyter": { "outputs_hidden": false }, "pycharm": { "name": "#%% md\n" } }, "source": [ "## Further topics\n", "\n", "Further features are available, among them:\n", "\n", "* Model Selection\n", "* Hierarchical Optimization of scaling/noise parameters\n", "* Categorical Data" ] } ], "metadata": { "kernelspec": { "display_name": "dev_venv_11", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.11.12" } }, "nbformat": 4, "nbformat_minor": 4 }