fenicsx.neutronics package
Submodules
fenicsx.neutronics.neutr_diff module
- class fenicsx.neutronics.neutr_diff.steady_neutron_diff(domain: dolfinx.mesh.Mesh, ct: dolfinx.cpp.mesh.MeshTags_int32, ft: dolfinx.cpp.mesh.MeshTags_int32, physical_param: dict, regions_markers: list, void_bound_marker: int, albedo: ndarray | None = None, coupling: str | None = None)[source]
Bases:
object
” This class implements the solution of the steady multigroup neutron diffusion equation through the \(k\)-eigenvalue method.
- Parameters:
domain (dolfinx.mesh.Mesh) – Mesh imported.
ct (dolfinx.cpp.mesh.MeshTags_int32) – Cell tags of the regions in the domain.
ft (dolfinx.cpp.mesh.MeshTags_int32) – Face tags of the boundaries.
physical_param (dict) – Dictionary with the main physical parameters at reference condition.
regions_markers (dict) – List of the region markers inside the domain.
void_bound_marker (int) – Integer describing the index of the void boundary.
albedo (np.ndarray, optional (Default = None)) – If not None, it, denoted as \(\gamma\), must be an array of \(G\), as the energy groups, to impose the following BC \(-D_g\nabla \phi_g\cdot \mathbf{n} = \gamma \phi_g\).
coupling (str, optional (Default = None)) – If not None, it indicates the type of coupling with the temperature: log and sqrt are available.
- assembleForm(direct: bool = True)[source]
This function assembles the variational forms and the correspondent linear algebra structures.
- Parameters:
direct (bool, optional (Default=True)) – If True, a direct linear solver is used (Gauss Elimination), otherwise an iterative one is adopted (preconditioned GMRES with ILU).
- solve(temperature: dolfinx.fem.Function | None = None, power: float = 1, nu: float = 2.41, Ef: float = 3.2e-11, tol: float = 1e-10, maxIter: int = 200, LL: int = 50, verbose: bool = False)[source]
This function implements the \(k\)-eigenvalue method: starting from an initial non-trivial guess of eigenvalue-eigenvector pair \(\left(k_{\text{eff}}^{(0)},\,\boldsymbol{\phi}^{(0)}\right)\), at each iteration \(n\geq 1\) the following linear system is solved
\[\mathbb{A}\boldsymbol{\phi}^{(n)} = \frac{1}{k_{\text{eff}}^{(n-1)}} \mathbb{B}\boldsymbol{\phi}^{(n-1)}\]The convergence check is made on the relative difference for the eigenvalue as
\[\frac{|k_{\text{eff}}^{(n)} - k_{\text{eff}}^{(n-1)}|}{k_{\text{eff}}^{(n)}} < \epsilon\]- Parameters:
temperature (dolfinx.fem.Function, optional (Default = None)) – Temperature field as a dolfinx.fem.Function.
power (float, optional (Default = 1)) – Power level \(P\) used to normalise the flux as \(P = \displaystyle\int_V E_f\cdot \sum_{g=1}^G \Sigma_{f,g}\phi_g \, dV\).
nu (float, optional (Default = 2.41)) – Number of neutrons emitted by the fission event.
Ef (float, optional (Default = 1.6e-19 * 200e6)) – Energy released by the fission event.
tol (float, optional (Default = 1e-10)) – Tolerance of the inverse power method, \(\epsilon\).
maxIter (int, optional (Default = 200)) – Maximum iteration for the inverse power method.
LL (int, optional (Default = 50)) – Frequency for printing the output.
verbose (bool, optional (Default = False)) – If True, the output is printed every LL iterations.
- Returns:
normalised_fluxes (list[dolfinx.fem.Function]) – List of dolfinx.fem.Function with the solution of the inverse power method.
k_eff_ (float) – Converged eigenvalue.
- class fenicsx.neutronics.neutr_diff.transient_neutron_diff(domain: dolfinx.mesh.Mesh, ct: dolfinx.cpp.mesh.MeshTags_int32, ft: dolfinx.cpp.mesh.MeshTags_int32, physical_param: dict, regions_markers: list, void_bound_marker: int, albedo: ndarray | None = None, coupling: str | None = None)[source]
Bases:
object
” This class implements the solution of the transient multigroup neutron diffusion equation coupled with the precursors equations with a semi-implicit first-order advancement scheme.
- Parameters:
domain (dolfinx.mesh.Mesh) – Mesh imported.
ct (dolfinx.cpp.mesh.MeshTags_int32) – Cell tags of the regions in the domain.
ft (dolfinx.cpp.mesh.MeshTags_int32) – Face tags of the boundaries.
physical_param (dict) – Dictionary with the main physical parameters at reference condition.
regions_markers (dict) – List of the region markers inside the domain.
void_bound_marker (int) – Integer describing the index of the void boundary.
albedo (np.ndarray, optional (Default = None)) – If not None, it, denoted as \(\gamma\), must be an array of \(G\), as the energy groups, to impose the following BC \(-D_g\nabla \phi_g\cdot \mathbf{n} = \gamma \phi_g\).
coupling (str, optional (Default = None)) – If not None, it indicates the type of coupling with the temperature: log and sqrt are available.
- advance(t: float, sigma_ag_transient_value: list, temperature: dolfinx.fem.Function | None = None, return_prec: bool = False)[source]
This function implements a single step forward in time.
- Parameters:
t (float) – Time to advance to, say \(t_{n+1}\).
sigma_ag_transient_value (list) – List of callable function dependent on time used to impose variation of the reactivity.
temperature (dolfinx.fem.Function, optional (Default = None)) – Temperature field as a dolfinx.fem.Function.
return_prec (bool, optional (Default = False)) – If True, the precursors are returned in the solution.
- Returns:
power (float) – Power at time \(t_{n+1}\), defined as \(P(t) = \displaystyle\int_V E_f\cdot \sum_{g=1}^G \Sigma_{f,g}\phi_g \, dV\).
solution (list[dolfinx.fem.Function]) – List of dolfinx.fem.Function with the output of this time step.
- assembleForm(phi_ss: list[dolfinx.fem.Function], dt: float, nu: float = 2.41, Ef: float = 3.2e-11, direct=True)[source]
This function assembles the variational forms and the correspondent linear algebra structures.
- Parameters:
phi_ss (list[dolfinx.fem.Function) – List of the steady state group fluxes.
dt (float) – Time step size.
nu (float, optional (Default = 2.41)) – Number of neutrons emitted by the fission event.
Ef (float, optional (Default = 1.6e-19 * 200e6)) – Energy released by the fission event.
direct (bool, optional (Default=True)) – If True, a direct linear solver is used (Gauss Elimination), otherwise an iterative one is adopted (preconditioned GMRES with ILU).
fenicsx.neutronics.thermal module
- class fenicsx.neutronics.thermal.steady_thermal_diffusion(domain: dolfinx.mesh.Mesh, ct: dolfinx.cpp.mesh.MeshTags_int32, ft: dolfinx.cpp.mesh.MeshTags_int32, physical_param: dict, regions_markers: list, void_bound_marker: int, h: float | None = None, TD: float = 300.0, coupling=None)[source]
Bases:
object
This class implements the solution of the steady thermal diffusion equation with internal power generation made by neutrons.
- Parameters:
domain (dolfinx.mesh.Mesh) – Mesh imported.
ct (dolfinx.cpp.mesh.MeshTags_int32) – Cell tags of the regions in the domain.
ft (dolfinx.cpp.mesh.MeshTags_int32) – Face tags of the boundaries.
physical_param (dict) – Dictionary with the main physical parameters at reference condition.
regions_markers (dict) – List of the region markers inside the domain.
void_bound_marker (int) – Integer describing the index of the void boundary onto which BC is imposed (Dirichlet or Robin).
h (float, optional (Default = None)) – If not None, it, denoted as \(h\), represents the heat transfer coefficients in Newton’s cooling law.
TD (float, optional (Default = None)) – If not None, it, denoted as \(T_D\), represents either the bulk temperature in Newton’s cooling law or the Dirichlet condition to impose.
coupling (str, optional (Default = None)) – If not None, it indicates the type of coupling with the temperature: log and sqrt are available.
- assembleForm(direct=True)[source]
This function assembles the variational forms and the correspondent linear algebra structures.
- Parameters:
direct (bool, optional (Default=True)) – If True, a direct linear solver is used (Gauss Elimination), otherwise an iterative one is adopted (preconditioned GMRES with ILU).
- solve(phi: list[dolfinx.fem.Function], k_eff: float, temperature: dolfinx.fem.Function | None = None)[source]
This function solves the steady equation.
- Parameters:
phi (list[dolfinx.fem.Function]) – List of dolfinx.fem.Function for the internal power generation, \(q''' = E_f\cdot \sum_{g=1}^G \Sigma_{f,g}\phi_g\).
k_eff (float) – Effective multiplication factor.
temperature (dolfinx.fem.Function, optional (Default = None)) – Temperature field as a dolfinx.fem.Function used to calculate the cross sections.
- Returns:
solution – Output of the linear system.
- Return type:
dolfinx.fem.Function
- class fenicsx.neutronics.thermal.transient_thermal_diffusion(domain: dolfinx.mesh.Mesh, ct: dolfinx.cpp.mesh.MeshTags_int32, ft: dolfinx.cpp.mesh.MeshTags_int32, physical_param: dict, regions_markers: list, void_bound_marker: int, h: float | None = None, TD: float = 300.0, coupling=None)[source]
Bases:
object
This class implements the solution of the transient thermal diffusion equation with internal power generation made by neutrons.
- Parameters:
domain (dolfinx.mesh.Mesh) – Mesh imported.
ct (dolfinx.cpp.mesh.MeshTags_int32) – Cell tags of the regions in the domain.
ft (dolfinx.cpp.mesh.MeshTags_int32) – Face tags of the boundaries.
physical_param (dict) – Dictionary with the main physical parameters at reference condition.
regions_markers (dict) – List of the region markers inside the domain.
void_bound_marker (int) – Integer describing the index of the void boundary onto which BC is imposed (Dirichlet or Robin).
h (float, optional (Default = None)) – If not None, it, denoted as \(h\), represents the heat transfer coefficients in Newton’s cooling law.
TD (float, optional (Default = None)) – If not None, it, denoted as \(T_D\), represents either the bulk temperature in Newton’s cooling law or the Dirichlet condition to impose.
coupling (str, optional (Default = None)) – If not None, it indicates the type of coupling with the temperature: log and sqrt are available.
- advance(phi: list)[source]
This function advances in time with a single step.
- Parameters:
phi (list[dolfinx.fem.Function]) – List of dolfinx.fem.Function for the internal power generation, \(q''' = E_f\cdot \sum_{g=1}^G \Sigma_{f,g}\phi_g``\)
- Returns:
solution – Solution at the new time step..
- Return type:
dolfinx.fem.Function
- assembleForm(T_ss: dolfinx.fem.Function, dt: float, direct=True)[source]
This function assembles the variational forms and the correspondent linear algebra structures.
- Parameters:
T_ss (dolfinx.fem.Function) – Steady solution, i.e. the initial condition.
dt (float) – Time step size.
direct (bool, optional (Default=True)) – If True, a direct linear solver is used (Gauss Elimination), otherwise an iterative one is adopted (preconditioned GMRES with ILU).
Module contents
models/fenicsx/neutronics.
OFELIA: FEniCSx models for MultiPhysics.