DPN
The DPN solver is a moment-based ("double-Pn") angular discretization, used as an alternative to SN. It shares the same overall interface as SN but discretizes the angular variable on a polynomial basis instead of a quadrature.
Structure
Radiant.DPN — Type
DPNStructure used to define the moment-based ("double-Pn") angular discretization method associated with the transport of a particle.
Mandatory field(s)
particle::Particle: particle for which the discretization methods is defined.solver_type::String: type of solver for the transport calculations.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
legendre_order::Int64 = 64: maximum order of the Legendre expansion for the differential cross-sections.polynomial_basis::String: polynomial basis ("legendre"in 1D,"spherical-harmonics"in 2D or 3D by default).angular_fokker_planck::String = "galerkin": type of discretization for the angular Fokker-Planck operation.convergence_criterion::Float64 = 1e-7: convergence criterion of in-group iterations.maximum_iteration::Int64 = 300: maximum number of in-group iterations.acceleration::String = "none": acceleration method for the in-group iterations.
Setters
Radiant.set_particle — Method
set_particle(this::DPN,particle::Particle)To set the particle for which the transport discretization method is for.
Input Argument(s)
this::DPN: discretization method.particle::Particle: particle.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_particle(electron)Radiant.set_solver_type — Method
set_solver_type(this::DPN,solver_type::String)To set the solver for the particle transport.
Input Argument(s)
this::DPN: discretization method.solver_type::String: solver type, which can take the following values:solver_type = "BTE": Boltzmann transport equation.solver_type = "BFP": Boltzmann Fokker-Planck equation.solver_type = "BCSD": Boltzmann-CSD equation.solver_type = "FP": Fokker-Planck equation.solver_type = "CSD": Continuous slowing-down only equation.solver_type = "BFP-EF": Boltzmann Fokker-Planck without elastic.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_solver_type("BFP")Radiant.set_legendre_order — Method
set_legendre_order(this::DPN,legendre_order::Int64)To set the maximum order of the Legendre expansion of the differential cross-sections.
Input Argument(s)
this::DPN: discretization method.legendre_order::Int64: maximum order of the Legendre expansion of the differential cross-sections.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_legendre_order(7)Radiant.set_polynomial_basis — Method
set_polynomial_basis(this::DPN,basis::String)To set the polynomial basis used for the angular discretization.
Input Argument(s)
this::DPN: discretization method.basis::String: polynomial basis, which can take the following values:basis = "legendre": Legendre polynomials (1D only).basis = "spherical-harmonics": spherical harmonics (1D, 2D or 3D).
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_polynomial_basis("legendre")Radiant.set_angular_fokker_planck — Method
set_angular_fokker_planck(this::DPN,angular_fokker_planck::String)To set the discretization method for the angular Fokker-Planck term. With the DPN solver, only the "galerkin" discretization is currently supported.
Input Argument(s)
this::DPN: discretization method.angular_fokker_planck::String: discretization method for the angular Fokker-Planck term, which can take the following value:angular_fokker_planck = "galerkin": Galerkin moment-based discretization.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_angular_fokker_planck("galerkin")Radiant.set_scheme — Method
set_scheme(this::DPN,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::DPN: discretization method.axis::String: variable of the derivative for which the scheme is applied, which can take 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": energy axis (discretization of the continuous slowing-down term).
scheme_type::String: type of scheme to be applied, which can take 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 = DPN()
julia> m.set_scheme("x","DD",1)
julia> m.set_scheme("E","DG",2)Radiant.set_is_full_coupling — Method
set_is_full_coupling(this::DPN,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::DPN: discretization method.isFC::Bool: boolean indicating if the high-order moments are fully coupled or not.
Output Argument(s)
N/A
Radiant.set_acceleration — Method
set_acceleration(this::DPN,acceleration::String)To set the acceleration method for the in-group iteration process.
Input Argument(s)
this::DPN: discretization method.acceleration::String: acceleration method, which takes the following values:acceleration = "none": none.acceleration = "livolant": livolant acceleration method.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_acceleration("livolant")Radiant.set_convergence_criterion — Method
set_convergence_criterion(this::DPN,convergence_criterion::Float64)To set the convergence criterion for the in-group iterations.
Input Argument(s)
this::DPN: discretization method.convergence_criterion::Float64: convergence criterion.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_convergence_criterion(1e-5)Radiant.set_maximum_iteration — Method
set_maximum_iteration(this::DPN,maximum_iteration::Int64)To set the maximum number of in-group iterations.
Input Argument(s)
this::DPN: discretization method.maximum_iteration::Int64: maximum number of in-group iterations.
Output Argument(s)
N/A
Examples
julia> m = DPN()
julia> m.set_maximum_iteration(50)Getters
Radiant.get_particle — Method
get_particle(this::DPN)Get the particle associated with the discretization methods.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
particle::Particle: particle.
Radiant.get_solver_type — Method
get_solver_type(this::DPN)Get the type of solver for transport calculations.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
solver::Int64: type of solver for transport calculations.isCSD::Bool: indicate if continuous slowing-down term is used or not.
Radiant.get_legendre_order — Method
get_legendre_order(this::DPN)Get the Legendre truncation order.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
legendre_order::Int64: Legendre truncation order.
Radiant.get_polynomial_basis — Method
get_polynomial_basis(this::DPN,Ndims::Int64)Get the polynomial basis used for the angular discretization. If the basis was not set explicitly, the default basis is selected based on the geometry dimension ("legendre" in 1D, "spherical-harmonics" otherwise).
Input Argument(s)
this::DPN: discretization method.Ndims::Int64: geometry dimension.
Output Argument(s)
polynomial_basis::String: polynomial basis.
Radiant.get_angular_fokker_planck — Method
get_angular_fokker_planck(this::DPN)Get the type of angular discretization for the Fokker-Planck operator.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
angular_fokker_planck::String: type of angular discretization for the Fokker-Planck operator.
Radiant.get_schemes — Method
get_schemes(this::DPN,geometry::Geometry,isFC::Bool)Get the space and/or energy schemes informations.
Input Argument(s)
this::DPN: 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.
Radiant.get_is_full_coupling — Method
get_is_full_coupling(this::DPN)Get, for multidimensional high-order schemes, if the high-order moments are fully coupled or not.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
isFC::Bool: boolean indicating if the high-order moments are fully coupled or not.
Radiant.get_acceleration — Method
get_acceleration(this::DPN)Get the acceleration method for in-group iteration convergence.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
acceleration::String: acceleration method.
Radiant.get_convergence_criterion — Method
get_convergence_criterion(this::DPN)Get the convergence criterion for in-group iteration convergence.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
convergence_criterion::Float64: convergence criterion.
Radiant.get_maximum_iteration — Method
get_maximum_iteration(this::DPN)Get the maximum number of iterations for in-group iteration convergence.
Input Argument(s)
this::DPN: discretization method.
Output Argument(s)
maximum_iteration::Int64: maximum number of iterations.