10 #include "sll_assert.h"
17 sll_f_horner_1d_eval, &
18 sll_f_horner_2d_eval, &
41 sll_pure
function sll_f_horner_1d_eval(coeffs, x, deriv)
result(res)
42 real(
wp),
intent(in) :: coeffs(:)
43 real(
wp),
intent(in) :: x
44 integer,
intent(in) :: deriv
50 n1 = ubound(coeffs, 1)
52 fmmjdr = real(n1 - 1 - deriv,
wp)
54 do i1 = 1, n1 - 1 - deriv
55 res = (res/fmmjdr)*x + coeffs(n1 - i1)
56 fmmjdr = fmmjdr - 1.0_wp
59 end function sll_f_horner_1d_eval
70 sll_pure
function sll_f_horner_2d_eval(coeffs, x, deriv)
result(res)
71 real(
wp),
intent(in) :: coeffs(:, :)
72 real(
wp),
intent(in) :: x(:)
73 integer,
intent(in) :: deriv(:)
76 real(
wp) :: coeffi(ubound(coeffs, 2))
79 n2 = ubound(coeffs, 2)
82 coeffi(i2) = sll_f_horner_1d_eval(coeffs(:, i2), x(1), deriv(1))
85 res = sll_f_horner_1d_eval(coeffi, x(2), deriv(2))
87 end function sll_f_horner_2d_eval
98 sll_pure
function sll_f_horner_3d_eval(coeffs, x, deriv)
result(res)
99 real(
wp),
intent(in) :: coeffs(:, :, :)
100 real(
wp),
intent(in) :: x(:)
101 integer,
intent(in) :: deriv(:)
104 real(
wp) :: coeffi(ubound(coeffs, 3))
107 n3 = ubound(coeffs, 3)
110 coeffi(i3) = sll_f_horner_2d_eval(coeffs(:, :, i3), x(1:2), deriv(1:2))
113 res = sll_f_horner_1d_eval(coeffi, x(3), deriv(3))
115 end function sll_f_horner_3d_eval
module for evaluation of a polynomial in pp form using Horner's algorithm
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)