Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Parallel 3D quasi-neutrality solver on "extruded" 2D polar mesh.
This module is a 3D wrapper around the 2D quasi-neutral solver in polar coordinates, with simple cycle over the 3rd direction.
Derived types and interfaces | |
type | sll_t_qn_solver_3d_polar_par |
Class for the 3D quasi-neutral solver in polar coordinates It is basically a wrapper around a 2D solver, which is reused within a cycle over the x3 coordinate. More... | |
Functions/Subroutines | |
subroutine, public | sll_s_qn_solver_3d_polar_par_init (solver, layout_r, layout_a, rmin, rmax, nr, ntheta, bc_rmin, bc_rmax, rho_m0, b_magn, lambda, use_zonal_flow, epsilon_0, rgrid) |
Initialize the 3D (wrapper) solver. More... | |
subroutine, public | sll_s_qn_solver_3d_polar_par_solve (solver, rhs, phi) |
Solve the quasi-neutrality equation and get the electrostatic potential. More... | |
subroutine, public | sll_s_qn_solver_3d_polar_par_free (solver) |
Delete contents (local storage) of quasi-neutrality solver. More... | |
subroutine, public sll_m_qn_solver_3d_polar_par::sll_s_qn_solver_3d_polar_par_free | ( | type(sll_t_qn_solver_3d_polar_par), intent(inout) | solver | ) |
Delete contents (local storage) of quasi-neutrality solver.
Definition at line 119 of file sll_m_qn_solver_3d_polar_par.F90.
subroutine, public sll_m_qn_solver_3d_polar_par::sll_s_qn_solver_3d_polar_par_init | ( | type(sll_t_qn_solver_3d_polar_par), intent(inout) | solver, |
type(sll_t_layout_2d), pointer | layout_r, | ||
type(sll_t_layout_2d), pointer | layout_a, | ||
real(f64), intent(in) | rmin, | ||
real(f64), intent(in) | rmax, | ||
integer(i32), intent(in) | nr, | ||
integer(i32), intent(in) | ntheta, | ||
integer(i32), intent(in) | bc_rmin, | ||
integer(i32), intent(in) | bc_rmax, | ||
real(f64), dimension(:), intent(in) | rho_m0, | ||
real(f64), dimension(:), intent(in) | b_magn, | ||
real(f64), dimension(:), intent(in), optional | lambda, | ||
logical, intent(in), optional | use_zonal_flow, | ||
real(f64), intent(in), optional | epsilon_0, | ||
real(f64), dimension(:), intent(in), optional, target | rgrid | ||
) |
Initialize the 3D (wrapper) solver.
[in,out] | solver | Poisson solver class |
layout_r | sequential in r direction | |
layout_a | sequential in theta direction | |
[in] | nr | number of cells radial |
[in] | ntheta | number of cells angular |
[in] | bc_rmin | boundary condition at r_min |
[in] | bc_rmax | boundary condition at r_max |
[in] | rho_m0 | radial profile: total mass density of equilibrium |
[in] | b_magn | radial profile: intensity of magnetic field |
[in] | lambda | radial profile: electron Debye length |
[in] | use_zonal_flow | if .false. set flux average to zero |
[in] | epsilon_0 | override default: vacuum permittivity |
[in] | rgrid | grid points along r |
Definition at line 49 of file sll_m_qn_solver_3d_polar_par.F90.
subroutine, public sll_m_qn_solver_3d_polar_par::sll_s_qn_solver_3d_polar_par_solve | ( | type(sll_t_qn_solver_3d_polar_par), intent(inout) | solver, |
real(f64), dimension(:, :, :), intent(in) | rhs, | ||
real(f64), dimension(:, :, :), intent(out) | phi | ||
) |
Solve the quasi-neutrality equation and get the electrostatic potential.
[in,out] | solver | 3D solver |
[in] | rhs | Charge density |
[out] | phi | Potential |
Definition at line 102 of file sll_m_qn_solver_3d_polar_par.F90.