Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_finite_difference_interpolator_1d.F90
Go to the documentation of this file.
1 
18 #include "sll_working_precision.h"
19 #include "sll_memory.h"
20 #include "sll_assert.h"
22  implicit none
23  private
24 
27  sll_int32 :: num_points
28  sll_real64 :: delta ! cell size, distance between data points
29  sll_real64 :: r_delta ! reciprocal of the cell size
30  contains
31  procedure, pass(interpolator) :: initialize => &
33  procedure :: compute_interpolants => null_fd_1d_array_msg
34  procedure :: interpolate_from_interpolant_value => null_fd_1d_arg_msg
35  ! FIXME: define is we want this to explicitly state that we are taking
36  ! derivatives with respect to eta or whatnot.
37 ! procedure :: interpolate_from_interpolant_derivative_eta1 => null_interp_one_arg_msg
38  ! yes, we are using the null functions more than once.
39  procedure :: interpolate_from_interpolant_derivative_eta1 => null_fd_1d_arg_msg
40  procedure :: interpolate_from_interpolant_array => null_fd_1d_array_sub
41  procedure :: interpolate_pointer_values => null_fd_1d_ptr_sub
42  procedure :: interpolate_from_interpolant_derivatives_eta1 => null_fd_1d_array_sub
43  procedure :: interpolate_pointer_derivatives => &
45  procedure, pass :: interpolate_array => null_fd_1d_interp_array
46  procedure, pass :: reconstruct_array => null_fd_1d_reconstruct
47  procedure, pass :: interpolate_array_disp => null_fd_1d_array_disp
49 
50 contains
51 
53  interpolator, &
54  num_points, &
55  delta)
56 
57  class(sll_finite_difference_interpolator_1d), intent(inout) :: interpolator
58  sll_int32, intent(in) :: num_points
59  sll_real64, intent(in) :: delta ! cell spacing
60  interpolator%num_points = num_points
61  interpolator%delta = delta
62  interpolator%r_delta = 1.0_f64/delta
64 
66  interpolator, &
67  num_pts, &
68  vals_to_interpolate, &
69  output)
70 
71  class(sll_finite_difference_interpolator_1d), intent(in) :: interpolator
72  sll_int32, intent(in) :: num_pts
73  ! the following input array must be a pointer so that this function
74  ! is usable with splitting methods in higher dimensional data.
75  sll_real64, dimension(:), pointer :: vals_to_interpolate
76  sll_real64, dimension(:), pointer :: output
77  sll_int32 :: i
78  sll_real64 :: r_delta
79  sll_real64, dimension(:), pointer :: dptr
80 
81  sll_assert(size(vals_to_interpolate) >= num_pts)
82  sll_assert(size(output) >= num_pts)
83 
84  r_delta = interpolator%r_delta
85 
86  ! Compute derivative in first point with a forward scheme (-3/2, 2, -1/2)
87  output(1) = r_delta*(-1.5_f64*vals_to_interpolate(1) + &
88  2.0_f64*vals_to_interpolate(2) - &
89  0.5_f64*vals_to_interpolate(3))
90 
91  ! Compute derivative in the bulk of the array with a centered scheme.
92  do i = 2, num_pts - 1
93  output(i) = 0.5_f64*r_delta*(vals_to_interpolate(i + 1) - &
94  vals_to_interpolate(i - 1))
95  end do
96 
97  ! Compute derivative in the last point with a backward scheme (1/2, -2, 3/2)
98  output(num_pts) = r_delta*(0.5_f64*vals_to_interpolate(num_pts - 2) - &
99  2.0_f64*vals_to_interpolate(num_pts - 1) + &
100  1.5_f64*vals_to_interpolate(num_pts))
101 
103 
104  ! Define the functions that are not needed in this particular interpolator
105  ! with the use of the null-function build macros.
106 
107  define_null_interp_one_arg_msg(sll_finite_difference_interpolator_1d, null_fd_1d_arg_msg)
108 
109  define_null_interp_1d_array(sll_finite_difference_interpolator_1d, null_fd_1d_interp_array)
110 
111  define_null_interp_1d_array_msg(sll_finite_difference_interpolator_1d, null_fd_1d_array_msg)
112 
113  define_null_interp_1d_array_sub(sll_finite_difference_interpolator_1d, null_fd_1d_array_sub)
114 
115  define_null_interp_1d_pointer_sub(sll_finite_difference_interpolator_1d, null_fd_1d_ptr_sub)
116 
117  define_null_reconstruct_1d_array(sll_finite_difference_interpolator_1d, null_fd_1d_reconstruct)
118 
119  define_null_interpolate_1d_disp(sll_finite_difference_interpolator_1d, null_fd_1d_array_disp)
120 
Finite differences implementation of sll_c_interpolator_1d.
subroutine finite_difference_1d_interp_derivatives(interpolator, num_pts, vals_to_interpolate, output)
subroutine initialize_finite_difference_1d_interp(interpolator, num_points, delta)
Module for 1D interpolation and reconstruction.
Abstract class for 1D interpolation and reconstruction.
    Report Typos and Errors