{
"cells": [
{
"cell_type": "markdown",
"id": "fd77aec4",
"metadata": {},
"source": [
"# Model\n",
"\n",
"As shown previously, when `Simulation.from_config()` is called, a `SimulationState` object is created. This is done via the class method `SimulationState.from_config()`. This model object contains important information about the shell structure, density, abundance, radiative temperature, and dilution factor throughout the supernova.\n",
"\n",
"Throughout this notebook, we show various configuration inputs into the TARDIS model and the resulting model. In interactive mode, these parameters can be varied to explore how the model changes. Editing configuration parameters in the notebook is explained [here](../../io/configuration/tutorial_read_configuration.ipynb).\n",
"\n",
"This notebook is based on the [built-in TARDIS models](../../io/configuration/components/models/index.rst#built-in-structure-density-and-abundance), but these parameters can also be input via [custom model configurations](../../io/configuration/components/models/index.rst#custom-model-configurations), the [CSVY model](../../io/configuration/components/models/index.rst#csvy-model) or the [custom abundance widget](../../io/visualization/how_to_abundance_widget.ipynb).\n",
"\n",
"## Shell Structure\n",
"\n",
"TARDIS models supernovae as expanding homologously. This means that at the beginning of the explosion, the\n",
"supernova starts at a single point and proceeds to expand radially outward such that the ratio of the velocity of\n",
"the ejecta to the distance from the ejecta to the supernova's center is uniform throughout the supernova. As an\n",
"example, if the outer edge of the ejecta moves outward at some velocity $v_\\mathrm{boundary\\_outer}$, the\n",
"velocity of the ejecta half way between the outer edge and the center would be\n",
"$\\frac{v_\\mathrm{boundary\\_outer}}{2}$. The animation below demonstrates this type of expansion.\n",
"\n",
"TARDIS simulates radiative transfer between an inner boundary (the photosphere) and an outer\n",
"boundary (the outer edge of the supernova ejecta). The velocity of the inner boundary\n",
"$v_\\mathrm{boundary\\_inner}$ and the velocity of the outer boundary $v_\\mathrm{boundary\\_outer}$ are\n",
"supplied in the configuration file, as well as the time after the explosion for\n",
"which TARDIS is calculating the spectrum ($t_\\mathrm{explosion}$). The radii of the inner and outer boundaries\n",
"are therefore calculated by $r_\\mathrm{boundary\\_inner}=v_\\mathrm{boundary\\_inner}*t_\\mathrm{explosion}$ and\n",
"$r_\\mathrm{boundary\\_outer}=v_\\mathrm{boundary\\_outer}*t_\\mathrm{explosion}$. Plasma at a distance $r$\n",
"from the center of the supernova would then be traveling outward at a speed $v=\\frac{r}{r_\\mathrm{boundary\\_outer}}v_\\mathrm{boundary\\_outer} = \\frac{r}{t_\\mathrm{explosion}}$. This is\n",
"also shown in the animation.\n",
"\n",
"Additionally, TARDIS divides the space between the inner and outer computational boundaries into cells -- radial\n",
"shells for which the state/condition of the ejecta is (spatially) constant. In the animation, 6 cells are shown,\n",
"being divided by the light blue lines. As TARDIS is a time-independent code which calculates the spectra at an\n",
"instant in time, the radii of the boundaries (either of the computational domain or of the cells) do not change\n",
"throughout the simulation.\n",
"\n",
"\n",
"\n",
"\n",
"We now demonstrate how the shell structure works using the TARDIS code."
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "f48c4843",
"metadata": {},
"outputs": [],
"source": [
"# We import the necessary packages\n",
"from tardis.io.configuration.config_reader import Configuration\n",
"from tardis.model import SimulationState\n",
"from tardis.io.atom_data import download_atom_data\n",
"from tardis.io.atom_data.base import AtomData\n",
"from astropy import units as u\n",
"import matplotlib.pyplot as plt\n",
"import copy\n",
"\n",
"# We download the atomic data needed to run the simulation\n",
"download_atom_data('kurucz_cd23_chianti_H_He')\n",
"atom_data = AtomData.from_hdf('kurucz_cd23_chianti_H_He.h5')"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "8b3a14a9",
"metadata": {},
"outputs": [],
"source": [
"# Set up a base configuration to use throughout the notebook\n",
"base_config = Configuration.from_yaml(\"tardis_example.yml\")"
]
},