Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_qn_2d_polar_splines_solver.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
19 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20 #include "sll_memory.h"
21 #include "sll_working_precision.h"
22 
23  use sll_m_qn_2d_base, only: &
25 
26  use sll_m_qn_2d_polar, only: &
32 
33  implicit none
34 
35  public :: &
39 
40  private
41 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
42 
44 
45  type(sll_t_qn_2d_polar) :: quasineutral
46 
47  contains
48  procedure, pass(qn) :: init => &
50  procedure, pass(qn) :: precompute_qn => &
52  procedure, pass(qn) :: solve_qn => &
54 
56 
57 contains
58 
60  eta_min, &
61  eta_max, &
62  Nc, &
63  N_points, &
64  lambda, &
65  T_i) &
66  result(qn)
67 
68  type(sll_t_qn_2d_polar_splines_solver), pointer :: qn
69  sll_real64, intent(in) :: eta_min(2)
70  sll_real64, intent(in) :: eta_max(2)
71  sll_int32, intent(in) :: nc(2)
72  sll_int32, intent(in) :: n_points
73  sll_real64, dimension(:), intent(in) :: lambda
74  sll_real64, dimension(:), intent(in) :: t_i
75  sll_int32 :: ierr
76 
77  sll_allocate(qn, ierr)
78  call qn%init( &
79  eta_min, &
80  eta_max, &
81  nc, &
82  n_points, &
83  lambda, &
84  t_i)
85 
87 
89  qn, &
90  eta_min, &
91  eta_max, &
92  Nc, &
93  N_points, &
94  lambda, &
95  T_i)
96 
97  class(sll_t_qn_2d_polar_splines_solver), intent(out) :: qn
98  sll_real64, intent(in) :: eta_min(2)
99  sll_real64, intent(in) :: eta_max(2)
100  sll_int32, intent(in) :: nc(2)
101  sll_int32, intent(in) :: n_points
102  sll_real64, dimension(:), intent(in) :: lambda
103  sll_real64, dimension(:), intent(in) :: t_i
104 
105  call sll_s_qn_2d_polar_init(qn%quasineutral, eta_min, eta_max, nc, n_points, lambda, t_i)
106 
108 
109  subroutine precompute_qn_2d_polar_splines(qn, mu_points, mu_weights, N_mu)
110  class(sll_t_qn_2d_polar_splines_solver), target :: qn
111  sll_int32, intent(in) :: n_mu
112  sll_real64, dimension(1:N_mu), intent(in) :: mu_points
113  sll_real64, dimension(1:N_mu), intent(in) :: mu_weights
114  sll_real64, dimension(1:N_mu) :: rho_points
115  sll_int32 :: i
116 
117  do i = 1, n_mu
118  rho_points(i) = sqrt(2._f64*mu_points(i))
119  end do
120 
121  call sll_s_precompute_gyroaverage_index(qn%quasineutral, rho_points, n_mu)
122  call sll_s_precompute_inverse_qn_matrix_polar_splines(qn%quasineutral, mu_points, mu_weights, n_mu)
123 
124  end subroutine precompute_qn_2d_polar_splines
125 
126  subroutine solve_qn_2d_polar_splines(qn, phi)
127  class(sll_t_qn_2d_polar_splines_solver), target :: qn
128  sll_real64, dimension(:, :), intent(inout) :: phi
129 
130  call sll_s_qn_2d_polar_solve(qn%quasineutral, phi)
131 
132  end subroutine solve_qn_2d_polar_splines
133 
type(sll_t_qn_2d_polar_splines_solver) function, pointer, public sll_f_new_qn_2d_polar_splines_solver(eta_min, eta_max, Nc, N_points, lambda, T_i)
subroutine, public sll_s_qn_2d_polar_splines_solver_init(qn, eta_min, eta_max, Nc, N_points, lambda, T_i)
subroutine precompute_qn_2d_polar_splines(qn, mu_points, mu_weights, N_mu)
subroutine, public sll_s_precompute_inverse_qn_matrix_polar_splines(quasineutral, mu_points, mu_weights, N_mu)
subroutine, public sll_s_qn_2d_polar_init(this, eta_min, eta_max, Nc, N_points, lambda, T_i)
subroutine, public sll_s_precompute_gyroaverage_index(quasineutral, rho, N_rho)
subroutine, public sll_s_qn_2d_polar_solve(quasineutral, phi)
    Report Typos and Errors