Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_bsplines_base.F90
Go to the documentation of this file.
1 
5 
7 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_assert.h"
9 
10  use sll_m_working_precision, only: f64
11 
12  implicit none
13 
14  public :: &
16 
17  private
18 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
19 
21  integer, parameter :: wp = f64
22 
24  type, abstract :: sll_c_bsplines
25 
26  integer :: degree
27  logical :: periodic
28  logical :: uniform
29  logical :: radial
30  integer :: ncells
31  integer :: nbasis
32 
33  real(wp) :: xmin
34  real(wp) :: xmax
35 
36  real(wp), allocatable :: knots(:) ! Only used by non-uniform B-splines
37 
38  contains
39  procedure(i_fun_find_cell), deferred :: find_cell
40  procedure(i_sub_eval_basis), deferred :: eval_basis
41  procedure(i_sub_eval_deriv), deferred :: eval_deriv
42  procedure(i_sub_eval_basis_and_n_derivs), deferred :: eval_basis_and_n_derivs
43  procedure(i_sub_free), deferred :: free
44 
45  end type sll_c_bsplines
46 
47  abstract interface
48 
49  !---------------------------------------------------------------------------
54  !---------------------------------------------------------------------------
55  sll_pure function i_fun_find_cell(self, x) result(icell)
56  import sll_c_bsplines, wp
57  class(sll_c_bsplines), intent(in) :: self
58  real(wp), intent(in) :: x
59  integer :: icell
60  end function i_fun_find_cell
61 
62  !---------------------------------------------------------------------------
70  !---------------------------------------------------------------------------
71  sll_pure subroutine i_sub_eval_basis(self, x, values, jmin)
72  import sll_c_bsplines, wp
73  class(sll_c_bsplines), intent(in) :: self
74  real(wp), intent(in) :: x
75  real(wp), intent(out) :: values(:)
76  integer, intent(out) :: jmin
77  end subroutine i_sub_eval_basis
78 
79  !---------------------------------------------------------------------------
87  !---------------------------------------------------------------------------
88  sll_pure subroutine i_sub_eval_deriv(self, x, derivs, jmin)
89  import sll_c_bsplines, wp
90  class(sll_c_bsplines), intent(in) :: self
91  real(wp), intent(in) :: x
92  real(wp), intent(out) :: derivs(:)
93  integer, intent(out) :: jmin
94  end subroutine i_sub_eval_deriv
95 
96  !---------------------------------------------------------------------------
105  !---------------------------------------------------------------------------
106  sll_pure subroutine i_sub_eval_basis_and_n_derivs(self, x, n, derivs, jmin)
107  import sll_c_bsplines, wp
108  class(sll_c_bsplines), intent(in) :: self
109  real(wp), intent(in) :: x
110  integer, intent(in) :: n
111  real(wp), intent(out) :: derivs(:, :)
112  integer, intent(out) :: jmin
113  end subroutine i_sub_eval_basis_and_n_derivs
114 
115  !---------------------------------------------------------------------------
118  !---------------------------------------------------------------------------
119  subroutine i_sub_free(self)
120  import sll_c_bsplines
121  class(sll_c_bsplines), intent(inout) :: self
122  end subroutine i_sub_free
123 
124  end interface
125 
126 end module sll_m_bsplines_base
Find which grid cell contains the given point.
Abstract class for B-splines of arbitrary degree.
integer, parameter wp
Working precision.
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
    Report Typos and Errors