Files
OpenMC_NEA_Training_nov25/course/live_notebooks/CriticalitySearch.ipynb
Stafie Alex PSI 850de66b07 first files
2025-12-02 11:57:33 +01:00

431 lines
42 KiB
Plaintext

{
"cells": [
{
"cell_type": "markdown",
"id": "237392ad-3011-416d-9d6f-5f1da349c7ac",
"metadata": {},
"source": [
"# Criticality Searches\n",
"\n",
"Criticality searches, or the adjustment of some free parameter of your problem (e.g. boron concentation, control rod insertion point), is a routine part of neutronics. In this example, we'll perform a critical search based on fuel enrichment. OpenMC performs this criticality search using the methodology from Price and Roskoff (https://linkinghub.elsevier.com/retrieve/pii/S014919702300166X), a GRsecant (generalized regressive secant) method. This method takes into account the statistical uncertainty in $k$ and also dynamically adjusts the number of batches run (using more batches as one gets closer to criticality)."
]
},
{
"cell_type": "code",
"execution_count": 1,
"id": "dccc5eea-6753-4de2-933d-df86cb3ede37",
"metadata": {},
"outputs": [],
"source": [
"import openmc\n",
"\n",
"model = openmc.Model()"
]
},
{
"cell_type": "code",
"execution_count": 2,
"id": "3a733c71-dff3-4d7e-b1fd-6ea9ba2ed904",
"metadata": {},
"outputs": [],
"source": [
"uo2 = openmc.Material()\n",
"uo2.add_element('U', 1.0, enrichment=3.5)\n",
"uo2.add_element('O', 2.0)\n",
"uo2.set_density('g/cm3', 10.97)\n",
"\n",
"zirconium = openmc.Material()\n",
"zirconium.add_element('Zr', 1.0)\n",
"zirconium.set_density('g/cm3', 6.55)\n",
"\n",
"water = openmc.model.borated_water(1500, temperature=600, temp_unit='K', pressure=15, press_unit='MPa')"
]
},
{
"cell_type": "code",
"execution_count": 6,
"id": "ad31d83c-a146-4cae-a601-361327456bf8",
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"Material\n",
"\tID =\t1\n",
"\tName =\t\n",
"\tTemperature =\tNone\n",
"\tDensity =\t10.97 [g/cm3]\n",
"\tVolume =\tNone [cm^3]\n",
"\tDepletable =\tTrue\n",
"\tS(a,b) Tables \n",
"\tNuclides \n",
"\tU234 =\t0.0003166930253944235 [ao]\n",
"\tU235 =\t0.03543164439454172 [ao]\n",
"\tU238 =\t0.964089368630351 [ao]\n",
"\tU236 =\t0.00016229394971280895 [ao]\n",
"\tO16 =\t1.999242 [ao]\n",
"\tO17 =\t0.000758 [ao]"
]
},
"execution_count": 6,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"uo2"
]
},
{
"cell_type": "code",
"execution_count": 3,
"id": "eca994f7-fccc-4878-9e95-45e871091fae",
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"Material\n",
"\tID =\t3\n",
"\tName =\t\n",
"\tTemperature =\t600\n",
"\tDensity =\t0.6603788831247912 [g/cc]\n",
"\tVolume =\tNone [cm^3]\n",
"\tDepletable =\tFalse\n",
"\tS(a,b) Tables \n",
"\tS(a,b) =\t('c_H_in_H2O', 1.0)\n",
"\tNuclides \n",
"\tH1 =\t0.11085782613285239 [ao]\n",
"\tH2 =\t1.726768711152118e-05 [ao]\n",
"\tO16 =\t0.05541653607970307 [ao]\n",
"\tO17 =\t2.101083027888316e-05 [ao]\n",
"\tB10 =\t2.749767045580906e-05 [ao]\n",
"\tB11 =\t0.00011123931468954441 [ao]"
]
},
"execution_count": 3,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"water"
]
},
{
"cell_type": "code",
"execution_count": 4,
"id": "d1500d4e-b872-4b80-9ea4-231a1848e93e",
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"<Axes: xlabel='x [cm]', ylabel='y [cm]'>"
]
},
"execution_count": 4,
"metadata": {},
"output_type": "execute_result"
},
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAASUAAAEMCAYAAACCxKCJAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAHT5JREFUeJzt3X1UVGUeB/DvoDKANoMGMpCIGia+guLqDusRXTlLqa2es2fXt4RMIdzcDemosJVkbpHpmltR2dlF29bMLN/WCo+SZiWLiYyRmScMxdRBzZhBU1R49o+Os07MwNxh7swz8P2cM+c4l+e593cv3K/PvTPzjEYIIUBEJIkAXxdARHQ7hhIRSYWhRERSYSgRkVQYSkQkFYYSEUmFoUREUmEoEZFUGEpEJBWGEhFJxa9Caf/+/bj//vsRFRUFjUaDbdu2tdpn3759GDFiBLRaLWJjY7F+/XrV6yQi9/lVKF25cgXx8fEoLCx0qX11dTUmTZqE8ePHw2QyITs7G/PmzcOuXbtUrpSI3KXx1w/kajQabN26FVOnTnXaZsmSJXj//ffx5Zdf2pZNnz4ddXV1KC4u9kKVRKRUZ18XoKbS0lKkpKTYLUtNTUV2drbTPg0NDWhoaLA9b2pqwqVLl3DnnXdCo9GoVSqRXxJCoL6+HlFRUQgI8MyFV7sOJbPZjIiICLtlERERsFqtuHr1KoKDg5v1KSgowLJly7xVIlG7cPr0afTq1csj62rXoeSOvLw85OTk2J5bLBb07t0bp0+fhk6n82FlRPKxWq2Ijo7GHXfc4bF1tutQMhgMqK2ttVtWW1sLnU7ncJQEAFqtFlqtttlynU7HUCJywpO3Nvzq1TeljEYjSkpK7Jbt3r0bRqPRRxURUWv8KpQuX74Mk8kEk8kE4KeX/E0mE2pqagD8dOmVlpZma5+VlYVvv/0Wixcvxtdff41XXnkF77zzDhYuXOiL8onIFcKP7N27VwBo9khPTxdCCJGeni6Sk5Ob9UlISBCBgYGiX79+Yt26dYq2abFYBABhsVg8sxNE7Yga54ffvk/JW6xWK/R6PSwWC+8pEf2MGueHX12+EVH7x1AiIqkwlIhIKgwlIpIKQ4mIpMJQIiKpMJSISCoMJSKSCkOJiKTCUCIiqTCUiEgqDCUikgpDiYikwlAiIqkwlIhIKgwlIpIKQ4mIpMJQIiKpMJSISCoMJSKSCkOJiKTCUCIiqTCUiEgqDCUikgpDiYikwlAiIqkwlIhIKgwlIpIKQ4mIpMJQIiKpMJSISCoMJSKSCkOJiKTS2dcFUNt12qnxdQmkosbJwtcleBVHSkQkFYYSEUmFoUREUmEoEZFUGEpEJBW/C6XCwkL06dMHQUFBGD16NA4ePOi07fr166HRaOweQUFBXqyWiJTyq1DatGkTcnJykJ+fj8OHDyM+Ph6pqak4f/680z46nQ7nzp2zPU6dOuXFiolIKb8KpdWrVyMjIwNz5szBoEGD8NprryEkJARFRUVO+2g0GhgMBtsjIiLCixUTkVJ+E0rXr19HeXk5UlJSbMsCAgKQkpKC0tJSp/0uX76MmJgYREdHY8qUKTh69GiL22loaIDVarV7EJH3+E0oXbx4EY2Njc1GOhERETCbzQ77DBgwAEVFRdi+fTv+/e9/o6mpCUlJSfjuu++cbqegoAB6vd72iI6O9uh+EFHL/CaU3GE0GpGWloaEhAQkJydjy5YtCA8Px9q1a532ycvLg8VisT1Onz7txYqJyG8++xYWFoZOnTqhtrbWbnltbS0MBoNL6+jSpQuGDx+Oqqoqp220Wi20Wm2baiUi9/nNSCkwMBCJiYkoKSmxLWtqakJJSQmMRqNL62hsbERlZSUiIyPVKpOI2shvRkoAkJOTg/T0dIwcORKjRo3CmjVrcOXKFcyZMwcAkJaWhrvuugsFBQUAgKeffhq//OUvERsbi7q6OqxcuRKnTp3CvHnzfLkbRNQCvwqladOm4cKFC1i6dCnMZjMSEhJQXFxsu/ldU1ODgID/D/5++OEHZGRkwGw2o3v37khMTMSBAwcwaNAgX+0CuWjRdOeX2C1Z+Xashyshb9MIITrWZC0KWa1W6PV6WCwW6HQ6X5fjkL/Pp+RuALnK34NK5vmU1Dg//GqkRO2H2kHkbFv+HlAdAUOJvMqbYdTS9hlO8mIokep8HUSOcPQkL4YSqUbGMHKEoye5MJRIFZ4IpLHv5ypqv3/Sc23a3qLpVQwmCTCUyKPcDSOlAeTqOpQGFUdNvsdQIo9RGkieCCIl21ASUBw1+Q5DidpMxjBqabuuhhNHTb7BUKI2URJIvgqjn3MnnBhM3sNQIre5GkiyhNHPKQknBpP3MJTILa4Ekqxh9HOuhhODyTv8ZuoSkkd7CqTbuVKzv7z3yp8xlEiR9hpItzCYfI+Xb+Sy1k5Gfw6j27lyOcdLOfVwpEQu6SiBdLvW9okjJnUwlKhVHTGQbmEweR9DidqkPQfSLR1hH2XCUKIWtTQS6Egna0v7ytGSZzGUyCmebK7jsfIchhK5pSONkm7piPvsCwwlcoiXbY7xMk59DCVqhoHUMgaTuhhKRCQVl97RnZOTo3jFTzzxBHr06KG4H/kWR0muGft+rtN3fPPd3m3jUiitWbMGRqMRgYGBLq30008/xYIFCxhKRKSYy59927p1K3r27OlS2zvuuMPtgkhOHCU119Joidzn0j2ldevWQa/Xu7zStWvXIiIiwu2iyDd4k9ZzeCzd51IopaenQ6vVurzSmTNnomvXrm4XRXLhKMk5HhvPa9Orb5cvX4bVarV7kH/i/+yex2PqHsWhVF1djUmTJqFr167Q6/Xo3r07unfvjtDQUHTv3l2NGsmHOBJoHY+RZyme5O2BBx6AEAJFRUWIiIiARqNRoy4i6qAUh9KRI0dQXl6OAQMGqFEP+QAvM9TD9ywpp/jy7Re/+AVOnz6tRi0kGV6WuI7HynMUj5T+8Y9/ICsrC2fOnMGQIUPQpUsXu58PGzbMY8URUcejOJQuXLiAEydOYM6cObZlGo0GQghoNBo0NjZ6tEAi6lgUh9JDDz2E4cOHY+PGjbzRTUQepziUTp06hR07diA2ljfv2gPe5FYfb3Yro/hG969//WscOXJEjVpIIrxxqxyPmWcoHindf//9WLhwISorKzF06NBmN7p/+9vfeqw4Iup4FIdSVlYWAODpp59u9jNv3OguLCzEypUrYTabER8fj5deegmjRo1y2n7z5s148skncfLkSfTv3x8rVqzAxIkTVa2RiNyn+PKtqanJ6UPtQNq0aRNycnKQn5+Pw4cPIz4+HqmpqTh//rzD9gcOHMCMGTMwd+5cVFRUYOrUqZg6dSq+/PJLVeskIvf51XS4q1evRkZGBubMmYNBgwbhtddeQ0hICIqKihy2//vf/457770XixYtwsCBA7F8+XKMGDECL7/8spcrJyJXKQ6lP//5z3jxxRebLX/55ZeRnZ3tiZocun79OsrLy5GSkmJbFhAQgJSUFJSWljrsU1paatceAFJTU522B4CGhgbOfEDkQ4pD6b333sOvfvWrZsuTkpLw7rvveqQoRy5evIjGxsZmk8dFRETAbDY77GM2mxW1B4CCggLo9XrbIzo6uu3FE5HLFIfS999/73AWSp1Oh4sXL3qkKF/Ky8uDxWKxPfg5PyLvUhxKsbGxKC4ubrb8ww8/RL9+/TxSlCNhYWHo1KkTamtr7ZbX1tbCYDA47GMwGBS1BwCtVgudTmf3ICLvURxKOTk5WLx4MfLz8/Hxxx/j448/xtKlS5Gbm4uFCxeqUSMAIDAwEImJiSgpKbEta2pqQklJCYxGo8M+RqPRrj0A7N6922l7IvI9tz771tDQgGeeeQbLly8HAPTp0wevvvoq0tLSPF7g7XJycpCeno6RI0di1KhRWLNmDa5cuWL7cHBaWhruuusuFBQUAAAeffRRJCcn429/+xsmTZqEt99+G4cOHcLrr7+uap1E5D7FoQQA8+fPx/z583HhwgUEBwejW7dunq7LoWnTpuHChQtYunQpzGYzEhISUFxcbLuZXVNTg4CA/w/+kpKS8NZbb+GJJ57AX/7yF/Tv3x/btm3DkCFDvFIvESmnEUIIXxchM6vVCr1eD4vFIu39pU473Z+pgd+I6zktfQdcWz6Q2zhZ3lNUjfPDpXtKI0aMwA8//ODySseMGYMzZ864XRR5Dz+9rj4eY2VcunwzmUw4cuSIy1/DbTKZ0NDQ0KbCiKhjcvme0oQJE+DqlR4nfiMid7kUStXV1YpX3KtXL8V9iIhcCqWYmBi16yAJ7Z/0HG92u6ilm9ykjF/NEkDq4I1Y9fDYKsdQIiKpMJSoRbwsaR2PkWcxlAgALzPUwGPqHsWhlJ6ejv3796tRC0mKIwHneGw8T3EoWSwWpKSkoH///nj22Wf5zu12hP+zew6PpfsUh9K2bdtw5swZzJ8/H5s2bUKfPn1w33334d1338WNGzfUqJEkwBFBczwm6nDrnlJ4eDhycnJw5MgRlJWVITY2FrNnz0ZUVBQWLlyIb775xtN1ElEH0aYb3efOncPu3buxe/dudOrUCRMnTkRlZSUGDRqEF154wVM1khe1dNnBkcH/qTUjALkRSjdu3MB7772HyZMnIyYmBps3b0Z2djbOnj2LN954A3v27ME777zj8MsqiYhaoziUIiMjkZGRgZiYGBw8eBCHDh1CVlaW3Vwq48ePR2hoqCfrJC/iaKllHCWpS3EovfDCCzh79iwKCwuRkJDgsE1oaKhbH+IleTCYHGMgqU9xKM2ePRtBQUFq1EJ+pCMGU0fcZ1/gO7rJKf7P7zoeK89hKFGLeBn3E162eQ9DidqkIwRTR9hHmTCUqFWtjQTa80nb2r5xlOR5DCVySUcMJgaSb7j1ZZTUMa18O7bF74m7dRL7+xS6rgQsA0k9HCmRIq6cjP48amIg+R5DiRRrr8HEQJIDL9/ILa1dygH+cznnaoAykLyDoURucyWYAHnDSclojoHkPQwlahNXgwmQJ5yUXloykLyLoURtduuklT2cGEb+gaFEHqNk1ATYh4RaAeXuDXcGku8wlMijlI6abnEUHkqDyhOv+DGMfI+hRKpQOmpyxNtvK2AgyYGhRKpxd9TkbQwjuTCUSHW3n/SyBBSDSF4MJfIqX4+eGEbyYyiRT3hz9MQg8i8MJfI5R6HhblAxgPwfQ4mkxHDpuPxmloBLly5h1qxZ0Ol0CA0Nxdy5c3H58uUW+4wbNw4ajcbukZWV5aWKicgdfjNSmjVrlu1rwm/cuIE5c+YgMzMTb731Vov9MjIy7L6tNyQkRO1SiagN/CKUjh07huLiYnz++ecYOXIkAOCll17CxIkTsWrVKkRFRTntGxISAoPB4K1SiaiN/OLyrbS0FKGhobZAAoCUlBQEBASgrKysxb4bNmxAWFgYhgwZgry8PPz4448ttm9oaIDVarV7EJH3+MVIyWw2o2fPnnbLOnfujB49esBsNjvtN3PmTMTExCAqKgpffPEFlixZguPHj2PLli1O+xQUFGDZsmUeq52IlPFpKOXm5mLFihUttjl27Jjb68/MzLT9e+jQoYiMjMSECRNw4sQJ3H333Q775OXlIScnx/bcarUiOjra7RqISBmfhtJjjz2GBx98sMU2/fr1g8FgwPnz5+2W37x5E5cuXVJ0v2j06NEAgKqqKqehpNVqodVqXV4nEXmWT0MpPDwc4eHhrbYzGo2oq6tDeXk5EhMTAQAfffQRmpqabEHjCpPJBACIjIx0q14iUp9f3OgeOHAg7r33XmRkZODgwYP47LPPsGDBAkyfPt32ytuZM2cQFxeHgwcPAgBOnDiB5cuXo7y8HCdPnsSOHTuQlpaGsWPHYtiwYb7cHSJqgV+EEvDTq2hxcXGYMGECJk6ciDFjxuD111+3/fzGjRs4fvy47dW1wMBA7NmzB7/5zW8QFxeHxx57DL/73e/wn//8x1e7QEQu0AghhK+LkJnVaoVer4fFYoFOp/N1OQ512qnxdQmkosbJ8p6iapwffjNSIqKOgaFERFJhKBGRVPziHd3UMpnvORApxZESEUmFoUREUmEoEZFUGEpEJBWGEhFJhaFERFJhKBGRVBhKRCQVhhIRSYWhRERSYSgRkVQYSkQkFYYSEUmFoUREUmEoEZFUGEpEJBWGEhFJhaFERFJhKBGRVBhKRCQVhhIRSYWhRERSYSgRkVQYSkQkFYYSEUmFoUREUmEoEZFUGEpEJBWGEhFJhaFERFJhKBGRVBhKRCQVhhIRSYWhRERSYSgRkVT8JpSeeeYZJCUlISQkBKGhoS71EUJg6dKliIyMRHBwMFJSUvDNN9+oWygRtYnfhNL169fx+9//HvPnz3e5z/PPP48XX3wRr732GsrKytC1a1ekpqbi2rVrKlZKRG0i/My6deuEXq9vtV1TU5MwGAxi5cqVtmV1dXVCq9WKjRs3urw9i8UiAAiLxeJOuUTtmhrnR2dfh6JaqqurYTabkZKSYlum1+sxevRolJaWYvr06Q77NTQ0oKGhwfbcYrEAAKxWq7oFE/mhW+eFEMJj62y3oWQ2mwEAERERdssjIiJsP3OkoKAAy5Yta7Y8OjraswUStSPff/899Hq9R9bl01DKzc3FihUrWmxz7NgxxMXFeakiIC8vDzk5ObbndXV1iImJQU1NjccOui9YrVZER0fj9OnT0Ol0vi7HbdwPuVgsFvTu3Rs9evTw2Dp9GkqPPfYYHnzwwRbb9OvXz611GwwGAEBtbS0iIyNty2tra5GQkOC0n1arhVarbbZcr9f79R/PLTqdjvshkfayHwEBnnvNzKehFB4ejvDwcFXW3bdvXxgMBpSUlNhCyGq1oqysTNEreETkXX7zloCamhqYTCbU1NSgsbERJpMJJpMJly9ftrWJi4vD1q1bAQAajQbZ2dn461//ih07dqCyshJpaWmIiorC1KlTfbQXRNQav7nRvXTpUrzxxhu258OHDwcA7N27F+PGjQMAHD9+3PZqGQAsXrwYV65cQWZmJurq6jBmzBgUFxcjKCjI5e1qtVrk5+c7vKTzJ9wPuXA/nNMIT76WR0TURn5z+UZEHQNDiYikwlAiIqkwlIhIKgwlB9rLNCmXLl3CrFmzoNPpEBoairlz59q9hcKRcePGQaPR2D2ysrK8VPFPCgsL0adPHwQFBWH06NE4ePBgi+03b96MuLg4BAUFYejQofjggw+8VGnLlOzH+vXrmx13Ja8Sq2H//v24//77ERUVBY1Gg23btrXaZ9++fRgxYgS0Wi1iY2Oxfv16xdtlKDnQXqZJmTVrFo4ePYrdu3dj586d2L9/PzIzM1vtl5GRgXPnztkezz//vBeq/cmmTZuQk5OD/Px8HD58GPHx8UhNTcX58+cdtj9w4ABmzJiBuXPnoqKiAlOnTsXUqVPx5Zdfeq1mR5TuB/DTu7tvP+6nTp3yYsXNXblyBfHx8SgsLHSpfXV1NSZNmoTx48fDZDIhOzsb8+bNw65du5Rt2GPzDbRD3p4mxZO++uorAUB8/vnntmUffvih0Gg04syZM077JScni0cffdQLFTo2atQo8cgjj9ieNzY2iqioKFFQUOCw/R/+8AcxadIku2WjR48WDz/8sKp1tkbpfrj6t+YrAMTWrVtbbLN48WIxePBgu2XTpk0TqampirbFkZIHtDZNii+UlpYiNDQUI0eOtC1LSUlBQEAAysrKWuy7YcMGhIWFYciQIcjLy8OPP/6odrkAfhqhlpeX2x3HgIAApKSkOD2OpaWldu0BIDU11WfHHXBvPwDg8uXLiImJQXR0NKZMmYKjR496o1yP8dTvwm/e0S0zd6dJUZPZbEbPnj3tlnXu3Bk9evRosaaZM2ciJiYGUVFR+OKLL7BkyRIcP34cW7ZsUbtkXLx4EY2NjQ6P49dff+2wj9lsluq4A+7tx4ABA1BUVIRhw4bBYrFg1apVSEpKwtGjR9GrVy9vlN1mzn4XVqsVV69eRXBwsEvr6TAjpdzc3GY3En/+cPYHIxO19yMzMxOpqakYOnQoZs2ahX/961/YunUrTpw44cG9oJ8zGo1IS0tDQkICkpOTsWXLFoSHh2Pt2rW+Ls3rOsxIScZpUtzh6n4YDIZmN1Vv3ryJS5cu2ep1xejRowEAVVVVuPvuuxXXq0RYWBg6deqE2tpau+W1tbVOazYYDIrae4M7+/FzXbp0wfDhw1FVVaVGiapw9rvQ6XQuj5KADhRK7WWaFFf3w2g0oq6uDuXl5UhMTAQAfPTRR2hqarIFjStMJhMA2IWtWgIDA5GYmIiSkhLbTA5NTU0oKSnBggULHPYxGo0oKSlBdna2bdnu3bthNBpVr9cZd/bj5xobG1FZWYmJEyeqWKlnGY3GZm/HcOt3ofQufEdw6tQpUVFRIZYtWya6desmKioqREVFhaivr7e1GTBggNiyZYvt+XPPPSdCQ0PF9u3bxRdffCGmTJki+vbtK65eveqLXRBCCHHvvfeK4cOHi7KyMvHpp5+K/v37ixkzZth+/t1334kBAwaIsrIyIYQQVVVV4umnnxaHDh0S1dXVYvv27aJfv35i7NixXqv57bffFlqtVqxfv1589dVXIjMzU4SGhgqz2SyEEGL27NkiNzfX1v6zzz4TnTt3FqtWrRLHjh0T+fn5okuXLqKystJrNTuidD+WLVsmdu3aJU6cOCHKy8vF9OnTRVBQkDh69KivdkHU19fb/vYBiNWrV4uKigpx6tQpIYQQubm5Yvbs2bb23377rQgJCRGLFi0Sx44dE4WFhaJTp06iuLhY0XYZSg6kp6cLAM0ee/futbUBINatW2d73tTUJJ588kkREREhtFqtmDBhgjh+/Lj3i7/N999/L2bMmCG6desmdDqdmDNnjl2wVldX2+1XTU2NGDt2rOjRo4fQarUiNjZWLFq0yOvf5PLSSy+J3r17i8DAQDFq1Cjx3//+1/az5ORkkZ6ebtf+nXfeEffcc48IDAwUgwcPFu+//75X63VGyX5kZ2fb2kZERIiJEyeKw4cP+6Dq/9u7d6/D8+BW3enp6SI5OblZn4SEBBEYGCj69etnd464ilOXEJFUOsyrb0TkHxhKRCQVhhIRSYWhRERSYSgRkVQYSkQkFYYSEUmFoUREUmEokVROnjxpm+3A0x9m/rnbp6C9/bNz5FsMJZLSnj17UFJSouo2pk2bhnPnzvn0w7vUXIeZJYD8y5133ok777xT1W0EBwcjODgYgYGBqm6HlOFIiVRz4cIFGAwGPPvss7ZlBw4cQGBgoFujoKKiIgwePBharRaRkZF204BoNBqsXbsWkydPRkhICAYOHIjS0lJUVVVh3Lhx6Nq1K5KSkjhZnR9gKJFqwsPDUVRUhKeeegqHDh1CfX09Zs+ejQULFmDChAmK1vXqq6/ikUceQWZmJiorK7Fjxw7ExsbatVm+fDnS0tJgMpkQFxeHmTNn4uGHH0ZeXh4OHToEIYTL8xmRD7VxdgOiVv3xj38U99xzj5g5c6YYOnSouHbtmtO2t6ZTqaiosFseFRUlHn/8caf9AIgnnnjC9ry0tFQAEP/85z9tyzZu3CiCgoKa9fX1N7iQPY6USHWrVq3CzZs3sXnzZmzYsAFarVZR//Pnz+Ps2bOtjq6GDRtm+/etCeyHDh1qt+zatWuwWq2Ktk/exVAi1Z04cQJnz55FU1MTTp48qbi/q/M7d+nSxfZvjUbjdFlTU5PiGsh7GEqkquvXr+OBBx7AtGnTsHz5csybN6/Fb4l15I477kCfPn1Uf4sAyYFvCSBVPf7447BYLHjxxRfRrVs3fPDBB3jooYewc+dORet56qmnkJWVhZ49e+K+++5DfX09PvvsM/zpT39SqXLyFY6USDX79u3DmjVr8Oabb0Kn0yEgIABvvvkmPvnkE7z66quK1pWeno41a9bglVdeweDBgzF58mR88803KlVOvsQ5ukkqJ0+eRN++fVFRUaH6x0xuGTduHBISErBmzRqvbI9axpESSSkpKQlJSUmqbmPDhg3o1q0bPvnkE1W3Q8pwpERSuXnzpu0VOq1Wi+joaNW2VV9fb/tG19DQUISFham2LXIdQ4mIpMLLNyKSCkOJiKTCUCIiqTCUiEgqDCUikgpDiYikwlAiIqkwlIhIKv8DThObj2JbR1QAAAAASUVORK5CYII=",
"text/plain": [
"<Figure size 258.065x259.74 with 1 Axes>"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"box = openmc.model.RectangularPrism(1.2, 1.2, boundary_type='reflective')\n",
"top = openmc.ZPlane(z0=300.0, boundary_type='vacuum')\n",
"bottom = openmc.ZPlane(z0=0.0, boundary_type='vacuum')\n",
"\n",
"fuel_surf = openmc.ZCylinder(r=0.39)\n",
"clad_surf = openmc.ZCylinder(r=0.45)\n",
"\n",
"fuel = openmc.Cell(region=-fuel_surf & -box, fill=uo2)\n",
"clad = openmc.Cell(region=+fuel_surf & -clad_surf & -box, fill=zirconium)\n",
"h2o = openmc.Cell(region=+clad_surf & -box, fill=water)\n",
"\n",
"model.geometry = openmc.Geometry(openmc.Universe(cells=[fuel, clad, h2o]))\n",
"\n",
"model.geometry.root_universe.plot(width=(2, 2))"
]
},
{
"cell_type": "code",
"execution_count": 5,
"id": "61952a0e-4f29-42ec-afdc-689a4f715c36",
"metadata": {},
"outputs": [],
"source": [
"model.settings.particles = 1000\n",
"model.settings.inactive = 10\n",
"model.settings.batches = 50"
]
},
{
"cell_type": "code",
"execution_count": 7,
"id": "0780ae81-3a39-4bc3-9ae5-71895c6f858b",
"metadata": {},
"outputs": [],
"source": [
"def change_enrichment(wo):\n",
" uo2.remove_element('U')\n",
" uo2.add_element('U', 1.0, enrichment=wo)"
]
},
{
"cell_type": "markdown",
"id": "f4d61e7d-c71b-4bff-8b7f-89d9afe680f1",
"metadata": {},
"source": [
"When running a criticality search, we can specify the desired maximum uncertainty in $k$ with `sigma_final` and the tolerance to which we want $k$ to match unity with `k_tol`."
]
},
{
"cell_type": "code",
"execution_count": 9,
"id": "8bbd54dc-8e3c-42e9-be32-6d76be788656",
"metadata": {
"scrolled": true
},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Iteration 1: batches=50, x=1, keff=0.80469+/-0.00354\n",
"Iteration 2: batches=50, x=5, keff=1.29765+/-0.00437\n",
"Iteration 3: batches=30, x=2.58479, keff=1.11761+/-0.00486\n",
"Iteration 4: batches=30, x=2.3531, keff=1.09803+/-0.00459\n",
"Iteration 5: batches=30, x=2.17815, keff=1.08116+/-0.00583\n",
"Iteration 6: batches=30, x=1.061, keff=0.82630+/-0.00817\n",
"Iteration 7: batches=30, x=1.87998, keff=1.03062+/-0.00446\n",
"Iteration 8: batches=44, x=1.80779, keff=1.01285+/-0.00527\n",
"Iteration 9: batches=106, x=1.77392, keff=1.00482+/-0.00283\n",
"Iteration 10: batches=253, x=1.75563, keff=1.00593+/-0.00173\n",
"Iteration 11: batches=201, x=1.73126, keff=0.99987+/-0.00195\n",
"Iteration 12: batches=987, x=1.72707, keff=1.00016+/-0.00083\n"
]
}
],
"source": [
"result = model.keff_search(change_enrichment, x0=1.0, x1=5.0, output=True, k_tol=0.00100, sigma_final=0.001)"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "062283ed-3854-49d2-8fea-d825dc82a8a5",
"metadata": {},
"outputs": [],
"source": []
},
{
"cell_type": "markdown",
"id": "27c3f865-81c5-42b1-8456-406a730d213e",
"metadata": {},
"source": [
"## Geometry\n",
"\n",
"Criticality searches work for any change to a mutable object - in this next example, we'll take a sphere of uranium enclosed in a metal aluminum shell of 1 cm thickness, and determine the critical radius."
]
},
{
"cell_type": "code",
"execution_count": 10,
"id": "a237cca8-aa86-48ab-a400-8de0ed9570fa",
"metadata": {},
"outputs": [],
"source": [
"model = openmc.Model()\n",
"\n",
"fuel = openmc.Material()\n",
"fuel.add_nuclide('U235', 1.0)\n",
"fuel.set_density('g/cm3', 19.1)\n",
"\n",
"al = openmc.Material()\n",
"al.add_element('Al', 1.0)\n",
"al.set_density('g/cm3', 2.7)\n",
"dr = 1.0\n",
"\n",
"in_sphere = openmc.Sphere(r=10.0)\n",
"out_sphere = openmc.Sphere(r=in_sphere.r + dr, boundary_type='vacuum')\n",
"\n",
"fc = openmc.Cell(region=-in_sphere, fill=fuel)\n",
"oc = openmc.Cell(region=+in_sphere & -out_sphere, fill=al)\n",
"\n",
"model.geometry = openmc.Geometry(openmc.Universe(cells=[fc, oc]))"
]
},
{
"cell_type": "code",
"execution_count": 13,
"id": "6a259fd7-a0ee-4443-a8cd-d41d4f7c5abc",
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"<Axes: xlabel='x [cm]', ylabel='y [cm]'>"
]
},
"execution_count": 13,
"metadata": {},
"output_type": "execute_result"
},
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAARYAAAEMCAYAAAABAJmyAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAHLlJREFUeJzt3XtQVPfdBvBnMbJykftl2QkqBJViDGhSKTR544URqdqx7dhUq6Jp2mgkFtG00njDG15ay5gYTKYqjuOoTTI6ad6YqTJabUUdTbCOYzqCGIyyiFpYcIZF4bx/+LKysAvnLGf3XPb5zOyMe/bs7nc37JPv71x+xyAIggAiIhn5KV0AEekPg4WIZMdgISLZMViISHYMFiKSHYOFiGTHYCEi2TFYiEh2DBYikh2DhYhkp6lgOX36NKZPnw6z2QyDwYCjR486PD5//nwYDAaH25QpU5QplsiHaSpYHj58iNTUVOzcudPlOlOmTEFdXZ39dvDgQS9WSEQA8IzSBUiRk5ODnJycXtcxGo0wmUxeqoiInNFUsIhx6tQpxMTEIDw8HBMnTsSGDRsQGRnpcn2bzQabzWa/39HRgQcPHiAyMhIGg8EbJRNphiAIaG5uhtlshp+f6wGProJlypQp+OlPf4qEhARUV1fjD3/4A3JyclBRUYEBAwY4fU5xcTGKioq8XCmRtt26dQvPPvusy8cNWp2PxWAw4MiRI5gxY4bLdW7cuIHnnnsOJ06cwKRJk5yu071jaWpqwpAhQ3Dr1i2EhITIXTaRplmtVsTHx6OxsRGhoaEu19NVx9JdYmIioqKiUFVV5TJYjEYjjEZjj+UhISEMFiIX+tpMoKm9QlJ99913uH//PuLi4pQuhcinaKpjaWlpQVVVlf1+TU0NKisrERERgYiICBQVFeFnP/sZTCYTqqur8bvf/Q5JSUnIzs5WsGoi36OpYLl48SImTJhgv19QUAAAyM3NRWlpKf79739j3759aGxshNlsxuTJk7F+/XqnQx0i8hzNbrz1FKvVitDQUDQ1NXEbC1E3Yn8fut7GQkTKYLAQkewYLEQkOwYLEcmOwUJEsmOwEJHsGCxEJDsGCxHJjsFCRLJjsBCR7BgsRCQ7BgsRyY7BQkSyY7AQkewYLEQkOwYLEcmOwUJEsmOwEJHsGCxEJDtNTaZN+uGpqZZ5WVx1YLCQxzkLkcLBNzzyXsXNiT2WMWy8j8FCsuseJJ4KEWecvVf3sGHQeJ6mtrGcPn0a06dPh9lshsFgwNGjRx0eFwQBq1evRlxcHAICApCVlYXr168rU6yPEQTBfiscfMPhprTu9XStlTxDU8Hy8OFDpKamYufOnU4f37p1K3bs2IFdu3bh/PnzCAoKQnZ2NlpbW71cqW9wFSZqx5DxPM1esMxgMODIkSOYMWMGgCd/5GazGcuWLcPy5csBAE1NTYiNjUVZWRl+8YtfiHpdXrCsd13/XLQQIlJ0HTJxuOSc2N+Hbrax1NTUwGKxICsry74sNDQU6enpqKiocBksNpsNNpvNft9qtXq8Vi3qDBTPh4nY/8/J/8Pv+tk6Q4YB4x7dBIvFYgEAxMbGOiyPjY21P+ZMcXExioqKPFqblnkmUFyHx//8b2HfzxaAM9OKe1mj/2HQ+XkZMO7RTbC4q7Cw0H5xeeBJxxIfH69gReogX6D0DBEx4dEbg8H1a7gOHfeCgQHjHt0Ei8lkAgDU19cjLi7Ovry+vh5paWkun2c0GmE0Gj1dnmbIEyhPw6S/ISKVs9DpGTbSw4EBI41ugiUhIQEmkwnl5eX2ILFarTh//jwWLVqkbHEa0bl3x81n2//l7TDpS9ew6W/IdA0YhotrmgqWlpYWVFVV2e/X1NSgsrISERERGDJkCPLz87FhwwYMHz4cCQkJWLVqFcxms33PETnnfpei3jBxRa6QKRx8g91LLzS1u/nUqVOYMGFCj+W5ubkoKyuDIAhYs2YNPvroIzQ2NuLll1/GBx98gBEjRoh+D1/b3exel/LkT0YrYSKGY8hICwpf6l7E/j40FSze4CvB4l6Xor9A6c7dgPGV7oXB4iZfCBbpXYr+A6W7/gSMnsOFweImPQeL9C7F9wKlO3cCRs/dC4PFTXoNFne6FF8OlO6eBoxvdy9ifx+aOgmR3CMtVAQwVHoyGIBXPi9E5/cjRudJjr6IHUs3eutYpIYKA6VvUrsXPXUu7FhIQqiwS5FCavfii50LO5Zu9NKxSAkVBor7pHQveuhc2LH4MIaK9zh2L73zpc6FwaIz4kKFQx85SRka+Uq4cCjUjZaHQmJDhYHiOWKHRlodFnEo5GMYKuogdmik986FwaIDDBV1YbgwWDSPoaJOvh4uDBYNY6iomy+HC4NFoxgq2iAlXPSEwaJbDBW1EBsueupaGCwa1He3wlBRm6fh4pqehkQMFo3p34TXpDzfGBIxWHSH3Ypa+dKQiMGiIRwCaZ+YcNHDkIjBohEMFf0Qu71FyxgsRIrR75CIwaIB7Fb0R+9DIl0Fy9q1a2EwGBxuycnJSpfVLwwV/dLzkEhTl1gVY9SoUThx4oT9/jPP6O4jku4I6G2aBUEQNDfFgq46FuBJkJhMJvstKiqq1/VtNhusVqvDTS3Yreif2CGR1uguWK5fvw6z2YzExET88pe/RG1tba/rFxcXIzQ01H6Lj4/3UqX9xVDRCzFDIq1ta9FVsKSnp6OsrAxffvklSktLUVNTg1deeQXNzc0un1NYWIimpib77datW16s2DUeYeuL9NO16HpqysbGRgwdOhTbt2/Hr371K1HPUcvUlL0HC7sVPXoyreVml4+rYTpLTk0JICwsDCNGjEBVVZXSpUjCbsWX6aNr0XWwtLS0oLq6GnFxcUqXIiN2K3qlp20tugqW5cuX4x//+Adu3ryJs2fP4ic/+QkGDBiAWbNmKV2aaOxWSA9di66C5bvvvsOsWbMwcuRI/PznP0dkZCTOnTuH6OhopUuTCbsVvRPTtWiBro4eO3TokNIlEMnE9UFzWjhgTlcdi9ZxTxABfXctWhgOMViISHYMFiLVcr0RV+17hxgsKsFhEHWl9eEQg4WIZMdgIVI1bQ6HRJ0rVFBQIPmFV65ciYiICLeKUpIS5wpxGESu9HX+0OaW57xYjfjfh6jjWEpKSpCRkQF/f39Rb/7Pf/4TeXl5mgwWIuo/0QfIHTlyBDExMaLWHTx4sNsFEVF3vc8wp0aitrHs3bsXoaGhol/0ww8/RGxsrNtFUScOg3xdX3uH1LqdRVSw5Obmwmg0in7R2bNnIygoyO2ifAlPOqT+UOvfTr/OFWppaUFHR4fDMiUnRyIidZC8u7mmpgZTp05FUFAQQkNDER4ejvDwcISFhSE8PNwTNRKRxkjuWObMmQNBELBnzx7Exsaq/ixLIn3Q1gZcycFy+fJlXLp0CSNHjvREPWTHDbf0ROcGXFfHs6hxGgXJQ6Hvf//7qpnJXuu44ZbkoMa/Ickdy1/+8hcsXLgQt2/fxvPPP4+BAwc6PP7CCy/IVhwRaZPkYGloaEB1dTUWLFhgX2YwGOztWHt7u6wFEpH2SA6W119/HWPGjMHBgwe58ZaInJIcLN9++y0+++wzJCUleaIeItIByRtvJ06ciMuXL3uiFiLqlToP33dGcscyffp0LF26FFeuXMHo0aN7bLz98Y9/LFtxvou7mslRX7uc1UZysCxcuBAAsG7duh6PqWXj7c6dO7Ft2zZYLBakpqbivffew7hx45Qui8hnSB4KdXR0uLypIVQOHz6MgoICrFmzBl999RVSU1ORnZ2Nu3fvKl0akc/Q3dSU27dvx69//WssWLAAKSkp2LVrFwIDA7Fnzx6lSyPyGZKDZcmSJdixY0eP5e+//z7y8/PlqMltbW1tuHTpErKysuzL/Pz8kJWVhYqKCqfPsdlssFqtDjci6h/JwfLpp5/ihz/8YY/lmZmZ+OSTT2Qpyl337t1De3t7j0mmYmNjYbFYnD6nuLgYoaGh9lt8fLw3SiXSNcnBcv/+faezyYWEhODevXuyFOVNhYWFaGpqst94HhRR/0kOlqSkJHz55Zc9lh87dgyJiYmyFOWuqKgoDBgwAPX19Q7L6+vrYTKZnD7HaDQiJCTE4UZE/SN5d3NBQQHy8vLQ0NCAiRMnAgDKy8vxpz/9CSUlJXLXJ4m/vz9efPFFlJeXY8aMGQCe7MUqLy9HXl6eorUR+RK3zhWy2WzYuHEj1q9fDwAYNmwYSktLMW/ePNkLlKqgoAC5ubl46aWXMG7cOJSUlODhw4cOJ01qgSA8OSiKSItEXbDMlYaGBgQEBCA4OFjOmvrt/ffftx8gl5aWhh07diA9PV3Uc715wTJeqIzEUsuFy2S9YJkr0dHR/Xm6x+Tl5XHoQ6QgURtvx44di//+97+iX/Tll1/G7du33S6KiLRNVMdSWVmJy5cvi75kamVlJWw2W78KIyLtEj0UmjRpkuirrnHyJyLfJipYampqJL/ws88+K/k5RKQPooJl6NChnq7DJxkMBhQ3J7rcM8RdziSGt/YISaG7s5v1w4Az04qVLoJUoK9dzWrEYCEi2TFYiEh2DBYikp3kYMnNzcXp06c9UYtP6tyA64r7J1yQL1DjhlvAjWBpampCVlYWhg8fjk2bNvEIW4/iBlxfp8UNt4AbwXL06FHcvn0bixYtwuHDhzFs2DDk5OTgk08+waNHjzxRIxFpjFvbWKKjo1FQUIDLly/j/PnzSEpKwty5c2E2m7F06VJcv35d7jp9GodDpDX92nhbV1eH48eP4/jx4xgwYAB+9KMf4cqVK0hJScGf//xnuWrUvd63s3A45KvUMlWCOyQHy6NHj/Dpp59i2rRpGDp0KD7++GPk5+fjzp072LdvH06cOIG//vWvTi9oRkS+QfJ8LHFxcejo6MCsWbNw4cIFpKWl9VhnwoQJCAsLk6E86sTD+0lLJM8gt3//fsycORODBg3yVE2K8uYMct1xRjnqpNZhkMdmkJs7d26/CiMi/eORtxrCvUOkFQwWFeHeIQLUOwySgsGiMexaSAsYLJrCrkXvtHoIf3e6CpZhw4bBYDA43DZv1tZ/pL5OSgTYtfgyLQyDAJ0FCwCsW7cOdXV19tvbb7+tdEkyY9eiV3rpVoB+XrBMjQYPHuzyAvBa0ddcuAAPmPNFWulWAB12LJs3b0ZkZCTGjBmDbdu24fHjx72ub7PZYLVaHW7qx65Fb/TUrQA661iWLFmCsWPHIiIiAmfPnkVhYSHq6uqwfft2l88pLi5GUVGRF6sUh10LdaWlbgXo50XhvWHFihXYsmVLr+tcu3YNycnJPZbv2bMHb775JlpaWmA0Gp0+12azOVy10Wq1Ij4+XpFD+rvr/RB/ABDwyueFDBeNE9OtqCVYxB7Sr/pgaWhowP3793tdJzExEf7+/j2WX716Fc8//zy++eYbjBw5UtT7KXmukDNiwoXnEGmXlkIF8OC5Qt4WHR2N6Ohot55bWVkJPz8/xMTEyFyVunBIRGqjm423FRUVKCkpweXLl3Hjxg0cOHAAS5cuxZw5cxAeHq50eW7r+7iWJxty1d13kjNa61ak0E2wGI1GHDp0CK+++ipGjRqFjRs3YunSpfjoo4+ULq3fxIYLaYeeQwXQwFBIrLFjx+LcuXNKl6EoDolILXTTsegdh0T6ofduBWCwaArDRft8IVQABosOcXuLWunt6NreMFg0RszZzwDPgNYqPXQrAINFkzgk0h5fGQJ1YrDoFsNFLXxpCNSJwaJR4oZEDBeliQ0VPXUrAINF0xgu6uaroQIwWDSP4aJOvhwqAINFFxgu6uLroQIwWHRDargwYOTX+b36eqgAGpiPxdvUNh+LVH3P32Jfk5NEyUjKnh8th4rY3wc7Fp0RewAdh0by8ZVQkYLBokPuhAsDRjopQx/Ad0IF4FCoB60PhboSPywCODSSRupBb3oJFQ6FSELnArB7EUdqlwLoJ1SkYMfSjZ46lk7SOheA3Ytz7hyar7dQYcdCdtI6F4DdiyN3uhRAf6EiBTuWbvTYsXTlbvcC+N60l52/DAbKU7q5rpC36T1YAHfCBfClgHE3UAB9hwrAYHGbLwRLJwaMIwZK3xgsbvKlYAHcDRega8AA2g2Zrn/97s6Z4iuhAjBY3OZrwdLJ/YABtBYycoQJ4FuB0kl3e4U2btyIzMxMBAYGIiwszOk6tbW1mDp1KgIDAxETE4N33nkHjx8/9m6hGiV9z5HDs3Fm2ub/vz3dm6Sm/2V1relprQwVT9HMBcva2towc+ZMZGRkYPfu3T0eb29vx9SpU2EymXD27FnU1dVh3rx5GDhwIDZt2qRAxdpjMBiwueW5fnYvhi4/WMdO5un7uF2iKM4CTa6pIRko4mhuKFRWVob8/Hw0NjY6LD927BimTZuGO3fuIDY2FgCwa9cu/P73v0dDQwP8/f2dvp7NZoPNZrPft1qtiI+P97mhkDP9CxiXr+o0bDqJDZ3e/mo9Mb8sA+UJsUMhzXQsfamoqMDo0aPtoQIA2dnZWLRoEa5evYoxY8Y4fV5xcTGKioq8VaamyNPB9HjVXn74vYdOV96anJqB4h7dBIvFYnEIFQD2+xaLxeXzCgsLUVBQYL/f2bHQU54JGKfvpJrZ7Bko/aNosKxYsQJbtmzpdZ1r164hOTnZYzUYjUYYjUaPvb6edAYM4KlhkrIYJvJRNFiWLVuG+fPn97pOYqK4PRUmkwkXLlxwWFZfX29/jOSll5BhmHiGosESHR2N6OhoWV4rIyMDGzduxN27dxETEwMAOH78OEJCQpCSkiLLe5BzWgsZhonnaWYbS21tLR48eIDa2lq0t7ejsrISAJCUlITg4GBMnjwZKSkpmDt3LrZu3QqLxYKVK1di8eLFHOp4UdeQAdQRNAwS79PM7ub58+dj3759PZafPHkS48ePBwB8++23WLRoEU6dOoWgoCDk5uZi8+bNeOYZ8fnpq0feepMnw4Yh4lk8pN9NDBYi13R3SD8RaQeDhYhkx2AhItkxWIhIdgwWIpIdg4WIZMdgISLZMViISHYMFiKSHYOFiGTHYCEi2TFYiEh2DBYikh2DhYhkx2AhItkxWIhIdgwWIpIdg4WIZMdgISLZMViISHYMFiKSHYOFiGSnmWDZuHEjMjMzERgYiLCwMKfrGAyGHrdDhw55t1Ai0s6VENva2jBz5kxkZGRg9+7dLtfbu3cvpkyZYr/vKoSIyHM0EyxFRUUAgLKysl7XCwsL40XgiRSmmWARa/HixXjjjTeQmJiIhQsXYsGCBTAYDC7Xt9lssNls9vtNTU0AnlzxjYgcdf4u+rqAqq6CZd26dZg4cSICAwPx97//HW+99RZaWlqwZMkSl88pLi62d0NdxcfHe7JUIk1rbm5GaGioy8cVvXbzihUrsGXLll7XuXbtGpKTk+33y8rKkJ+fj8bGxj5ff/Xq1di7dy9u3brlcp3uHUtHRwcePHiAyMjIXjsdT7JarYiPj8etW7d8/vrR/C6eUMv3IAgCmpubYTab4efnet+Poh3LsmXLMH/+/F7XSUxMdPv109PTsX79ethsNhiNRqfrGI3GHo+pZYNvSEiIT/+YuuJ38YQavofeOpVOigZLdHQ0oqOjPfb6lZWVCA8PdxkqROQZmtnGUltbiwcPHqC2thbt7e2orKwEACQlJSE4OBh/+9vfUF9fjx/84AcYNGgQjh8/jk2bNmH58uXKFk7kiwSNyM3NFQD0uJ08eVIQBEE4duyYkJaWJgQHBwtBQUFCamqqsGvXLqG9vV3Zwt3Q2toqrFmzRmhtbVW6FMXxu3hCa9+DohtviUifNHNIPxFpB4OFiGTHYCEi2TFYiEh2DBaVETM9RG1tLaZOnYrAwEDExMTgnXfewePHj71bqBfs3LkTw4YNw6BBg5Ceno4LFy4oXZLHnT59GtOnT4fZbIbBYMDRo0cdHhcEAatXr0ZcXBwCAgKQlZWF69evK1NsLxgsKtM5PcSiRYucPt7e3o6pU6eira0NZ8+exb59+1BWVobVq1d7uVLPOnz4MAoKCrBmzRp89dVXSE1NRXZ2Nu7evat0aR718OFDpKamYufOnU4f37p1K3bs2IFdu3bh/PnzCAoKQnZ2NlpbW71caR8U3t1NLuzdu1cIDQ3tsfyLL74Q/Pz8BIvFYl9WWloqhISECDabzYsVeta4ceOExYsX2++3t7cLZrNZKC4uVrAq7wIgHDlyxH6/o6NDMJlMwrZt2+zLGhsbBaPRKBw8eFCBCl1jx6IxFRUVGD16NGJjY+3LsrOzYbVacfXqVQUrk09bWxsuXbqErKws+zI/Pz9kZWWhoqJCwcqUVVNTA4vF4vC9hIaGIj09XXXfC4NFYywWi0OoALDft1gsSpQku3v37qG9vd3p59TLZ3RH52fXwvfCYPGCFStWOJ2Pt+vtm2++UbpMItlo5iRELZNzegiTydRj70h9fb39MT2IiorCgAED7J+rU319vW4+ozs6P3t9fT3i4uLsy+vr65GWlqZQVc6xY/GC6OhoJCcn93rz9/cX9VoZGRm4cuWKw96R48ePIyQkBCkpKZ76CF7l7++PF198EeXl5fZlHR0dKC8vR0ZGhoKVKSshIQEmk8nhe7FarTh//rzqvhd2LCrT1/QQkydPRkpKCubOnYutW7fCYrFg5cqVWLx4sa7mnSkoKEBubi5eeukljBs3DiUlJXj48CEWLFigdGke1dLSgqqqKvv9mpoaVFZWIiIiAkOGDEF+fj42bNiA4cOHIyEhAatWrYLZbMaMGTOUK9oZpXdLkaO+pocQBEG4efOmkJOTIwQEBAhRUVHCsmXLhEePHilXtIe89957wpAhQwR/f39h3Lhxwrlz55QuyeNOnjzp9L9/bm6uIAhPdjmvWrVKiI2NFYxGozBp0iThP//5j7JFO8FpE4hIdtzGQkSyY7AQkewYLEQkOwYLEcmOwUJEsmOwEJHsGCxEJDsGCxHJjsFCirl586b97G5Pn0RXVlZmf6/8/HyPvhcxWEgFTpw44XBinSe89tprqKurU93JenrFkxBJcZGRkYiMjPToewQEBCAgIED0WeTUP+xYSBYNDQ0wmUzYtGmTfdnZs2fh7+/vVjeyZ88ejBo1CkajEXFxccjLy7M/ZjAY8OGHH2LatGkIDAzE9773PVRUVKCqqgrjx49HUFAQMjMzUV1dLctnI+kYLCSL6Oho7NmzB2vXrsXFixfR3NyMuXPnIi8vD5MmTZL0WqWlpVi8eDF+85vf4MqVK/jss8+QlJTksM769esxb948VFZWIjk5GbNnz8abb76JwsJCXLx4EYIgOIQReZnCZ1eTzrz11lvCiBEjhNmzZwujR48WWltbXa5bU1MjABC+/vprh+Vms1l49913XT4PgLBy5Ur7/YqKCgGAsHv3bvuygwcPCoMGDerx3FdffVX47W9/K/4DkVvYsZCs/vjHP+Lx48f4+OOPceDAAcmTT929exd37tzps8t54YUX7P/unFx69OjRDstaW1thtVolvT/Jg8FCsqqursadO3fQ0dGBmzdvSn5+QECAqPUGDhxo/7fBYHC5rKOjQ3IN1H8MFpJNW1sb5syZg9deew3r16/HG2+8IfnKhYMHD8awYcM8vvuZPIu7m0k27777LpqamrBjxw4EBwfjiy++wOuvv47PP/9c0uusXbsWCxcuRExMDHJyctDc3Ix//etfePvttz1UOcmNHQvJ4tSpUygpKcH+/fsREhICPz8/7N+/H2fOnEFpaamk18rNzUVJSQk++OADjBo1CtOmTVPlhc/JNc55S4q5efMmEhIS8PXXX3vtujjjx49HWloaSkpKvPJ+voodCykuMzMTmZmZHn2PAwcOIDg4GGfOnPHo+9AT7FhIMY8fP7bvOTIajYiPj/fYezU3N9uvrBgWFoaoqCiPvRcxWIjIAzgUIiLZMViISHYMFiKSHYOFiGTHYCEi2TFYiEh2DBYikh2DhYhk939kSvtiZZl3IQAAAABJRU5ErkJggg==",
"text/plain": [
"<Figure size 258.065x259.74 with 1 Axes>"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"model.plot(width=(30, 30))"
]
},
{
"cell_type": "markdown",
"id": "689982e2-4d85-4014-9ed2-8960746b5eb9",
"metadata": {},
"source": [
"In order to change the radius of the fuel sphere, we want to modify the `in_sphere.r` parameter. Note that we have to also update the `out_sphere.r` parameter! When we use `r=in_sphere.r + dr`, that only applies a value (not a reference to `in_sphere.r`), so it also needs to be updated by us."
]
},
{
"cell_type": "code",
"execution_count": 14,
"id": "40c467d1-d431-4edb-8f8b-098bccaa2945",
"metadata": {},
"outputs": [],
"source": [
"def change_r(r):\n",
" in_sphere.r = r\n",
" out_sphere.r = in_sphere.r + dr"
]
},
{
"cell_type": "code",
"execution_count": 22,
"id": "c765ce8c-eac4-4300-b1b3-6a907f7f777b",
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Iteration 1: batches=100, x=12, keff=1.34877+/-0.00093\n",
"Iteration 2: batches=100, x=15, keff=1.53733+/-0.00077\n",
"Iteration 3: batches=30, x=9.6331, keff=1.15725+/-0.00257\n",
"Iteration 4: batches=30, x=9.84362, keff=1.17560+/-0.00210\n",
"Iteration 5: batches=30, x=9.95447, keff=1.18665+/-0.00232\n",
"Iteration 6: batches=30, x=10.1907, keff=1.20553+/-0.00153\n",
"Iteration 7: batches=37, x=10.1236, keff=1.20029+/-0.00154\n",
"Iteration 8: batches=129, x=10.1225, keff=1.19880+/-0.00075\n",
"Iteration 9: batches=112, x=10.1314, keff=1.19972+/-0.00076\n"
]
}
],
"source": [
"model.settings.particles = 10000\n",
"model.settings.inactive = 10\n",
"model.settings.batches = 100\n",
"\n",
"result = model.keff_search(change_r, x0=12, x1=15, target=1.2, output=True, k_tol = 0.00100, sigma_final = 0.00100)"
]
},
{
"cell_type": "code",
"execution_count": 19,
"id": "0a28df86-cdd4-4a3e-8dc2-e3055d055b7c",
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"<Axes: xlabel='x [cm]', ylabel='y [cm]'>"
]
},
"execution_count": 19,
"metadata": {},
"output_type": "execute_result"
},
{
"data": {
"image/png": "iVBORw0KGgoAAAANSUhEUgAAARYAAAEMCAYAAAABAJmyAAAAOnRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjEwLjMsIGh0dHBzOi8vbWF0cGxvdGxpYi5vcmcvZiW1igAAAAlwSFlzAAAPYQAAD2EBqD+naQAAHOFJREFUeJzt3X1QFPcdP/D3YeTkQY6ng+MmqBBUijGgSbXQ5BcfGJGqHduOtVoVTdNGIrGIppXGJ3zCh9YypgaTqYrjZNQmGZ00v+hUGa22Io4mWMcxHSEYjHKIWjhwhkNhf3/448LBHewee7e3e+/XzM3k9vb2Pnfh3n6+373d1QmCIICISEYBShdARNrDYCEi2TFYiEh2DBYikh2DhYhkx2AhItkxWIhIdgwWIpIdg4WIZMdgISLZqSpYzp07h1mzZsFsNkOn0+H48eMOjy9evBg6nc7hNn36dGWKJfJjqgqWR48eITU1FXv27HG5zvTp01FfX2+/HT582IsVEhEAPKN0AVJkZ2cjOzu7z3X0ej1MJpOXKiIiZ1QVLGKcPXsWMTExiIiIwJQpU7B582ZERUW5XN9ms8Fms9nvd3Z24uHDh4iKioJOp/NGyUSqIQgCWlpaYDabERDgesCjqWCZPn06fvrTnyIhIQE1NTX4wx/+gOzsbFRUVGDQoEFOn1NcXIyioiIvV0qkbrdv38azzz7r8nGdWs/HotPpcOzYMcyePdvlOl9//TWee+45nD59GlOnTnW6Ts+Opbm5GcOGDcPt27cRFhYmd9lEqma1WhEfH4+mpiYYDAaX62mqY+kpMTER0dHRqK6udhkser0eer2+1/KwsDAGC5EL/U0TqGqvkFTffvstHjx4gLi4OKVLIfIrqupYWltbUV1dbb9fW1uLqqoqREZGIjIyEkVFRfjZz34Gk8mEmpoa/O53v0NSUhKysrIUrJrI/6gqWC5fvozJkyfb7xcUFAAAcnJyUFpaiv/85z84ePAgmpqaYDabMW3aNGzatMnpUIeIPEe1k7eeYrVaYTAY0NzczDkWoh7Efj80PcdCRMpgsBCR7BgsRCQ7BgsRyY7BQkSyY7AQkewYLEQkOwYLEcmOwUJEsmOwEJHsGCxEJDsGCxHJjsFCRLJjsBCR7BgsRCQ7BgsRyY7BQkSyY7AQkexUdc5b0h65zozKq1b6FgYLeYWrACkc+rUs2y9uSXS6nIGjDAYLeUTPIJErQFxxtf2egcOg8Q4GC8nC20EiVs86GDTeoarJ23PnzmHWrFkwm83Q6XQ4fvy4w+OCIGDdunWIi4tDUFAQMjMzcfPmTWWK9QOCINhvhUO/drj5qp51dn8PJB9VBcujR4+QmpqKPXv2OH18x44d2L17N/bu3YvKykqEhIQgKysLbW1tXq5U25yFiVo5CxkaONVesEyn0+HYsWOYPXs2gKd/7GazGStXrsSqVasAAM3NzYiNjUVZWRl+8YtfiNouL1jmXPc/EzUHiRjdh0scKjkS+/3QzBxLbW0tLBYLMjMz7csMBgMmTpyIiooKl8Fis9lgs9ns961Wq8drVZOuQJE/TDz179nAg6D7e+0KGQaMNJoJFovFAgCIjY11WB4bG2t/zJni4mIUFRV5tDY1kjdQeofI//m/hTJst8erCMD5mcVOHnE/FLrePwNGGs0Ei7sKCwvtF5cHnnYs8fHxClakvK75kwFsweGeJ0LEGZ2u92s5Dxvp4dA9YBgu/dNMsJhMJgBAQ0MD4uLi7MsbGhqQlpbm8nl6vR56vd7T5anCwLqU78LEW0EiRs+w6R000kKicOjX7F5E0EywJCQkwGQyoby83B4kVqsVlZWVyM3NVbY4HydHoPhSmPSle9A4hoz4kODwqH+qCpbW1lZUV1fb79fW1qKqqgqRkZEYNmwY8vPzsXnzZowcORIJCQlYu3YtzGazfc8R9ebesMc3uxOpukLG3S6GwyPXVLW7+ezZs5g8eXKv5Tk5OSgrK4MgCFi/fj0++OADNDU14eWXX8Z7772HUaNGiX4Nf9nd7F6Xoq7uxB3udjH+0r2I/X6oKli8wR+CRXqXov1A6WkgAaPlcGGwuEnrwSItVPwvUHpyJ2C0HC4MFjdpNVikD30Evw6Unr4LGPHhAmhvaCT2+6GqY4XIPd2P6xGxNhgqvel0wCufFaLr8+lP9+OP/BE7lh601rFIHfowUPrnTveilc6FHQtJCBV2KVK4073427/f7Fh60ErHIiVUGCjuk9K9aKFzYcfixxgq3uPYvfTNnzoXBovGiAsVDn3kJGVo5C/hwqFQD2oeCokNFQaK54gdGql1WMShkJ9hqPgGsUMjrXcuDBYNYKj4FoYLg0X1GCq+yd/DhcGiYgwV3+bP4cJgUSmGijpICRctYbBoFkPFV4gNFy11LQwWFeq/W2Go+Box4aKlIRGDRWUYKur1Xbi4ppVwYbCoyMAvy0G+QfvzLQwWTWG34uv8Zb6FwaISHAJphz/MtzBYVIChoj1i51vUisFCpChtDokYLD6O3Yp2iR0SqZGmgmXDhg3Q6XQOt+TkZKXL8iCGitqJGRKpsWtR1SVWxRgzZgxOnz5tv//MM+p9i9y97E8EuDqHS9eF6NV0/hZNdSzA0yAxmUz2W3R0dJ/r22w2WK1Wh5sv4BDIf2hxSKS5YLl58ybMZjMSExPxy1/+EnV1dX2uX1xcDIPBYL/Fx8d7qdKBYKhojdaGRJoKlokTJ6KsrAwnT55EaWkpamtr8corr6ClpcXlcwoLC9Hc3Gy/3b5924sVO8chkD/TRtei6XPeNjU1Yfjw4di1axd+9atfiXqOL5zztu9gYbeiZU/PmbvN5eNKz7XwnLcAwsPDMWrUKFRXVytdimjsVkgLXYumg6W1tRU1NTWIi4tTuhSZsFvROjFzLWqgqWBZtWoV/vnPf+LWrVu4cOECfvKTn2DQoEGYN2+e0qURSeS6a1HD7IWmguXbb7/FvHnzMHr0aPz85z9HVFQULl68CKPRqHRponBuhYD+uxY1DIfU++sxJ44cOaJ0CUReIQiCT/9gTlMdi5px0pZ6U+8kLoNFFTgM8jdqn8RlsBD5NHVO4jJYfAAnbckZNU/iMliISHYMFiKf57tDHldE7W4uKCiQvOE1a9YgMjJS8vOoOw6D/F3XcMjV8UO+uttZVLCUlJQgPT0dgYGBojb6r3/9C3l5eQwWEbibmQaicOjX2Nb6nNJl9CL6B3LHjh1DTEyMqHWHDh3qdkFEpH6i5lgOHDgAg8EgeqPvv/8+YmNj3S6KiHpS1zyLqGDJycmBXq8XvdH58+cjJCTE7aII4PwKdelvt7Mv/p5lQHuFWltbffJ8sWrB+RWSgy/+DUkOltraWsyYMQMhISEwGAyIiIhAREQEwsPDERER4YkaiUhlJB/dvGDBAgiCgP379yM2NtYnd3URkbIkB8vVq1dx5coVjB492hP1EJFLrq895GskD4W+//3v+8SZ7LWNE7fkSG1HO0vuWP76179i6dKluHPnDp5//nkMHjzY4fEXXnhBtuKISBxf+wWu5GBpbGxETU0NlixZYl+m0+nsb6yjo0PWArWKe4RITr72C1zJwfLaa69h3LhxOHz4MCdvicgpycHyzTff4NNPP0VSUpIn6iEiDZA8eTtlyhRcvXrVE7UQkUZI7lhmzZqFFStW4Nq1axg7dmyvydsf//jHshVHROokOViWLl0KANi4cWOvx3xl8nbPnj3YuXMnLBYLUlNT8e6772LChAlKl0XkNyQPhTo7O13efCFUjh49ioKCAqxfvx5ffPEFUlNTkZWVhXv37ildGpHf0NypKXft2oVf//rXWLJkCVJSUrB3714EBwdj//79SpdG5DckB8vy5cuxe/fuXsv/8pe/ID8/X46a3Nbe3o4rV64gMzPTviwgIACZmZmoqKhw+hybzcYjtIlkJjlYPvnkE/zwhz/stTwjIwMff/yxLEW56/79++jo6Oh1kqnY2FhYLBanzykuLobBYLDf4uPjvVEqkaZJDpYHDx44PZtcWFgY7t+/L0tR3lRYWIjm5mb7jcdBEQ2c5GBJSkrCyZMney0/ceIEEhMTZSnKXdHR0Rg0aBAaGhocljc0NMBkMjl9jl6vR1hYmMONiAZG8u7mgoIC5OXlobGxEVOmTAEAlJeX409/+hNKSkrkrk+SwMBAvPjiiygvL8fs2bMBPN2LVV5ejry8PEVrI/Inbh0rZLPZsGXLFmzatAkAMGLECJSWlmLRokWyFyhVQUEBcnJy8NJLL2HChAkoKSnBo0ePHA6aJCLPkhwsAJCbm4vc3Fw0NjYiKCgIoaGhctfltrlz56KxsRHr1q2DxWJBWloaTp48yasGEHmRW8HSxWg0ylWHrPLy8jj0IVKQqMnb8ePH43//+5/ojb788su4c+eO20URkbqJ6liqqqpw9epV0ZdMraqqgs1mG1BhRKReoodCU6dOFX1hJJ78qX86nQ7FLYk8ixzJwpfOHgeIDJba2lrJG3722WclP4e+IwhPT6BMpEaigmX48OGeroMc6HB+ZjHP1E92ggCcn7lN6TJE09zRzUSkPAYLEcmOwUJEsmOwKKhrz5ArInfCkZ/ztT1CgBvBkpOTg3PnznmiFnLwdAKXSG0Tt4AbwdLc3IzMzEyMHDkSW7du5S9siagXycFy/Phx3LlzB7m5uTh69ChGjBiB7OxsfPzxx3j8+LEnaiQilXFrjsVoNKKgoABXr15FZWUlkpKSsHDhQpjNZqxYsQI3b96Uu07N4jwLDYQvzq8AA5y8ra+vx6lTp3Dq1CkMGjQIP/rRj3Dt2jWkpKTgz3/+s1w1+jHOs/g7Nc6vAG4Ey+PHj/HJJ59g5syZGD58OD766CPk5+fj7t27OHjwIE6fPo2//e1vTi9oRkT+QfL5WOLi4tDZ2Yl58+bh0qVLSEtL67XO5MmTER4eLkN5BPC4IVIfnSD2kOX/79ChQ5gzZw6GDBniqZoUZbVaYTAY0Nzc7NUTawuC0MeRzgKPG/JD/Q2DlJhfEfv9kNyxLFy4cECFkXvYtZCa8Je3qsBJXH+j1knbLgwWH9Hfbmei7nx1N3MXBouK8DctpBYMFh/Sd9fC4ZC/8MVJW6k0FSwjRoyATqdzuG3bpt5xqjPsWkgNNBUsALBx40bU19fbb2+99ZbSJcmIXYvWqX3StsuALljmi4YOHeryAvBqIObs/dz17L/UMAwCNNixbNu2DVFRURg3bhx27tyJJ0+e9Lm+zWaD1Wp1uPk2di1apZVuBdBYx7J8+XKMHz8ekZGRuHDhAgoLC1FfX49du3a5fE5xcTGKioq8WGX/2LWQM2rpVgA3ftLvbatXr8b27dv7XOfGjRtITk7utXz//v1444030NraCr1e7/S5NpvN4aqNVqsV8fHxXv9Jf099/8QfAAS88lkhw0UjxHQrvhAsYn/S7/PB0tjYiAcPHvS5TmJiIgIDA3stv379Op5//nl89dVXGD16tKjXU+pYIWfEhAuPIVI/tYQK4MFjhbzNaDTCaDS69dyqqioEBAQgJiZG5qq8g0MiAnwnVKTw+WARq6KiApWVlZg8eTKGDh2KiooKrFixAgsWLEBERITS5XnI04lcDonUS0sTtt1pZq+QXq/HkSNH8Oqrr2LMmDHYsmULVqxYgQ8++EDp0gak/2OIuJdIrdQ0BJJKMx3L+PHjcfHiRaXLUAyHRORLNNOxaJnYrsW3p+GpOy13KwCDRTUYLtqh9VABGCwaw/kWX6fVydqeGCwqIvZkUOxa1E3t3QrAYFEdDonUyx+GQF0YLCrEcFEffwoVgMGiYQwXX+Ev8yrdMVhUStx8C8NFaWJDRUvdCsBgUTWGi2/z11ABGCyqx3DxTf4cKgCDRRMYLr7F30MFYLBohtRwYcDIr+tz9fdQAVRwoidv86UTPbmj/5ND2dfk6RZkJGXPj5pDRez3gx2Lxoi/VCuHRnLxl1CRgsGiQe6ECwNGOilDH8B/QgXgUKgXtQ+FuhM/LAI4NJJG6o/etBIqHAqRhM4FYPcijtQuBdBOqEjBjqUHLXUsXaR1LgC7F+fc+Wm+1kKFHQvZ6XQ6bGt9jt2Lm9ztUrQWKlKwY+lBix1Ld+52L4D/nVO365vh711Kd5q5YJm3aT1YAHfCBfCngHE3UABthwrAYHGbPwRLl4EGDKCdkOn+LWCguMZgcZM/BQvgbrjYn636LmYg3UkXfwkVQIOTt1u2bEFGRgaCg4MRHh7udJ26ujrMmDEDwcHBiImJwdtvv40nT554t1CVkT6x6/BsnJ+5zWGiVw3/THWv9Wn97oWKv0/Q9kU1Fyxrb2/HnDlzkJ6ejn379vV6vKOjAzNmzIDJZMKFCxdQX1+PRYsWYfDgwdi6dasCFatLV8C418Houn05HYdKT7ctS4lu6xl2Az2bG8Okf6obCpWVlSE/Px9NTU0Oy0+cOIGZM2fi7t27iI2NBQDs3bsXv//979HY2IjAwECn27PZbLDZbPb7VqsV8fHxfjMUcmZgw6NeW+sVNIDnwsbZX7Ocp4X091AROxRSTcfSn4qKCowdO9YeKgCQlZWF3NxcXL9+HePGjXP6vOLiYhQVFXmrTFUYWPfSa2tOvtjOw0YOnjq3rL8HilSaCRaLxeIQKgDs9y0Wi8vnFRYWoqCgwH6/q2MhuQPGYcuqObk0A8U9igbL6tWrsX379j7XuXHjBpKTkz1Wg16vh16v99j2taArYAC5h0m+iWEycIoGy8qVK7F48eI+10lMFLe3wmQy4dKlSw7LGhoa7I+RPDzXxSiPgSIfRYPFaDTCaDTKsq309HRs2bIF9+7dQ0xMDADg1KlTCAsLQ0pKiiyvQd/RShfDMPEM1cyx1NXV4eHDh6irq0NHRweqqqoAAElJSQgNDcW0adOQkpKChQsXYseOHbBYLFizZg2WLVvGoY6HdQ8ZwLeDhkHiHarZ3bx48WIcPHiw1/IzZ85g0qRJAIBvvvkGubm5OHv2LEJCQpCTk4Nt27bhmWfE56e//fLWG5QMGgaJvPiTfjcxWLxHzsBhgHiH3/2OhdSn5xCKtEM1xwoRkXowWIhIdgwWIpIdg4WIZMdgISLZMViISHYMFiKSHYOFiGTHYCEi2TFYiEh2DBYikh2DhYhkx2AhItkxWIhIdgwWIpIdg4WIZMdgISLZMViISHYMFiKSHYOFiGTHYCEi2akmWLZs2YKMjAwEBwcjPDzc6To6na7X7ciRI94tlIjUc/mP9vZ2zJkzB+np6di3b5/L9Q4cOIDp06fb77sKISLyHNUES1FREQCgrKysz/XCw8N5EXgihakmWMRatmwZXn/9dSQmJmLp0qVYsmQJdDqdy/VtNhtsNpv9fnNzM4CnV3wjIkdd34v+LqCqqWDZuHEjpkyZguDgYPzjH//Am2++idbWVixfvtzlc4qLi+3dUHfx8fGeLJVI1VpaWmAwGFw+rui1m1evXo3t27f3uc6NGzeQnJxsv19WVob8/Hw0NTX1u/1169bhwIEDuH37tst1enYsnZ2dePjwIaKiovrsdDzJarUiPj4et2/f9vvrR/OzeMpXPgdBENDS0gKz2YyAANf7fhTtWFauXInFixf3uU5iYqLb2584cSI2bdoEm80GvV7vdB29Xt/rMV+Z8A0LC/PrL1N3/Cye8oXPoa9OpYuiwWI0GmE0Gj22/aqqKkRERLgMFSLyDNXMsdTV1eHhw4eoq6tDR0cHqqqqAABJSUkIDQ3F3//+dzQ0NOAHP/gBhgwZglOnTmHr1q1YtWqVsoUT+SNBJXJycgQAvW5nzpwRBEEQTpw4IaSlpQmhoaFCSEiIkJqaKuzdu1fo6OhQtnA3tLW1CevXrxfa2tqULkVx/CyeUtvnoOjkLRFpk2p+0k9E6sFgISLZMViISHYMFiKSHYPFx4g5PURdXR1mzJiB4OBgxMTE4O2338aTJ0+8W6gX7NmzByNGjMCQIUMwceJEXLp0SemSPO7cuXOYNWsWzGYzdDodjh8/7vC4IAhYt24d4uLiEBQUhMzMTNy8eVOZYvvAYPExXaeHyM3Ndfp4R0cHZsyYgfb2dly4cAEHDx5EWVkZ1q1b5+VKPevo0aMoKCjA+vXr8cUXXyA1NRVZWVm4d++e0qV51KNHj5Camoo9e/Y4fXzHjh3YvXs39u7di8rKSoSEhCArKwttbW1errQfCu/uJhcOHDggGAyGXss///xzISAgQLBYLPZlpaWlQlhYmGCz2bxYoWdNmDBBWLZsmf1+R0eHYDabheLiYgWr8i4AwrFjx+z3Ozs7BZPJJOzcudO+rKmpSdDr9cLhw4cVqNA1diwqU1FRgbFjxyI2Nta+LCsrC1arFdevX1ewMvm0t7fjypUryMzMtC8LCAhAZmYmKioqFKxMWbW1tbBYLA6fi8FgwMSJE33uc2GwqIzFYnEIFQD2+xaLRYmSZHf//n10dHQ4fZ9aeY/u6HrvavhcGCxesHr1aqfn4+1+++qrr5Quk0g2qjkIUc3kPD2EyWTqtXekoaHB/pgWREdHY9CgQfb31aWhoUEz79EdXe+9oaEBcXFx9uUNDQ1IS0tTqCrn2LF4gdFoRHJycp+3wMBAUdtKT0/HtWvXHPaOnDp1CmFhYUhJSfHUW/CqwMBAvPjiiygvL7cv6+zsRHl5OdLT0xWsTFkJCQkwmUwOn4vVakVlZaXPfS7sWHxMf6eHmDZtGlJSUrBw4ULs2LEDFosFa9aswbJlyzR13pmCggLk5OTgpZdewoQJE1BSUoJHjx5hyZIlSpfmUa2traiurrbfr62tRVVVFSIjIzFs2DDk5+dj8+bNGDlyJBISErB27VqYzWbMnj1buaKdUXq3FDnq7/QQgiAIt27dErKzs4WgoCAhOjpaWLlypfD48WPlivaQd999Vxg2bJgQGBgoTJgwQbh48aLSJXncmTNnnP7/z8nJEQTh6S7ntWvXCrGxsYJerxemTp0q/Pe//1W2aCd42gQikh3nWIhIdgwWIpIdg4WIZMdgISLZMViISHYMFiKSHYOFiGTHYCEi2TFYSDG3bt2yH93t6YPoysrK7K+Vn5/v0dciBgv5gNOnTzscWOcJc+fORX19vc8drKdVPAiRFBcVFYWoqCiPvkZQUBCCgoJEH0VOA8OOhWTR2NgIk8mErVu32pdduHABgYGBbnUj+/fvx5gxY6DX6xEXF4e8vDz7YzqdDu+//z5mzpyJ4OBgfO9730NFRQWqq6sxadIkhISEICMjAzU1NbK8N5KOwUKyMBqN2L9/PzZs2IDLly+jpaUFCxcuRF5eHqZOnSppW6WlpVi2bBl+85vf4Nq1a/j000+RlJTksM6mTZuwaNEiVFVVITk5GfPnz8cbb7yBwsJCXL58GYIgOIQReZnCR1eTxrz55pvCqFGjhPnz5wtjx44V2traXK5bW1srABC+/PJLh+Vms1l45513XD4PgLBmzRr7/YqKCgGAsG/fPvuyw4cPC0OGDOn13FdffVX47W9/K/4NkVvYsZCs/vjHP+LJkyf46KOP8OGHH0o++dS9e/dw9+7dfrucF154wf7fXSeXHjt2rMOytrY2WK1WSa9P8mCwkKxqampw9+5ddHZ24tatW5KfHxQUJGq9wYMH2/9bp9O5XNbZ2Sm5Bho4BgvJpr29HQsWLMDcuXOxadMmvP7665KvXDh06FCMGDHC47ufybO4u5lk884776C5uRm7d+9GaGgoPv/8c7z22mv47LPPJG1nw4YNWLp0KWJiYpCdnY2Wlhb8+9//xltvveWhyklu7FhIFmfPnkVJSQkOHTqEsLAwBAQE4NChQzh//jxKS0slbSsnJwclJSV47733MGbMGMycOdMnL3xOrvGct6SYW7duISEhAV9++aXXroszadIkpKWloaSkxCuv56/YsZDiMjIykJGR4dHX+PDDDxEaGorz58979HXoKXYspJgnT57Y9xzp9XrEx8d77LVaWlrsV1YMDw9HdHS0x16LGCxE5AEcChGR7BgsRCQ7BgsRyY7BQkSyY7AQkewYLEQkOwYLEcmOwUJEsvt/O0XBelCn8yIAAAAASUVORK5CYII=",
"text/plain": [
"<Figure size 258.065x259.74 with 1 Axes>"
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"model.plot(width=(30, 30))"
]
},
{
"cell_type": "code",
"execution_count": 20,
"id": "26e46f8c-0c60-4442-9c90-f2ba96bb8c89",
"metadata": {},
"outputs": [
{
"data": {
"text/plain": [
"10.113085130446276"
]
},
"execution_count": 20,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"in_sphere.r"
]
},
{
"cell_type": "code",
"execution_count": null,
"id": "43aad6d9-f5bf-462c-a24b-c9e4fd728281",
"metadata": {},
"outputs": [],
"source": []
}
],
"metadata": {
"kernelspec": {
"display_name": "Python 3 (ipykernel)",
"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.12.4"
}
},
"nbformat": 4,
"nbformat_minor": 5
}