Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_gyroaverage_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_gyroaverage_2d_base, only: &
25 
26  use sll_m_gyroaverage_2d_polar, only: &
35 
36  implicit none
37 
38  public :: &
40 
41  private
42 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
43 
45 
46  type(sll_t_plan_gyroaverage_polar), pointer :: gyro
47  sll_int32 :: splines_case
48  ! splines_case
49  ! 1 : splines
50  ! 2 : splines with invariance
51  ! 3 : splines pre-compute
52  ! 4 : splines pre-compute with FFT
53 
54  contains
55  procedure, pass(gyroaverage) :: initialize => &
57  procedure, pass(gyroaverage) :: compute_gyroaverage => &
59 
61 
62 contains
64  eta_min, &
65  eta_max, &
66  Nc, &
67  N_points, &
68  splines_case) &
69  result(gyroaverage)
70 
71  type(gyroaverage_2d_polar_splines_solver), pointer :: gyroaverage
72  sll_real64, intent(in) :: eta_min(2)
73  sll_real64, intent(in) :: eta_max(2)
74  sll_int32, intent(in) :: nc(2)
75  sll_int32 :: n_points
76  sll_int32, optional :: splines_case
77  sll_int32 :: ierr
78 
79  sll_allocate(gyroaverage, ierr)
81  gyroaverage, &
82  eta_min, &
83  eta_max, &
84  nc, &
85  n_points, &
86  splines_case)
87 
89 
91  gyroaverage, &
92  eta_min, &
93  eta_max, &
94  Nc, &
95  N_points, &
96  splines_case)
97  class(gyroaverage_2d_polar_splines_solver), intent(out) :: gyroaverage
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_int32, optional :: splines_case
103 
104  if (.not. (present(splines_case))) then
105  gyroaverage%splines_case = 1
106  else
107  gyroaverage%splines_case = splines_case
108  end if
109 
110  select case (gyroaverage%splines_case)
111  case (1, 2, 3, 4)
112  gyroaverage%gyro => sll_f_new_plan_gyroaverage_polar_splines( &
113  eta_min, &
114  eta_max, &
115  nc, &
116  n_points)
117  case default
118  print *, '#bad value of splines_case=', gyroaverage%splines_case
119  print *, '#not implemented'
120  print *, '#in initialize_gyroaverage_2d_polar_splines_solver'
121  stop
122  end select
123 
125 
126  subroutine compute_gyroaverage_2d_polar_splines(gyroaverage, larmor_rad, f)
127  class(gyroaverage_2d_polar_splines_solver), target :: gyroaverage
128  sll_real64, intent(in) :: larmor_rad
129  sll_real64, dimension(:, :), intent(inout) :: f
130 
131  select case (gyroaverage%splines_case)
132  case (1)
133  call sll_s_compute_gyroaverage_points_polar_spl(gyroaverage%gyro, f, larmor_rad)
134  case (2)
135  call sll_s_compute_gyroaverage_points_polar_with_invar_spl(gyroaverage%gyro, f, larmor_rad)
136  case (3)
137  call sll_s_pre_compute_gyroaverage_polar_spl(gyroaverage%gyro, larmor_rad)
138  call sll_s_compute_gyroaverage_pre_compute_polar_spl(gyroaverage%gyro, f)
139  case (4)
140  call sll_s_pre_compute_gyroaverage_polar_spl_fft(gyroaverage%gyro, larmor_rad)
142  case default
143  print *, '#bad value of splines_case=', gyroaverage%splines_case
144  print *, '#not implemented'
145  print *, 'compute_gyroaverage_2d_polar_splines'
146  stop
147  end select
148 
150 
subroutine compute_gyroaverage_2d_polar_splines(gyroaverage, larmor_rad, f)
type(gyroaverage_2d_polar_splines_solver) function, pointer, public sll_f_new_gyroaverage_2d_polar_splines_solver(eta_min, eta_max, Nc, N_points, splines_case)
subroutine initialize_gyroaverage_2d_polar_splines_solver(gyroaverage, eta_min, eta_max, Nc, N_points, splines_case)
type(sll_t_plan_gyroaverage_polar) function, pointer, public sll_f_new_plan_gyroaverage_polar_splines(eta_min, eta_max, Nc, N_points)
subroutine, public sll_s_compute_gyroaverage_points_polar_spl(gyro, f, rho)
subroutine, public sll_s_pre_compute_gyroaverage_polar_spl(gyro, rho)
subroutine, public sll_s_pre_compute_gyroaverage_polar_spl_fft(gyro, rho)
subroutine, public sll_s_compute_gyroaverage_pre_compute_polar_spl(gyro, f)
subroutine, public sll_s_compute_gyroaverage_points_polar_with_invar_spl(gyro, f, rho)
subroutine, public sll_s_compute_gyroaverage_pre_compute_polar_spl_fft(gyro, f)
    Report Typos and Errors