Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_spline_1d.F90
Go to the documentation of this file.
1 
20 
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 #include "sll_assert.h"
24 
25  use sll_m_working_precision, only: f64
27 
28  implicit none
29 
30  public :: &
32 
33  private
34 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
35 
37  integer, parameter :: wp = f64
38 
41 
42  real(wp), allocatable :: bcoef(:)
43  class(sll_c_bsplines), pointer, private :: bspl => null()
44 
45  contains
46 
47  procedure :: init => s_spline_1d__init
48  procedure :: free => s_spline_1d__free
49  procedure :: belongs_to_space => f_spline_1d__belongs_to_space
50  procedure :: eval => f_spline_1d__eval
51  procedure :: eval_deriv => f_spline_1d__eval_deriv
52  procedure :: eval_array => s_spline_1d__eval_array
53  procedure :: eval_array_deriv => s_spline_1d__eval_array_deriv
54 
55  end type sll_t_spline_1d
56 
57 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
58 contains
59 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 
61  !-----------------------------------------------------------------------------
65  !-----------------------------------------------------------------------------
66  subroutine s_spline_1d__init(self, bsplines)
67 
68  class(sll_t_spline_1d), intent(out) :: self
69  class(sll_c_bsplines), intent(in), target :: bsplines
70 
71  ! Store pointer to B-splines
72  self%bspl => bsplines
73 
74  ! Allocate array of spline coefficients: in case of periodic BCs, the last
75  ! p coefficients are a periodic copy of the first p ones
76  associate(n => bsplines%ncells, p => bsplines%degree)
77 
78  allocate (self%bcoef(1:n + p))
79 
80  end associate
81 
82  ! Set all coefficients to zero
83  self%bcoef = 0.0_f64
84 
85  end subroutine s_spline_1d__init
86 
87  !-----------------------------------------------------------------------------
90  !-----------------------------------------------------------------------------
91  subroutine s_spline_1d__free(self)
92 
93  class(sll_t_spline_1d), intent(inout) :: self
94 
95  deallocate (self%bcoef)
96  nullify (self%bspl)
97 
98  end subroutine s_spline_1d__free
99 
100  !-----------------------------------------------------------------------------
105  !-----------------------------------------------------------------------------
106  sll_pure function f_spline_1d__belongs_to_space(self, bsplines) result(in_space)
107 
108  class(sll_t_spline_1d), intent(in) :: self
109  class(sll_c_bsplines), intent(in), target :: bsplines
110  logical :: in_space
111 
112  in_space = associated(self%bspl, bsplines)
113 
114  end function f_spline_1d__belongs_to_space
115 
116  !-----------------------------------------------------------------------------
121  !-----------------------------------------------------------------------------
122  sll_pure function f_spline_1d__eval(self, x) result(y)
123 
124  class(sll_t_spline_1d), intent(in) :: self
125  real(wp), intent(in) :: x
126  real(wp) :: y
127 
128  real(wp) :: values(1:self%bspl%degree + 1)
129  integer :: jmin, jmax
130 
131  call self%bspl%eval_basis(x, values, jmin)
132 
133  jmax = jmin + self%bspl%degree
134 
135  y = dot_product(self%bcoef(jmin:jmax), values)
136 
137  end function f_spline_1d__eval
138 
139  !-----------------------------------------------------------------------------
144  !-----------------------------------------------------------------------------
145  sll_pure function f_spline_1d__eval_deriv(self, x) result(y)
146 
147  class(sll_t_spline_1d), intent(in) :: self
148  real(wp), intent(in) :: x
149  real(wp) :: y
150 
151  real(wp) :: derivs(1:self%bspl%degree + 1)
152  integer :: jmin, jmax
153 
154  call self%bspl%eval_deriv(x, derivs, jmin)
155 
156  jmax = jmin + self%bspl%degree
157 
158  y = dot_product(self%bcoef(jmin:jmax), derivs)
159 
160  end function f_spline_1d__eval_deriv
161 
162  !-----------------------------------------------------------------------------
167  !-----------------------------------------------------------------------------
168  sll_pure subroutine s_spline_1d__eval_array(self, x, y)
169 
170  class(sll_t_spline_1d), intent(in) :: self
171  real(wp), intent(in) :: x(:)
172  real(wp), intent(out) :: y(:)
173 
174  integer :: i
175 
176  sll_assert(size(x) == size(y))
177 
178  do i = 1, size(x)
179  y(i) = f_spline_1d__eval(self, x(i))
180  end do
181 
182  end subroutine s_spline_1d__eval_array
183 
184  !-----------------------------------------------------------------------------
189  !-----------------------------------------------------------------------------
190  sll_pure subroutine s_spline_1d__eval_array_deriv(self, x, y)
191 
192  class(sll_t_spline_1d), intent(in) :: self
193  real(wp), intent(in) :: x(:)
194  real(wp), intent(out) :: y(:)
195 
196  integer :: i
197 
198  sll_assert(size(x) == size(y))
199 
200  do i = 1, size(x)
201  y(i) = f_spline_1d__eval_deriv(self, x(i))
202  end do
203 
204  end subroutine s_spline_1d__eval_array_deriv
205 
206 end module sll_m_spline_1d
Abstract class for B-splines of arbitrary degree.
Module for 1D splines, linear combination of B-spline functions.
integer, parameter wp
Working precision.
subroutine s_spline_1d__free(self)
Destroy 1D spline (re-initialization is possible afterwards)
subroutine s_spline_1d__init(self, bsplines)
Initialize 1D spline object as element of span(B-splines)
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