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_hermite_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: &
34 
35  implicit none
36 
37  public :: &
39 
40  private
41 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
42 
44 
45  type(sll_t_plan_gyroaverage_polar), pointer :: gyro
46  sll_int32 :: hermite_case
47  ! hermite_case
48  ! 1 : hermite standard
49  ! 2 : hermite C^1
50  ! 3 : hermite C^1 with precomputation
51  ! 4 : hermite C^1 with invariance
52 
53  contains
54  procedure, pass(gyroaverage) :: initialize => &
56  procedure, pass(gyroaverage) :: compute_gyroaverage => &
58 
60 
61 contains
63  eta_min, &
64  eta_max, &
65  Nc, &
66  N_points, &
67  interp_degree, &
68  larmor_rad, &
69  hermite_case) &
70  result(gyroaverage)
71 
72  type(gyroaverage_2d_polar_hermite_solver), pointer :: gyroaverage
73  sll_real64, intent(in) :: eta_min(2)
74  sll_real64, intent(in) :: eta_max(2)
75  sll_int32, intent(in) :: nc(2)
76  sll_int32, intent(in) :: n_points
77  sll_int32, intent(in) :: interp_degree(2)
78  sll_real64, intent(in) :: larmor_rad
79  sll_int32, optional :: hermite_case
80  sll_int32 :: ierr
81 
82  sll_allocate(gyroaverage, ierr)
84  gyroaverage, &
85  eta_min, &
86  eta_max, &
87  nc, &
88  n_points, &
89  interp_degree, &
90  larmor_rad, &
91  hermite_case)
92 
94 
96  gyroaverage, &
97  eta_min, &
98  eta_max, &
99  Nc, &
100  N_points, &
101  interp_degree, &
102  larmor_rad, &
103  hermite_case)
104  class(gyroaverage_2d_polar_hermite_solver) :: gyroaverage
105  sll_real64, intent(in) :: eta_min(2)
106  sll_real64, intent(in) :: eta_max(2)
107  sll_int32, intent(in) :: nc(2)
108  sll_int32, intent(in) :: n_points
109  sll_int32, intent(in) :: interp_degree(2)
110  sll_real64, intent(in) :: larmor_rad
111  sll_int32 :: deriv_size
112  sll_int32, optional :: hermite_case
113 
114  if (.not. (present(hermite_case))) then
115  gyroaverage%hermite_case = 2
116  else
117  gyroaverage%hermite_case = hermite_case
118  end if
119 
120  deriv_size = 9
121  select case (gyroaverage%hermite_case)
122  case (2, 3, 4)
123  deriv_size = 4
124  end select
125 
126  select case (gyroaverage%hermite_case)
127  case (1, 2, 3, 4)
128  gyroaverage%gyro => sll_f_new_plan_gyroaverage_polar_hermite( &
129  eta_min, &
130  eta_max, &
131  nc, &
132  n_points, &
133  interp_degree, &
134  deriv_size)
135  case default
136  print *, '#bad value of hermite_case=', gyroaverage%hermite_case
137  print *, '#not implemented'
138  print *, '#in initialize_gyroaverage_2d_polar_hermite_solver'
139  stop
140  end select
141 
142  select case (gyroaverage%hermite_case)
143  case (3)
144  call sll_s_pre_compute_gyroaverage_polar_hermite_c1(gyroaverage%gyro, larmor_rad)
145  end select
146 
148 
149  subroutine compute_gyroaverage_2d_polar_hermite(gyroaverage, larmor_rad, f)
150  class(gyroaverage_2d_polar_hermite_solver), target :: gyroaverage
151  sll_real64, intent(in) :: larmor_rad
152  sll_real64, dimension(:, :), intent(inout) :: f
153 
154  select case (gyroaverage%hermite_case)
155  case (1)
156  call sll_s_compute_gyroaverage_points_polar_hermite(gyroaverage%gyro, f, larmor_rad)
157  case (2)
158  call sll_s_compute_gyroaverage_points_polar_hermite_c1(gyroaverage%gyro, f, larmor_rad)
159  case (3)
161  case (4)
162  call sll_s_compute_gyroaverage_points_polar_with_invar_hermite_c1(gyroaverage%gyro, f, larmor_rad)
163  case default
164  print *, '#bad value of hermite_case=', gyroaverage%hermite_case
165  print *, '#not implemented'
166  print *, 'compute_gyroaverage_2d_polar_hermite'
167  stop
168  end select
169 
171 
type(gyroaverage_2d_polar_hermite_solver) function, pointer, public sll_f_new_gyroaverage_2d_polar_hermite_solver(eta_min, eta_max, Nc, N_points, interp_degree, larmor_rad, hermite_case)
subroutine compute_gyroaverage_2d_polar_hermite(gyroaverage, larmor_rad, f)
subroutine initialize_gyroaverage_2d_polar_hermite_solver(gyroaverage, eta_min, eta_max, Nc, N_points, interp_degree, larmor_rad, hermite_case)
subroutine, public sll_s_compute_gyroaverage_points_polar_hermite(gyro, f, rho)
subroutine, public sll_s_compute_gyroaverage_points_polar_hermite_c1(gyro, f, rho)
type(sll_t_plan_gyroaverage_polar) function, pointer, public sll_f_new_plan_gyroaverage_polar_hermite(eta_min, eta_max, Nc, N_points, interp_degree, deriv_size)
subroutine, public sll_s_compute_gyroaverage_pre_compute_polar_hermite_c1(gyro, f)
subroutine, public sll_s_pre_compute_gyroaverage_polar_hermite_c1(gyro, rho)
subroutine, public sll_s_compute_gyroaverage_points_polar_with_invar_hermite_c1(gyro, f, rho)
    Report Typos and Errors