Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_qn_solver_3d_polar_par.F90
Go to the documentation of this file.
1 
9 
11 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
12 #include "sll_assert.h"
13 #include "sll_working_precision.h"
14 
15  use sll_m_remapper, only: &
16  sll_t_layout_2d
17 
18  use sll_m_qn_solver_2d_polar_par, only: &
23 
24  implicit none
25 
26  public :: &
31 
32  private
33 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 
39 
40  type(sll_t_qn_solver_2d_polar_par), private :: solver_2d
41 
43 
44 
45 contains
46 
47 
49  subroutine sll_s_qn_solver_3d_polar_par_init(solver, &
50  layout_r, &
51  layout_a, &
52  rmin, &
53  rmax, &
54  nr, &
55  ntheta, &
56  bc_rmin, &
57  bc_rmax, &
58  rho_m0, &
59  b_magn, &
60  lambda, &
61  use_zonal_flow, &
62  epsilon_0, &
63  rgrid)
64 
65  type(sll_t_qn_solver_3d_polar_par), intent(inout) :: solver
66  type(sll_t_layout_2d), pointer :: layout_r
67  type(sll_t_layout_2d), pointer :: layout_a
68  real(f64), intent(in) :: rmin
69  real(f64), intent(in) :: rmax
70  integer(i32), intent(in) :: nr
71  integer(i32), intent(in) :: ntheta
72  integer(i32), intent(in) :: bc_rmin
73  integer(i32), intent(in) :: bc_rmax
74  real(f64), intent(in) :: rho_m0(:)
75  real(f64), intent(in) :: b_magn(:)
76  real(f64), optional, intent(in) :: lambda(:)
77  logical, optional, intent(in) :: use_zonal_flow
78  real(f64), optional, intent(in) :: epsilon_0
79  real(f64), target, optional, intent(in) :: rgrid(:)
80 
81  ! Initialize the 2D solver
82  call sll_s_qn_solver_2d_polar_par_init(solver%solver_2d, &
83  layout_r, &
84  layout_a, &
85  rmin, &
86  rmax, &
87  nr, &
88  ntheta, &
89  bc_rmin, &
90  bc_rmax, &
91  rho_m0, &
92  b_magn, &
93  lambda, &
94  use_zonal_flow, &
95  epsilon_0, &
96  rgrid)
97 
99 
100  !=============================================================================
102  subroutine sll_s_qn_solver_3d_polar_par_solve(solver, rhs, phi)
103  type(sll_t_qn_solver_3d_polar_par), intent(inout) :: solver
104  real(f64), intent(in) :: rhs(:, :, :)
105  real(f64), intent(out) :: phi(:, :, :)
106 
107  integer(i32) :: i3
108 
109  sll_assert(all(shape(rhs) == shape(phi)))
110 
111  do i3 = lbound(rhs, 3), ubound(rhs, 3)
112  call sll_s_qn_solver_2d_polar_par_solve(solver%solver_2d, rhs(:, :, i3), phi(:, :, i3))
113  end do
114 
116 
117  !=============================================================================
120  type(sll_t_qn_solver_3d_polar_par), intent(inout) :: solver
121 
122  call sll_s_qn_solver_2d_polar_par_free(solver%solver_2d)
123 
124  end subroutine sll_s_qn_solver_3d_polar_par_free
125 
Parallel quasi-neutrality solver on 2D polar mesh; uses FFT in theta and 2nd-order FD in r.
subroutine, public sll_s_qn_solver_2d_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 solver.
subroutine, public sll_s_qn_solver_2d_polar_par_solve(solver, rho, phi)
Solve the quasi-neutrality equation and get the electrostatic potential.
subroutine, public sll_s_qn_solver_2d_polar_par_free(solver)
Delete contents (local storage) of quasi-neutrality solver.
Parallel 3D quasi-neutrality solver on "extruded" 2D polar mesh.
subroutine, public sll_s_qn_solver_3d_polar_par_solve(solver, rhs, phi)
Solve the quasi-neutrality equation and get the electrostatic potential.
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.
subroutine, public sll_s_qn_solver_3d_polar_par_free(solver)
Delete contents (local storage) of quasi-neutrality solver.
Class for the 3D quasi-neutral solver in polar coordinates It is basically a wrapper around a 2D solv...
    Report Typos and Errors