{ "cells": [ { "cell_type": "markdown", "id": "0", "metadata": {}, "source": [ "\n", "# Upgrade guide to 0.20\n", "\n", "With version 0.20, zfit prepares for a more stable and user-friendly interface. This guide will help you to upgrade your code to the new version and demonstrate the most significant changes. It is meant for people who are already familiar with zfit and want to upgrade their code.\n", "\n", "**Not all changes are everywhere reflected in the docs, [help with the docs is welcome](https://github.com/zfit/zfit/issues/556)** as well as [adding a more polished PDF](https://github.com/zfit/zfit/issues/512). (See [all issues with contributions wanted](https://github.com/zfit/zfit/issues?q=is%3Aissue+is%3Aopen+sort%3Aupdated-desc) or [reach out](https://github.com/zfit/zfit#contact) to us on Mattermost, Gitter, GitHub, or e-mail if you're interested **in becoming a contributor**, from beginners to advanced)." ] }, { "cell_type": "code", "execution_count": null, "id": "1", "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "# standard imports\n", "import zfit\n", "import zfit.z.numpy as znp # use this \"numpy-like\" for mathematical operations\n", "import zfit_physics as zphys # physics module, with new physics-inspired PDFs\n", "from zfit import z" ] }, { "cell_type": "code", "execution_count": null, "id": "2", "metadata": {}, "outputs": [], "source": [ "@z.function\n", "def maximum(x, y):\n", " return znp.maximum(x, y)\n", "\n", "\n", "# example usage of the numpy-like backend, use it if possible\n", "print(f\"sqrt with np:{np.sqrt(maximum(1., 2.))}, with z:{znp.sqrt(maximum(1., 2.))}\")\n", "print(f\"vectorized: {np.sqrt(maximum(znp.array([1., 2.]), znp.array([3., 4.])))}\")" ] }, { "cell_type": "markdown", "id": "3", "metadata": {}, "source": [ "## Parameters\n", "\n", "The largest news comes from parameters: the `NameAlreadyTakenError` is gone (!). Multiple parameters with the same name can now be created and co-exist. The only limit is: they must not exist within the same PDF/loss etc., as this would lead to ambiguities." ] }, { "cell_type": "code", "execution_count": null, "id": "4", "metadata": {}, "outputs": [], "source": [ "param1 = zfit.Parameter(\"param1\", 1, 0, 10)\n", "param1too = zfit.Parameter(\"param1\", 2, 0, 10)\n", "# no error!" ] }, { "cell_type": "markdown", "id": "5", "metadata": {}, "source": [ "## Labels\n", "\n", "Many objects, including parameters, can now have a label. This label is purely human-readable and can be used for plotting, printing, etc. It can contain arbitrary characters.\n", "\n", "The `name` of objects still exists and will in a future version probably be used for identification purpose (in conjunction with serialization). Therefore, use `label` for human-readable names and avoid `name` for that purpose." ] }, { "cell_type": "code", "execution_count": null, "id": "6", "metadata": {}, "outputs": [], "source": [ "param1again = zfit.Parameter(\"param1\", 3, 0, 10, label=r\"$param_1$ [GeV$^2$]\")" ] }, { "cell_type": "markdown", "id": "7", "metadata": {}, "source": [ "## Space\n", "\n", "As explained [in the github discussion thread](https://github.com/zfit/zfit/discussions/533) there are major improvements and changes.\n", "- multispaces (adding two spaces, having disjoint observables) are now deprecated and will be removed. The new `TruncatedPDF` supports multiple spaces as limits and achieves a similar, if not better, functionality.\n", "- the return type of `Space.limit` will be changed in the future. For forwards compatibility, use `Space.v1.limit` (or similar methods) instead of `Space.limit`. The old one is still available via `Space.v0.limit`.\n", "- new ways of creating spaces\n", "\n", "More [examples on how to use the new spaces](https://github.com/zfit/zfit/blob/develop/examples/spaces.py)." ] }, { "cell_type": "code", "execution_count": null, "id": "8", "metadata": {}, "outputs": [], "source": [ "obs1 = zfit.Space(\"obs1\", -1, 1) # no tuple needed anymore\n", "obs2 = zfit.Space(\"obs2\", lower=-1, upper=1, label=\"observable two\")\n", "\n", "# create a space with multiple observables\n", "obs12 = obs1 * obs2" ] }, { "cell_type": "code", "execution_count": null, "id": "9", "metadata": {}, "outputs": [], "source": [ "# limits are now as one would naively expect, and area has been renamed to volume (some are tensors, but that doesn't matter: they behave like numpy arrays)\n", "print(f\"old: limits={obs12.v0.limits}\")\n", "print(f\"new: limits={obs12.v1.limits}\") # we have just 1D arrays\n" ] }, { "cell_type": "code", "execution_count": null, "id": "10", "metadata": {}, "outputs": [], "source": [ "# this allows, for example, for a more intuitive way\n", "np.linspace(*obs12.v1.limits, num=7)" ] }, { "cell_type": "markdown", "id": "11", "metadata": {}, "source": [ "## Data\n", "\n", "Data handling has been significantly simplified and streamlined.\n", "- all places (or most) take directly `numpy` arrays, tensors or `pandas DataFrame` as input, using a `Data` object is not necessary anymore (but convenient, as it cuts the data to the space)\n", "- `Data` objects can now be created without the specific constructors (e.g. `zfit.Data.from_pandas`), but directly with the data. The constructors are still available for convenience and for more options.\n", "- `Data` objects are now stateless and offer instead `with`-methods to change the data. For example, `with_obs`, `with_weights` (can be `None` to have without weights) etc.\n", "- The `SamplerData` has been overhauld and has now a more public API including a `update_data` method that allows to change the data without creating a new object and without relying on a `create_sampler` method from a PDF.\n", "- `zfit.data.concat` has been added to concatenate data objects, similar to `pd.concat`." ] }, { "cell_type": "code", "execution_count": null, "id": "12", "metadata": {}, "outputs": [], "source": [ "data1_obs1 = zfit.Data(np.random.uniform(0, 1, 1000), obs=obs1)\n", "data1_obs2 = zfit.Data(np.random.uniform(0, 1, 1000), obs=obs2, label=\"My favourite $x$\")\n", "data2_obs1 = zfit.Data(np.random.normal(0, 1, 1000), obs=obs1)\n", "\n", "# similar like pd.concat\n", "data_obs12 = zfit.data.concat([data1_obs1, data1_obs2], axis=\"columns\")\n", "data_obs1 = zfit.data.concat([data1_obs1, data2_obs1], axis=\"index\")" ] }, { "cell_type": "code", "execution_count": null, "id": "13", "metadata": {}, "outputs": [], "source": [ "# data can be accessed with \"obs\" directly\n", "data_obs12[\"obs1\"] # returns a numpy array shape (n,)\n", "data_obs12[[\"obs2\", \"obs1\"]] # returns a numpy array shape (n, 2)" ] }, { "cell_type": "markdown", "id": "14", "metadata": {}, "source": [ "## Binning\n", "\n", "Using a binned space has now the effect of creating a binned version. This happens for `Data` and `PDF` objects." ] }, { "cell_type": "code", "execution_count": null, "id": "15", "metadata": {}, "outputs": [], "source": [ "obs1_binned = obs1.with_binning(12)\n", "data_binned = zfit.Data(np.random.normal(0, 0.2, 1000), obs=obs1_binned)" ] }, { "cell_type": "markdown", "id": "16", "metadata": {}, "source": [ "## PDFs\n", "\n", "- there are a plethora of new PDFs, mostly covering physics inspired use-cases. Amongst the interesting ones are a `GeneralizedCB`, a more general version of the `DoubleCB` that should be preferred in the future. A Voigt profile is available, Bernstein polynomials, QGauss, GaussExpTail, etc. and in [zfit-physics](https://zfit.readthedocs.io/en/latest/user_api/zfit.pdf.html#physics-pdfs) HEP specific PDFS , from `CMSShape`, `Cruijff`, `Novosibirsk` and more. \n", "- the `TruncatedPDF` has been added to allow for a more flexible way of truncating a PDF. Any PDF can be converted to a truncated version using `to_truncated` (which, by default, truncates to the limits of the space).\n", "- PDFs have a new `plot` method that allows for a quick plotting of the PDF (it takes an \"obs\" argument that allows to simply project it!). This is still experimental and may changes, the main purpose is to allow for a quick check of the PDF in interactive environments. The function is fully compatible with matplotlib and takes an `ax` argument, it also allows to pass through any keyword arguments to the plotting function." ] }, { "cell_type": "code", "execution_count": null, "id": "17", "metadata": {}, "outputs": [], "source": [ "# all the new PDFs\n", "print(zfit.pdf.__all__)" ] }, { "cell_type": "code", "execution_count": null, "id": "18", "metadata": {}, "outputs": [], "source": [ "print(zphys.pdf.__all__)" ] }, { "cell_type": "code", "execution_count": null, "id": "19", "metadata": { "jupyter": { "is_executing": true } }, "outputs": [], "source": [ "# create a PDF\n", "pdf = zfit.pdf.Gauss(\n", " mu=zfit.Parameter(\"mu\", 1.2), sigma=param1again, obs=obs1, extended=1000\n", ") # using an extended PDF, the truncated pdf automatically rescales\n", "pdf.plot.plotpdf() # quick plot\n", "# truncate it\n", "pdf_truncated = pdf.to_truncated(limits=(-0.7, 0.1))\n", "pdf_truncated.plot.plotpdf()" ] }, { "cell_type": "code", "execution_count": null, "id": "20", "metadata": { "jupyter": { "is_executing": true } }, "outputs": [], "source": [ "# a truncated PDF can also have multiple limits, replacing the MultiSpace concept\n", "range1 = zfit.Space('obs1', -0.9, -0.1, label=\"lower range\") # an excellent use-case of a label\n", "range2 = zfit.Space('obs1', 0.3, 0.5, label=\"upper range\")\n", "pdf_disjoint = pdf.to_truncated(limits=(range1, range2))\n", "pdf_disjoint.plot.plotpdf()" ] }, { "cell_type": "code", "execution_count": null, "id": "21", "metadata": { "jupyter": { "is_executing": true } }, "outputs": [], "source": [ "# binned pdfs from space work like the data\n", "gauss_binned = zfit.pdf.Gauss(mu=zfit.Parameter(\"mu\", 2.5), sigma=param1again, obs=obs1_binned, extended=1000)" ] }, { "cell_type": "code", "execution_count": null, "id": "22", "metadata": { "jupyter": { "is_executing": true } }, "outputs": [], "source": [ "# as mentioned before, PDFs can be evaluated directly on numpy arrays or pandas DataFrames\n", "pdf.pdf(data_obs1.to_pandas())" ] }, { "cell_type": "markdown", "id": "23", "metadata": {}, "source": [ "## Loss and minimizer\n", "\n", "They stay mostly the same (apart from improvements behind the scenes and bugfixes).\n", "Losses take now directly the data and the model, without the need of a `Data` object (if the data is already cut to the space).\n", "\n", "To use the automatic gradient, set `gradient=\"zfit\"` in the minimizer. This can speed up the minimization for larger fits.\n", "\n", "### Updated params\n", "\n", "The minimizer currently updates the parameter default values after each minimization. To keep this behavior, add `update_params()` call after the minimization.\n", "\n", "(experimentally, the update can be disabled with `zfit.run.experimental_disable_param_update(True)`, this will probably be the default in the future. Pay attention that using this experimental features most likely breaks current scripts. Feedback on this new feature is highly welcome!)" ] }, { "cell_type": "code", "execution_count": null, "id": "24", "metadata": { "jupyter": { "is_executing": true } }, "outputs": [], "source": [ "loss = zfit.loss.ExtendedUnbinnedNLL(model=pdf, data=data2_obs1.to_pandas())\n", "minimizer = zfit.minimize.Minuit(\n", " gradient=\"zfit\"\n", ") # to use the automatic gradient -> can fail, but can speed up the minimization" ] }, { "cell_type": "code", "execution_count": null, "id": "25", "metadata": { "jupyter": { "is_executing": true } }, "outputs": [], "source": [ "result = minimizer.minimize(loss).update_params()\n", "pdf.plot.plotpdf(full=False) # plot only the pdf, no labels" ] }, { "cell_type": "markdown", "id": "26", "metadata": {}, "source": [ "## Result\n", "\n", "The result has more usage for setting parameters, the `update_params` method is now available and can be used to set the parameters to the minimum as seen above.\n", "\n", "It can also be used in a context manager to temporarily set the parameters to the minimum and restore them afterwards.\n" ] }, { "cell_type": "code", "execution_count": null, "id": "27", "metadata": {}, "outputs": [], "source": [ "param1again.set_value(1.5)\n", "with result:\n", " print(f\"param1 set temporarily to {param1again}\")\n", "print(f\"param1 is now {param1again} again\")" ] }, { "cell_type": "code", "execution_count": null, "id": "28", "metadata": {}, "outputs": [], "source": [ "# or to set the parameters to the minimum\n", "zfit.param.set_values(result) # supports also a dict of {param: value}!\n", "print(param1again)" ] }, { "cell_type": "markdown", "id": "29", "metadata": {}, "source": [ "## Serialization\n", "\n", "The result can now be pickled like any object in zfit!\n", "(this was not possible before, only after calling `freeze` on the result)\n", " \n", "This works directly using `dill` (a library that extends `pickle`), but can fail if the garbage collector is not run. Therefore, zfit provides a slightly modified `dill` that can work as a drop-in replacement.\n", "\n", "The objects can be saved and loaded again and used as before. Ideal to store the minimum of a fit and use it later for statistical treatments, plotting, etc." ] }, { "cell_type": "code", "execution_count": null, "id": "30", "metadata": {}, "outputs": [], "source": [ "result_serialized = zfit.dill.dumps(result)\n", "result_deserialized = zfit.dill.loads(result_serialized)" ] }, { "cell_type": "code", "execution_count": null, "id": "31", "metadata": {}, "outputs": [], "source": [ "# the result can be used as before\n", "result_deserialized.hesse() # the default name is now \"hesse\" and not \"minuit_hesse\"\n", "result_deserialized.errors() # the default name is now \"errors\" and not \"minuit_minos\"\n", "print(result_deserialized)" ] }, { "cell_type": "markdown", "id": "32", "metadata": {}, "source": [ "## Parameters as arguments\n", "\n", "The values of the parameters can now be directly used as arguments in functions of PDFs/loss. Methods in the pdf also take the parameters as arguments.\n", "\n", "The name of the argument has to match the name of the parameter given in initialization (or can also be the parameter itself)." ] }, { "cell_type": "code", "execution_count": null, "id": "33", "metadata": {}, "outputs": [], "source": [ "from matplotlib import pyplot as plt\n", "\n", "x = znp.linspace(*obs1.v1.limits, 1000)\n", "plt.plot(x, pdf.pdf(x, params={\"param1\": 1.5}), label=\"param1=1.5\")\n", "plt.plot(x, pdf.pdf(x, params={param1again: 2.5}), label=\"param1=2.5\")\n", "plt.legend()" ] }, { "cell_type": "code", "execution_count": null, "id": "34", "metadata": {}, "outputs": [], "source": [ "import tqdm\n", "\n", "# similar for the loss\n", "param1dict = result_deserialized.params[param1again]\n", "param1min = param1dict[\"value\"]\n", "lower, upper = param1dict[\"errors\"][\"lower\"], param1dict[\"errors\"][\"upper\"]\n", "x = np.linspace(param1min + 2 * lower, param1min + 2 * upper, 50)\n", "y = []\n", "param1again.floating = False # not minimized\n", "for x_i in tqdm.tqdm(x):\n", " param1again.set_value(x_i)\n", " minimizer.minimize(loss).update_params() # set nuisance parameters to minimum\n", " y.append(loss.value())\n", "plt.plot(x, y)" ] }, { "cell_type": "code", "execution_count": null, "id": "35", "metadata": {}, "outputs": [], "source": [ "# a result can also be used as argument for PDFs or, here, for losses\n", "loss_before = loss.value()\n", "loss_min = loss.value(params=result_deserialized) # evaluate at minimum\n", "print(f\"loss before (random from before): {loss_before:.7} minimum value: {loss_min:.7} vs {result_deserialized.fmin:.7}\")" ] }, { "cell_type": "code", "execution_count": null, "id": "36", "metadata": {}, "outputs": [], "source": [ "# creating a PDF looks also different, but here we use the name of the parametrization and the axis (integers)\n", "\n", "\n", "class MyGauss2D(zfit.pdf.ZPDF):\n", " _PARAMS = (\"mu\", \"sigma\")\n", " _N_OBS = 2\n", "\n", " @zfit.supports() # this allows the params magic\n", " def _unnormalized_pdf(self, x, params):\n", " x0 = x[0] # this means \"axis 0\"\n", " x1 = x[1] # this means \"axis 1\"\n", " mu = params[\"mu\"]\n", " sigma = params[\"sigma\"]\n", " return znp.exp(-0.5 * ((x0 - mu) / sigma) ** 2) * x1 # fake, just for demonstration" ] }, { "cell_type": "code", "execution_count": null, "id": "37", "metadata": {}, "outputs": [], "source": [ "gauss2D = MyGauss2D(mu=0.8, sigma=param1again, obs=obs12, label=\"2D Gaussian$^2$\")\n", "gauss2D.plot.plotpdf(obs=\"obs1\") # we can project the 2D pdf to 1D" ] }, { "cell_type": "code", "execution_count": null, "id": "38", "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 2 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython2", "version": "2.7.6" } }, "nbformat": 4, "nbformat_minor": 5 }