
{
"cell_type": "markdown",
"id": "cee054e9",
"metadata": {},
"source": [
"In the cell below, we set up a model. We use the [specific structure](../../io/configuration/components/models/index.rst#specific-structure) where we supply $t_\\mathrm{explosion}$, the velocity of the inner and outer boundaries of the supernova (labeled `start` and `stop`), and the number of shells (labeled `num`). The shells are then evenly spaced between the inner and outer boundaries of the supernova. The time after the explosion, the inner and outer velocities, and the number of shells can be varied to get different shell structures. The `SimulationState` object stores information about the model in the following attributes: `velocity` shows the velocity of the shell boundaries, `v_inner` shows the velocities of the inner boundaries of each shell, `v_outer` shows the velocity of the outer boundaries of each shell, and `v_middle` shows the velocity of the middle of each shell. Similarly, `radius`, `r_inner`, `r_outer`, and `r_middle` show the radii of each shell boundary, the inner boundaries, the outer boundaries, and the middles of each shell, respectively. `v_inner_boundary` shows the velocity of the inner boundary of the computational domain, and `v_outer_boundary` shows the velocity of the outer boundary of the computational domain. Finally, `volume` shows the volume of each shell, calculated via the formula of the volume of a spherical shell: $V=\\frac{4}{3}\\pi (r_\\mathrm{outer}^3-r_\\mathrm{inner}^3)$."
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "17752fd2",
"metadata": {
"scrolled": false
},
"outputs": [],
"source": [
"shell_config = copy.deepcopy(base_config)\n",
"\n",
"shell_config.supernova.time_explosion = 10 * u.day\n",
"\n",
"# This line is necessary to indicate we are using a built-in shell structure.\n",
"# Do not change it.\n",
"shell_config.model.structure.type = 'specific'\n",
"\n",
"shell_config.model.structure.velocity.start = 1000 * u.km/u.s\n",
"shell_config.model.structure.velocity.stop = 2000 * u.km/u.s\n",
"shell_config.model.structure.velocity.num = 20\n",
"\n",
"shell_simulation_state = SimulationState.from_config(shell_config, atom_data=atom_data)\n",
"\n",
"print('velocity:\\n', shell_simulation_state.velocity)\n",
"print('v_inner:\\n', shell_simulation_state.v_inner)\n",
"print('v_outer:\\n', shell_simulation_state.v_outer)\n",
"print('v_middle:\\n', shell_simulation_state.v_middle)\n",
"print('v_inner_boundary:\\n', shell_simulation_state.v_inner_boundary)\n",
"print('v_outer_boundary:\\n', shell_simulation_state.v_outer_boundary)\n",
"print('radius:\\n', shell_simulation_state.radius)\n",
"print('r_inner:\\n', shell_simulation_state.r_inner)\n",
"print('r_outer:\\n', shell_simulation_state.r_outer)\n",
"print('r_middle:\\n', shell_simulation_state.r_middle)\n",
"print('volume:\\n', shell_simulation_state.volume)"
]
},
{
"cell_type": "markdown",
"id": "1ee56110",
"metadata": {},
"source": [
"Notice that `radius = velocity*time_explosion`, and similarly for `r_inner`, `r_outer`, and `r_middle`. You can get the radius of the photosphere via `v_inner_boundary*time_explosion` and outer edge of the supernova via `v_boundary_outer*time_explosion`.\n",
"\n",
"