{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# PWR Pin 2D Axial Symmetry - Solution" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import os\n", "import time\n", "\n", "import sys\n", "\n", "mesh_path = '../../../mesh/'\n", "benchmark_path = '../../../BenchmarkData/'\n", "sys.path.append('../../../models/')\n", "\n", "from ofelia.aux_obj import *\n", "\n", "import matplotlib.pyplot as plt\n", "from matplotlib import cm" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Setting OpenMC materials" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "The number of axial divisions is 8\n", "The pin power is set to 65.0 kW\n" ] } ], "source": [ "# from ofelia.constants import *\n", "from ofelia.materials import * # importing parameters and OpenMC materials\n", "\n", "path_res = './results_'+str(n_div)+'div_'+str(Power/1000)+'kW/'\n", "if not os.path.exists(path_res):\n", " os.system(\"mkdir \"+path_res)\n", "\n", "print('The number of axial divisions is '+ str(n_div))\n", "print('The pin power is set to ' +str(Power/1000)+' kW')\n", "\n", "xml_materials = updateXML(mat_dict, n_div)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Setting simulation parameters for FEniCSx" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "from ofelia.thermal import *\n", "from dolfinx.io import gmshio\n", "from mpi4py import MPI\n", "import gmsh\n", "from IPython.display import clear_output\n", "\n", "gdim = 2\n", "model_rank = 0\n", "mesh_comm = MPI.COMM_WORLD\n", "\n", "mesh_factor = 0.05\n", "\n", "# Initialize the gmsh module\n", "gmsh.initialize()\n", "\n", "# Load the .geo file\n", "gmsh.merge(mesh_path+'pin2D_pwr.geo')\n", "gmsh.model.geo.synchronize()\n", "\n", "gmsh.option.setNumber(\"Mesh.MeshSizeFactor\", mesh_factor)\n", "\n", "gmsh.model.mesh.generate(gdim)\n", "gmsh.model.mesh.optimize(\"Netgen\")\n", "\n", "# Domain\n", "domain, ct, ft = gmshio.model_to_mesh(gmsh.model, comm = mesh_comm, rank = model_rank, gdim = gdim )\n", "gmsh.finalize()\n", "\n", "clear_output()\n", "\n", "fdim = gdim - 1\n", "domain.topology.create_connectivity(fdim, gdim)\n", "\n", "region_markers = [1, 2, 3] # 1: fuel, 2: helium, 3: cladding\n", "\n", "# Boundary Markers\n", "water_out_marker = 10\n", "inl_marker = 20\n", "out_marker = 30\n", "\n", "# The thermal conducticity of the helium gap is increased to consider the radiation HT\n", "th_input = thermal_inputs(T_w_average, pressure)\n", "physical_param = {'th_cond': np.array([3, 0.25 * 3, 16.23]) / 100., # W/cm-K \n", " 'htc': th_input.compute_htc(pitch, clad_or, u_in)}\n", "\n", "# Defining and assembling the linear structures for the Thermal Calculations\n", "TH = thermal_solver(domain, ct, ft,\n", " physical_param, region_markers, water_out_marker)\n", "TH.assemble()\n", "\n", "Vnorms = norms(TH.funSpace, domain)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## RUN" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "##################################################################\n", "## ##\n", "## OOO FFFFF EEEEEE L III A ##\n", "## O O F E L I A A ##\n", "## O O F E L I AAAAA ##\n", "## O O FFFF EEEEEE L I A A ##\n", "## O O F E L I AAAAAAA ##\n", "## O O F E L I A A ##\n", "## OOO F EEEEEEE LLLLL III A A ##\n", "## ##\n", "##################################################################\n", "## ##\n", "## Openmc-FEnicsx for muLtiphysics tutorIAl ##\n", "## ------- ##\n", "## ERMETE Lab - Politecnico di Milano ##\n", "## ## \n", "## GitHub: https://github.com/ERMETE-Lab/MP-OFELIA ## \n", "## ##\n", "##################################################################\n", "\n", "\n", " -- Simulation Start --\n", "\n", "Tguess: fuel = 600.00 °C, water = 290.00 °C\n", "Tolerance for convergence: 5.00e-02\n", " \n", "Iter #00\n", " Running OpenMC\n", "\n", " Input Pin Power: 65.00 kW\n", " Pin Power calculated (MC): 67.0096 +- 1.75 kW\n", "\n", " keff: 1.380605 +- 0.001630\n", " alpha = 1.0000\n", " Next pop size = 1214 (neutrons/cycle)\n", "\n", " Running FEniCSx\n", " Temperature regional average error: 1.60e-01 | relative L2 error: 8.60e-01\n", " Average T_fuel: 741.78 °C | Average T_water: 307.00 °C\n", "\n", "Iter #01\n", " Running OpenMC\n", " keff: 1.360200 +- 0.001300\n", " alpha = 0.5357\n", " Next pop size = 1400 (neutrons/cycle)\n", "\n", " Running FEniCSx\n", " Temperature regional average error: 9.17e-02 | relative L2 error: 1.83e-01\n", " Power density: relative L2 error 3.985e-01\n", " Average T_fuel: 744.00 °C | Average T_water: 310.47 °C\n", "\n", "Iter #02\n", " Running OpenMC\n", " keff: 1.355326 +- 0.001248\n", " alpha = 0.4111\n", " Next pop size = 1824 (neutrons/cycle)\n", "\n", " Running FEniCSx\n", " Temperature regional average error: 1.27e-02 | relative L2 error: 2.56e-02\n", " Power density: relative L2 error 5.613e-02\n", " Average T_fuel: 743.73 °C | Average T_water: 309.98 °C\n", "\n", "Iter #03\n", " Running OpenMC\n", " keff: 1.354819 +- 0.001075\n", " alpha = 0.3352\n", " Next pop size = 2238 (neutrons/cycle)\n", "\n", " Running FEniCSx\n", " Temperature regional average error: 8.20e-03 | relative L2 error: 1.67e-02\n", " Power density: relative L2 error 3.630e-02\n", " Average T_fuel: 744.02 °C | Average T_water: 310.21 °C\n", "\n" ] } ], "source": [ "os.system(\"cat ../../../models/ofelia/header\")\n", "print(\"\\n -- Simulation Start --\\n\")\n", "\n", "# Picard settings initialization\n", "i = 0\n", "error = 100\n", "alpha = 1\n", "s1 = s1_val\n", "sum_s = s1\n", "\n", "# Guess temperature for the Picard Loop\n", "Tguess_w = Tin\n", "Tguess_f = 273.15 + 600\n", "\n", "Tguess = []\n", "Tguess.append(np.ones((n_div,)) * Tguess_f)\n", "Tguess.append(np.ones((n_div,)) * Tguess_w)\n", "\n", "# Defining class to convert power from OpenMC to FEniCSx\n", "ExtractPower = extract_power(n_div, Power, meshSize, pin_length, pin_r, J_to_eV, tallyDict)\n", "\n", "# The division in slices acts only on the active region\n", "slices = np.linspace(-l_active/2, l_active/2, n_div+1)\n", "\n", "# regLogic = parameters.assign(is_fuel + 2 * is_helium + 3 * is_clad + 4 * is_water + 5 * is_graph + 273.15)\n", "# TH.plot2D(regLogic, L, pitch, 'pictures/regions')\n", "\n", "#Define input folder path\n", "path = os.getcwd()\n", "\n", "phiE_list = [[], []] # mean, unc\n", "\n", "q3_relaxed_list = []\n", "q3_std_list = []\n", "z_list = []\n", "\n", "alpha_list = [1]\n", "pop_list = [s1]\n", "errors_q3_list = [1]\n", "\n", "T_old = Function(TH.funSpace)\n", "T_old.interpolate(lambda x: x[0] * 0. + Tin)\n", "\n", "## Set computational time\n", "comp_time = dict()\n", "start_time = time.time()\n", "start_cpu_time = time.process_time()\n", "\n", "# os.system('python ../../../models/ofelia/build_xml.py')\n", "from ofelia.build_xml import *\n", "RemoveFolders(path=os.getcwd())\n", "\n", "print('Tguess: fuel = {:.2f}'.format(np.mean(Tguess[0])-273.15)+' °C, water = {:.2f}'.format(np.mean(Tguess[1])-273.15)+' °C')\n", "print('Tolerance for convergence: {:.2e}'.format(tol))\n", "#print('Under relaxation constant: {:.2f}'.format(alpha))\n", "print(' ')\n", "\n", "while error > tol:\n", "\n", " print(f'Iter #{i+0:02}')\n", " input_folder = path + '/build_xml/it' + str(i)\n", " \n", " #Run the simulation\n", " print(' Running OpenMC')\n", " openmc.run(cwd=input_folder, output = False)\n", " \n", " #Rename the output with the info\n", " old_sp_name = 'statepoint.'+str(batches)+'.h5'\n", " old_sp_path = input_folder + '/' +old_sp_name\n", "\n", " new_sp_name = 'statepoint.'+str(batches)+'.it' + str(i) + '.h5'\n", " new_sp_path = input_folder + '/' +new_sp_name\n", "\n", " os.system('mv ' + old_sp_path + ' ' + new_sp_path)\n", "\n", " #Collecting the output (for simplicity, the k)\n", " sp=openmc.StatePoint(new_sp_path)\n", "\n", " # Get the energy flux from the iteration\n", " phiE_list, EnergyStairs = ExtractPower.getSpectrum(sp, phiE_list)\n", "\n", " # Get the axial power density from the iteration\n", " q3_unrelaxed, z, q3std, A = ExtractPower.eval(sp, i)\n", " q3_integral = np.trapz(q3_unrelaxed, z)\n", " uq3_integral = np.trapz(q3std, z)\n", "\n", " if i == 0:\n", " \n", " print('\\n Input Pin Power: {:.2f}'.format(Power/1000) + ' kW') \n", " print(' Pin Power calculated (MC): {:.4f}'.format(A*q3_integral/1000) + ' +- {:.2f}'.format(A*uq3_integral/1000) +' kW\\n')\n", "\n", " q3_relaxed_list.append(q3_unrelaxed)\n", " q3_relaxed = q3_unrelaxed\n", " q3_std_list.append(q3std)\n", " z_list.append(z)\n", "\n", " # Under Relaxation\n", " sn = (s1+ np.sqrt(s1*s1 + 4 * s1 * sum_s))/2\n", " sum_s=sn\n", " error = 1\n", "\n", " if i>0:\n", " sn = (s1+ np.sqrt(s1*s1 + 4 * s1 * sum_s))/2\n", " sum_s = sum_s + sn\n", " alpha = sn/sum_s\n", "\n", " alpha_list.append(alpha)\n", " pop_list.append(sn)\n", " \n", " q3_relaxed = q3_unrelaxed * alpha + q3_relaxed_list[i-1] * (1-alpha)\n", " \n", " errors_q3 = np.sqrt(np.trapz( (q3_relaxed - q3_relaxed_list[i-1])**2, z) / np.trapz( (q3_relaxed)**2, z) )\n", " errors_q3_list.append(errors_q3)\n", "\n", " q3_relaxed_list.append(q3_relaxed)\n", " q3_std_list.append(q3std)\n", "\n", " \n", " #Print the k value at the iteration\n", " k_eff = sp.keff\n", " keff = k_eff.n\n", " ukeff = k_eff.std_dev\n", " print(' keff: {:.6f}'.format(keff) + ' +- {:.6f}'.format(ukeff))\n", " print(' alpha = {:.4f}'.format(alpha)) \n", " print(' Next pop size = ' + str(round(sn))+' (neutrons/cycle)\\n') \n", "\n", " # np.savetxt(path + '/build_xml/it' + str(i) + '/q3_it' + str(i)+'.txt', np.vstack([q3, z]))\n", " # power_fun = importq3(path + '/build_xml/it' + str(i) + '/q3_it' + str(i)) # W/cm3\n", "\n", " print(' Running FEniCSx')\n", " power_fun, Tb_fun, res_bulk = th_input.mapping_q3_Tb(z, q3_relaxed, Tin, total_length, l_active, fuel_or)\n", " T_sol = TH.solve(power_fun, Tb_fun)\n", "\n", " # Computing regional average temperature\n", " average_T = []\n", " average_T.append(TH.computeSolidAverageT(region_markers[0], T_sol, slices))\n", "\n", " dimz_full = z.shape[0]\n", " if z.shape[0] < 100:\n", " dimz_full = 100\n", " else:\n", " dimz_full = z.shape[0]+100\n", " average_T.append(th_input.computeWaterAverageT(slices, dimz_full))\n", "\n", " # Computing regional errors\n", " errors_mat = np.zeros((len(average_T,)))\n", " for mat_ID in range(len(average_T)):\n", " errors_mat[mat_ID] = np.linalg.norm(average_T[mat_ID] - Tguess[mat_ID])/np.linalg.norm(Tguess[mat_ID])\n", "\n", " # Check convergence: maxIter vs error < tol?\n", " i += 1\n", " error = np.mean(errors_mat)\n", " error_L2 = Vnorms.L2norm(T_sol - T_old) / Vnorms.L2norm(T_old)\n", "\n", " Tguess = average_T.copy()\n", " T_old.x.array[:] = T_sol.x.array\n", " print(f' Temperature regional average error: {float(error):.2e} | relative L2 error: {float(error_L2):.2e}') \n", " if i > 1:\n", " print(f' Power density: relative L2 error {float(errors_q3):.3e}')\n", " error += errors_q3\n", " print(f' Average T_fuel: {float(np.mean(Tguess[0])-273.15):.2f} °C | Average T_water: {float(np.mean(Tguess[1])-273.15):.2f} °C\\n')\n", " \n", " error = errors_q3_list[i-1]\n", " if errors_q3_list[i-1] > tol:\n", " xml_materials.update(Tguess[0], Tguess[1], sn, i)\n", " \n", " if i >= maxIter:\n", " error = 0\n", " print(' ')\n", " print('###############################################################')\n", " print(' ')\n", " print('Warning: Maximum Number of iteration reached!')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Post-Processing and data storing" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Computational wall time = 2.9910 minutes\n", "Computational CPU time = 0.1665 minutes\n" ] } ], "source": [ "import pickle\n", "from dolfinx.io import XDMFFile\n", "\n", "plt.rcParams.update({\n", " \"text.usetex\": True,\n", " \"font.family\": \"serif\"\n", "})\n", "\n", "plt.rcParams['text.latex.preamble'] = r'\\usepackage{amssymb} \\usepackage{amsmath} \\usepackage{amsthm} \\usepackage{mathtools}'\n", "\n", "comp_time['Wall time'] = (time.time() - start_time) / 60\n", "comp_time['CPU time'] = (time.process_time() - start_cpu_time) / 60\n", "\n", "print('Computational wall time = {:.4f}'.format(comp_time['Wall time'])+' minutes')\n", "print('Computational CPU time = {:.4f}'.format(comp_time['CPU time'])+' minutes')\n", "\n", "# Temperature storage in .h5 format\n", "with XDMFFile(domain.comm, path_res+'T_sol.xdmf', \"w\") as loc:\n", " T_sol.name = \"T\"\n", " loc.write_mesh(domain)\n", " loc.write_function(T_sol)\n", "\n", "# Data store using numpy arrays\n", "res2d = TH.extract_2D_data(T_sol, total_length, clad_or)\n", "with open(path_res+'results.pin', 'wb') as f:\n", " pickle.dump([comp_time, res_bulk, Tguess, res2d], f)\n", " \n", "# Write axial power to txt\n", "zWrite = np.linspace(-total_length/2, total_length/2, int(1e3))\n", "q_to_write = np.zeros((len(zWrite), 2))\n", "q_to_write[:, 0] = zWrite[:].flatten()\n", "q_to_write[:, 1] = power_fun(zWrite, 0.).flatten()\n", "np.savetxt(path_res+'axial_powerDensity.txt', q_to_write, delimiter=',')\n", "\n", "# Write L2 error in taking the axial steps\n", "# step_errors = np.zeros((4,))\n", "# step_errors[0], step_errors[2] = TH.FuelAxError(T_sol, average_T, L_active, fuel_or, slices)\n", "# step_errors[1], step_errors[3] = TH.WaterAxError(T_sol, average_T, L_active, clad_or, pitch, slices)\n", "# np.savetxt('T_L2error_axsteps.txt', step_errors, delimiter=',')\n", "\n", "# Under relaxation plots\n", "fig = plt.figure (figsize = (6,4))\n", "plt.plot(range(len(pop_list)), pop_list, '--o', label=r'$S_n$')\n", "plt.ylabel('MC population (n/cycle)', fontsize = 20)\n", "plt.xlabel('Iteration $i$', fontsize = 20)\n", "plt.xticks(range(len(errors_q3_list)))\n", "plt.grid(which='major',linestyle='-')\n", "plt.grid(which='minor',linestyle='--')\n", "plt.tight_layout()\n", "#plt.legend(fontsize = 15)\n", "fig.savefig(path_res+'PopSizeIt.pdf', format='pdf', dpi=600, bbox_inches='tight')\n", "plt.close()\n", "\n", "\n", "fig = plt.figure (figsize = (6,4))\n", "plt.plot(range(len(alpha_list)), alpha_list, '--*', label=r'$\\alpha$')\n", "plt.ylabel('Under relaxation factor', fontsize = 20)\n", "plt.xticks(range(len(errors_q3_list)))\n", "plt.xlabel('Iteration $i$', fontsize = 20)\n", "plt.grid(which='major',linestyle='-')\n", "plt.grid(which='minor',linestyle='--')\n", "plt.tight_layout()\n", "#plt.legend(fontsize = 15)\n", "fig.savefig(path_res+'AlphaIt.pdf', format='pdf', dpi=600, bbox_inches='tight')\n", "plt.close()\n", "\n", "fig = plt.figure (figsize = (6,4))\n", "plt.semilogy(range(len(errors_q3_list)), errors_q3_list, '--x', label=r\"$ \\frac{\\lVert (q'''_i - q'''_{i-1}) \\rVert_{L^2}}{ \\lVert q'''_i \\rVert_{L^2}}$\")\n", "plt.axhline(y = tol, color = 'r', linestyle = '--', label=\"Threshold\")\n", "plt.xticks(range(len(errors_q3_list)))\n", "plt.xlabel('Iteration $i$', fontsize = 20)\n", "plt.ylabel(r'Relative $L^2$ Error', fontsize = 20)\n", "plt.grid(which='major',linestyle='-')\n", "plt.grid(which='minor',linestyle='--')\n", "plt.tight_layout()\n", "plt.legend(fontsize = 15)\n", "fig.savefig(path_res+'ErrorsIt.pdf', format='pdf', dpi=600, bbox_inches='tight')\n", "plt.close()\n", "\n", "\n", "# Plot flux spectrum\n", "fig,ax = plt.subplots(figsize = (10,6))\n", "iCol = np.linspace(0,1,len(phiE_list[0]))\n", "colors = [cm.jet(r) for r in iCol]\n", "\n", "for idx in range(len(phiE_list[0])):\n", " plt.stairs(phiE_list[0][idx], EnergyStairs, color = colors[idx], linewidth = 1.5, label = \"it \"+str(idx))\n", " plt.stairs(phiE_list[0][idx]+1*phiE_list[1][idx], EnergyStairs, baseline=phiE_list[0][idx]-1*phiE_list[1][idx], fill=True, alpha=0.2, facecolor = colors[idx])\n", "\n", "plt.xscale('log')\n", "plt.xlabel('Energy (eV)', fontsize = 20)\n", "plt.ylabel(r'Neutron fluence per unit lethargy $\\left(\\frac{\\text{n}}{\\text{cm}^2 \\, \\cdot \\, \\text{part}} \\right)$', fontsize = 16)\n", "plt.grid(linestyle = ':')\n", "plt.tight_layout()\n", "plt.legend(fontsize = 18, loc = 'upper left', framealpha = 1, ncol = 3)\n", "plt.tick_params(axis='both', which='both', labelsize=17)\n", "fig.savefig(path_res+'phiE.pdf', format='pdf', dpi=600, bbox_inches='tight')\n", "plt.close(fig)" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "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.10.12" } }, "nbformat": 4, "nbformat_minor": 2 }