SN
The SN struct describes the discrete-ordinates discretization method used for one particle. The alias Discrete_Ordinates is kept for backward compatibility.
Structure
Radiant.SN — Type
SNStructure used to define the discretization method associated with the transport of a particle.
Mandatory field(s)
particle::Particle: particle for which the discretization methods is definedsolver_type::String: type of solver for the transport calculations.quadrature_type::String: type of quadrature for the angular domain.quadrature_order::Int64: order of the quadrature for the angular domain.legendre_order::Int64: maximum order of the Legendre expansion for the differential cross-sections.scheme_type::Dict{String,String}: type of schemes for the spatial or energy discretization.scheme_order::Dict{String,Int64}: order of the expansion for the discretization schemes.
Optional field(s) - with default values
angular_fokker_planck::String = "finite-difference": type of discretization for the angular Fokker-Planck operation.angular_boltzmann::String = "galerkin-d": type of discretization for the Boltzmann operation.convergence_criterion::Float64 = 1e-7: convergence criterion of in-group iterations.maximum_iteration::Int64 = 300: maximum number of in-group iterations.acceleration::Int64 = "none": acceleration method for the in-group iterations.
Setters
Radiant.set_particle — Method
set_particle(this::SN,particle::Particle)To set the particle for which the transport discretization method is for.
Input Argument(s)
this::SN: discretization method.particle::Particle: particle.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_particle(electron) Radiant.set_solver_type — Method
set_solver_type(this::SN,solver_type::String)To set the solver for the particle transport.
Input Argument(s)
this::SN: discretization method.solver_type::String: solver type, which can takes the following values:solver_type = "BTE": Boltzmann transport equationsolver_type = "BFP": Boltzmann Fokker-Planck equationsolver_type = "BCSD": Boltzmann-CSD equationsolver_type = "FP": Fokker-Planck equationsolver_type = "CSD": Continuous slowing-down only equationsolver_type = "BFP-EF": Boltzmann Fokker-Planck without elastic
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_particle("BFP")Radiant.set_quadrature — Method
set_quadrature(this::SN,type::String,order::Int64,Qdims::Int64=0)To set the quadrature properties for the discretization of the angular domain.
Input Argument(s)
this::SN: discretization method.type::String: type of quadrature, which can takes the following values:type = "gauss-legendre": Gauss-Legendre quadrature (1D Cartesian geometry only).type = "gauss-lobatto": Gauss-Lobatto quadrature (1D Cartesian geometry only).type = "carlson": Carlson quadrature (2D or 3D Cartesian geometry only).type = "gauss-legendre-chebychev": product quadrature between Gauss-Legendre and Chebychev quadratures (2D or 3D Cartesian geometry only).type = "lebedev": Lebedev quadrature (2D or 3D Cartesian geometry only).
order::Int64: order of the quadrature, which is any integer greater than 2.Qdims::Int64: quadrature dimension, which can takes the following values:Qdims = 0: default value, Qdims = 1 with 1D quadrature, Qdims = 3 with quadrature over the unit sphere in 1D or in 2D geometry or in 3D geometry.Qdims = 1: 1D quadrature.Qdims = 2: quadrature over half the unit-sphere.Qdims = 3: quadrature over the unit sphere.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_quadrature("gauss-legendre",4)Radiant.set_legendre_order — Method
set_legendre_order(this::SN,legendre_order::Int64)To set the maximum order of the Legendre expansion of the differential cross-sections.
Input Argument(s)
this::SN: discretization method.legendre_order::Int64: maximum order of the Legendre expansion of the differential cross-sections.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_legendre_order(7)Radiant.set_angular_fokker_planck — Method
set_angular_fokker_planck(this::SN,angular_fokker_planck::String)To set the discretization method for the angular Fokker-Planck term.
Input Argument(s)
this::SN: discretization method.angular_fokker_planck::String: discretization method for the angular Fokker-Planck term, which can takes the following values:angular_fokker_planck = "finite-difference": finite difference discretization.angular_fokker_planck = "galerkin": galerkin moment-based discretization.angular_fokker_planck = "differential-quadrature": finite difference discretization (1D Cartesian geometry only).
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_angular_fokker_planck("differential-quadrature")Radiant.set_angular_boltzmann — Method
set_angular_boltzmann(this::SN,angular_boltzmann::String)To set the angular discretization method for the Boltzmann operator.
Input Argument(s)
this::SN: discretization method.angular_boltzmann::String: angular discretization method for the Boltzmann operator, which can takes the following values:angular_boltzmann = "standard": standard discrete ordinates (SN) method.angular_boltzmann = "galerkin-m": Galerkin method by inversion of the discrete-to-moment M matrix.angular_boltzmann = "galerkin-d": Galerkin method by inversion of the moment-to-discrete D matrix.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_angular_boltzmann("standard")Radiant.set_convergence_criterion — Method
set_convergence_criterion(this::SN,convergence_criterion::Float64)To set the convergence criterion for the in-group iterations.
Input Argument(s)
this::SN: discretization method.convergence_criterion::Float64: convergence criterion.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_convergence_criterion(1e-5)Radiant.set_maximum_iteration — Method
set_maximum_iteration(this::SN,maximum_iteration::Int64)To set the maximum number of in-group iterations.
Input Argument(s)
this::SN: discretization method.maximum_iteration::Int64: maximum number of in-group iterations.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_maximum_iteration(50)Radiant.set_scheme — Method
set_scheme(this::SN,axis::String,scheme_type::String,scheme_order::Int64)To set the type of discretization scheme for derivative along the specified spatial or energy axis.
Input Argument(s)
this::SN: discretization method.axis::String: variable of the derivative for which the scheme is applied, which can takes the following values:axis = "x": spatial x axis (discretization of the streaming term)axis = "y": spatial y axis (discretization of the streaming term)axis = "z": spatial z axis (discretization of the streaming term)axis = "E": spatial E axis (discretization of the continuous slowing-down term)
scheme_type::String: type of scheme to be applied, which can takes the following values:scheme_type = "DD": diamond difference scheme (any order)scheme_type = "DG": discontinuous Galerkin scheme (any order)scheme_type = "AWD": adaptive weighted scheme (1st and 2nd order only)
scheme_order::Int64: scheme order, which takes values greater than 1.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_scheme("x","DD",1)
julia> m.set_scheme("E","DG",2)Radiant.set_acceleration — Method
set_acceleration(this::SN,acceleration::String,parameter::Int64=0)To set the acceleration method for the in-group iteration process.
Input Argument(s)
this::SN: discretization method.acceleration::String: acceleration method, which takes the following valuesacceleration = "none": source iteration without accelerationacceleration = "livolant": Livolant two-point extrapolationacceleration = "anderson": depth-mAnderson accelerationacceleration = "gmres": matrix-free restarted GMRESacceleration = "bicgstab": matrix-free BiCGStab
parameter::Int64: optional tuning parameter for the chosen method. It is the Krylov subspace size before restart for"gmres"(default 30) and the memory depth for"anderson"(default 3); it is ignored by the other methods. A value of0keeps the current default.
Output Argument(s)
N/A
Examples
julia> m = SN()
julia> m.set_acceleration("gmres") # GMRES with the default restart of 30
julia> m.set_acceleration("gmres",50) # GMRES restarted every 50 Krylov vectors
julia> m.set_acceleration("anderson",5) # Anderson acceleration with memory depth 5Radiant.set_is_full_coupling — Method
set_is_full_coupling(this::SN,isFC::Bool)Set, for multidimensional high-order schemes, if the high-order moments are fully coupled or not. For example, with two linear schemes, the moments are either fully coupled [00,10,01,11] or not [00,10,01].
Input Argument(s)
this::SN: discretization method.isFC::Bool: boolean indicating if the high-orde moments are fully coupled or not.
Output Argument(s)
N/A
Getters
Radiant.get_particle — Method
get_particle(this::SN)Get the particle associated with the discretization methods.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
particle::Particle: particle.
Radiant.get_solver_type — Method
get_solver_type(this::SN)Get the type of solver for transport calculations.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
solver::String: type of solver for transport calculations.isCSD::Bool: indicate if continuous slowing-down term is used or not.
Radiant.get_quadrature_type — Method
get_quadrature_type(this::SN)Get the quadrature type.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
quadrature_type::String: quadrature type.
Radiant.get_quadrature_order — Method
get_quadrature_order(this::SN)Get the quadrature order.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
quadrature_order::Int64: quadrature order.
Radiant.get_quadrature_dimension — Method
get_quadrature_dimension(this::SN)Get the quadrature dimension.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
quadrature_dimension::Int64: quadrature dimension.
Radiant.get_legendre_order — Method
get_legendre_order(this::SN)Get the Legendre truncation order.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
legendre_order::Int64: Legendre truncation order.
Radiant.get_angular_boltzmann — Method
get_angular_boltzmann(this::SN)Get the type of angular discretization for the Boltzmann operator.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
angular_boltzmann::String: type of angular discretization for the Boltzmann operator.
Radiant.get_angular_fokker_planck — Method
get_angular_fokker_planck(this::SN)Get the type of angular discretization for the Fokker-Planck operator.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
angular_fokker_planck::String: type of angular discretization for the Fokker-Planck operator.
Radiant.get_convergence_criterion — Method
get_convergence_criterion(this::SN)Get the convergence criterion for in-group iteration convergence.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
convergence_criterion::Float64: convergence criterion.
Radiant.get_maximum_iteration — Method
get_maximum_iteration(this::SN)Get the maximum number of iterations for in-group iteration convergence.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
maximum_iteration::Int64: maximum number of iterations.
Radiant.get_acceleration — Method
get_acceleration(this::SN)Get the acceleration method for in-group iteration convergence.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
acceleration::String: acceleration method.
Radiant.get_gmres_restart — Method
get_gmres_restart(this::SN)Get the GMRES restart parameter for in-group iteration convergence.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
gmres_restart::Int64: Krylov subspace size before restart.
Radiant.get_anderson_depth — Method
get_anderson_depth(this::SN)Get the Anderson acceleration memory depth for in-group iteration convergence.
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
anderson_depth::Int64: Anderson memory depth.
Radiant.get_is_full_coupling — Method
get_is_full_coupling(this::SN)Get, for multidimensional high-order schemes, if the high-order moments are fully coupled or not. For example, with two linear schemes, the moments are either fully coupled [00,10,01,11] or not [00,10,01].
Input Argument(s)
this::SN: discretization method.
Output Argument(s)
isFC::Bool: boolean indicating if the high-orde moments are fully coupled or not.
Radiant.get_schemes — Method
get_schemes(this::SN,geometry::Geometry,isFC::Bool)Get the space and/or energy schemes informations.
Input Argument(s)
this::SN: discretization method.geometry::Geometry: geometry.isFC::Bool: boolean indicating if the high-order moments are fully coupled.
Output Argument(s)
schemes::Vector{String}: scheme types.𝒪::Vector{Int64}: order of the schemes.Nm::Vector{Int64}: numbers of moments.