Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_horner.F90
Go to the documentation of this file.
1 
4 !
7 !
9 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
10 #include "sll_assert.h"
11 
12  use sll_m_working_precision, only: f64
13 
14  implicit none
15 
16  public :: &
17  sll_f_horner_1d_eval, &
18  sll_f_horner_2d_eval, &
19  sll_f_horner_3d_eval
20 
21  private
22 
23 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
24 
26  integer, parameter :: wp = f64
27 
28 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
29 contains
30 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 
32  !-----------------------------------------------------------------------------
40  !-----------------------------------------------------------------------------
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
45  real(wp) :: res
46 
47  real(wp) :: fmmjdr
48  integer :: i1, n1
49 
50  n1 = ubound(coeffs, 1)
51 
52  fmmjdr = real(n1 - 1 - deriv, wp)
53  res = coeffs(n1)
54  do i1 = 1, n1 - 1 - deriv
55  res = (res/fmmjdr)*x + coeffs(n1 - i1)
56  fmmjdr = fmmjdr - 1.0_wp
57  end do
58 
59  end function sll_f_horner_1d_eval
60 
61  !-----------------------------------------------------------------------------
69  !-----------------------------------------------------------------------------
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(:)
74  real(wp) :: res
75 
76  real(wp) :: coeffi(ubound(coeffs, 2))
77  integer :: i2, n2
78 
79  n2 = ubound(coeffs, 2)
80 
81  do i2 = 1, n2
82  coeffi(i2) = sll_f_horner_1d_eval(coeffs(:, i2), x(1), deriv(1))
83  end do
84 
85  res = sll_f_horner_1d_eval(coeffi, x(2), deriv(2))
86 
87  end function sll_f_horner_2d_eval
88 
89  !-----------------------------------------------------------------------------
97  !-----------------------------------------------------------------------------
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(:)
102  real(wp) :: res
103 
104  real(wp) :: coeffi(ubound(coeffs, 3))
105  integer :: i3, n3
106 
107  n3 = ubound(coeffs, 3)
108 
109  do i3 = 1, n3
110  coeffi(i3) = sll_f_horner_2d_eval(coeffs(:, :, i3), x(1:2), deriv(1:2))
111  end do
112 
113  res = sll_f_horner_1d_eval(coeffi, x(3), deriv(3))
114 
115  end function sll_f_horner_3d_eval
116 
117 end module sll_m_horner
module for evaluation of a polynomial in pp form using Horner's algorithm
Definition: sll_m_horner.F90:8
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