GN
The GN solver is a moment-based angular discretization that subdivides the angular domain and uses a polynomial expansion on each subdivision. It generalizes the DPN solver and shares most of its interface.
Structure
Radiant.GN — Type
GNStructure used to define the moment-based angular discretization method that subdivides the angular domain and uses a polynomial expansion on each subdivision. It generalizes the DPN solver.
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 = 16: global Legendre order.legendre_order_local::Int64 = 0: per-subdivision local Legendre order.subdivision::Int64 = 1: number of angular subdivisions.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::GN,particle::Particle)To set the particle for which the transport discretization method is for.
Input Argument(s)
this::GN: discretization method.particle::Particle: particle.
Output Argument(s)
N/A
Examples
julia> m = GN()
julia> m.set_particle(electron)Radiant.set_solver_type — Method
set_solver_type(this::GN,solver_type::String)To set the solver for the particle transport.
Input Argument(s)
this::GN: 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 = GN()
julia> m.set_solver_type("BFP")Radiant.set_legendre_order — Method
set_legendre_order(this::GN,legendre_order::Int64,legendre_order_local::Int64)To set the global Legendre order and the per-subdivision local Legendre order.
Input Argument(s)
this::GN: discretization method.legendre_order::Int64: global Legendre order.legendre_order_local::Int64: per-subdivision local Legendre order.
Output Argument(s)
N/A
Examples
julia> m = GN()
julia> m.set_legendre_order(16,0)Radiant.set_subdivision — Method
set_subdivision(this::GN,subdivision::Int64)To set the number of angular subdivisions used by the GN solver.
Input Argument(s)
this::GN: discretization method.subdivision::Int64: number of angular subdivisions (must be at least 1).
Output Argument(s)
N/A
Examples
julia> m = GN()
julia> m.set_subdivision(4)Radiant.set_polynomial_basis — Method
set_polynomial_basis(this::GN,basis::String)To set the polynomial basis used for the angular discretization.
Input Argument(s)
this::GN: 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 = GN()
julia> m.set_polynomial_basis("legendre")Radiant.set_angular_fokker_planck — Method
set_angular_fokker_planck(this::GN,angular_fokker_planck::String)To set the discretization method for the angular Fokker-Planck term.
Input Argument(s)
this::GN: discretization method.angular_fokker_planck::String: discretization method for the angular Fokker-Planck term, which can take the following values:angular_fokker_planck = "galerkin": Galerkin moment-based discretization (diagonal on full-range spherical harmonics,ℳ[p,p] = -ℓ_p(ℓ_p+1)).angular_fokker_planck = "finite-difference": finite-volume (TPFA) discretization of the Laplace-Beltrami operator on the GN patch graph, using the natural connectivity of the tiling (no Voronoi needed).
Output Argument(s)
N/A
Examples
julia> m = GN()
julia> m.set_angular_fokker_planck("finite-difference")Radiant.set_scheme — Method
set_scheme(this::GN,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::GN: 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 = GN()
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::GN,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::GN: 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::GN,acceleration::String,parameter::Int64=0)To set the acceleration method for the in-group iteration process.
Input Argument(s)
this::GN: discretization method.acceleration::String: acceleration method, which takes the following values:acceleration = "none": source iteration without acceleration.acceleration = "livolant": Livolant two-point extrapolation.acceleration = "anderson": depth-mAnderson acceleration.acceleration = "gmres": matrix-free restarted GMRES.acceleration = "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 = GN()
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_convergence_criterion — Method
set_convergence_criterion(this::GN,convergence_criterion::Float64)To set the convergence criterion for the in-group iterations.
Input Argument(s)
this::GN: discretization method.convergence_criterion::Float64: convergence criterion.
Output Argument(s)
N/A
Examples
julia> m = GN()
julia> m.set_convergence_criterion(1e-5)Radiant.set_maximum_iteration — Method
set_maximum_iteration(this::GN,maximum_iteration::Int64)To set the maximum number of in-group iterations.
Input Argument(s)
this::GN: discretization method.maximum_iteration::Int64: maximum number of in-group iterations.
Output Argument(s)
N/A
Examples
julia> m = GN()
julia> m.set_maximum_iteration(50)Getters
Radiant.get_particle — Method
get_particle(this::GN)Get the particle associated with the discretization methods.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
particle::Particle: particle.
Radiant.get_solver_type — Method
get_solver_type(this::GN)Get the type of solver for transport calculations.
Input Argument(s)
this::GN: 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::GN)Get the global Legendre truncation order.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
legendre_order::Int64: Legendre truncation order.
Radiant.get_legendre_order_local — Method
get_legendre_order_local(this::GN)Get the per-subdivision local Legendre truncation order.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
legendre_order_local::Int64: per-subdivision local Legendre truncation order.
Radiant.get_subdivision — Method
get_subdivision(this::GN)Get the number of angular subdivisions.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
subdivision::Int64: number of angular subdivisions.
Radiant.get_polynomial_basis — Method
get_polynomial_basis(this::GN,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::GN: 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::GN)Get the type of angular discretization for the Fokker-Planck operator.
Input Argument(s)
this::GN: 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::GN,geometry::Geometry,isFC::Bool)Get the space and/or energy schemes informations.
Input Argument(s)
this::GN: 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::GN)Get, for multidimensional high-order schemes, if the high-order moments are fully coupled or not.
Input Argument(s)
this::GN: 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::GN)Get the acceleration method for in-group iteration convergence.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
acceleration::String: acceleration method.
Radiant.get_gmres_restart — Method
get_gmres_restart(this::GN)Get the GMRES restart parameter for in-group iteration convergence.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
gmres_restart::Int64: Krylov subspace size before restart.
Radiant.get_anderson_depth — Method
get_anderson_depth(this::GN)Get the Anderson acceleration memory depth for in-group iteration convergence.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
anderson_depth::Int64: Anderson memory depth.
Radiant.get_convergence_criterion — Method
get_convergence_criterion(this::GN)Get the convergence criterion for in-group iteration convergence.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
convergence_criterion::Float64: convergence criterion.
Radiant.get_maximum_iteration — Method
get_maximum_iteration(this::GN)Get the maximum number of iterations for in-group iteration convergence.
Input Argument(s)
this::GN: discretization method.
Output Argument(s)
maximum_iteration::Int64: maximum number of iterations.