23 #include "sll_assert.h"
42 real(
wp),
allocatable :: bcoef(:)
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
76 associate(n => bsplines%ncells, p => bsplines%degree)
78 allocate (self%bcoef(1:n + p))
95 deallocate (self%bcoef)
106 sll_pure
function f_spline_1d__belongs_to_space(self, bsplines)
result(in_space)
109 class(sll_c_bsplines),
intent(in),
target :: bsplines
112 in_space =
associated(self%bspl, bsplines)
114 end function f_spline_1d__belongs_to_space
122 sll_pure
function f_spline_1d__eval(self, x)
result(y)
125 real(wp),
intent(in) :: x
128 real(wp) :: values(1:self%bspl%degree + 1)
129 integer :: jmin, jmax
131 call self%bspl%eval_basis(x, values, jmin)
133 jmax = jmin + self%bspl%degree
135 y = dot_product(self%bcoef(jmin:jmax), values)
137 end function f_spline_1d__eval
145 sll_pure
function f_spline_1d__eval_deriv(self, x)
result(y)
148 real(wp),
intent(in) :: x
151 real(wp) :: derivs(1:self%bspl%degree + 1)
152 integer :: jmin, jmax
154 call self%bspl%eval_deriv(x, derivs, jmin)
156 jmax = jmin + self%bspl%degree
158 y = dot_product(self%bcoef(jmin:jmax), derivs)
160 end function f_spline_1d__eval_deriv
168 sll_pure
subroutine s_spline_1d__eval_array(self, x, y)
171 real(wp),
intent(in) :: x(:)
172 real(wp),
intent(out) :: y(:)
176 sll_assert(
size(x) ==
size(y))
179 y(i) = f_spline_1d__eval(self, x(i))
182 end subroutine s_spline_1d__eval_array
190 sll_pure
subroutine s_spline_1d__eval_array_deriv(self, x, y)
193 real(wp),
intent(in) :: x(:)
194 real(wp),
intent(out) :: y(:)
198 sll_assert(
size(x) ==
size(y))
201 y(i) = f_spline_1d__eval_deriv(self, x(i))
204 end subroutine s_spline_1d__eval_array_deriv
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)
Abstract type, B-splines.