Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_fornberg.F90
Go to the documentation of this file.
1 !Copyright (c) 2012-2013, Björn Ingvar Dahlgren
2 !All rights reserved.
3 !
4 !Redistribution and use in source and binary forms, with or without
5 !modification, are permitted provided that the following conditions
6 !are met:
7 !
8 !Redistributions of source code must retain the above copyright
9 !notice, this list of conditions and the following disclaimer.
10 !
11 !Redistributions in binary form must reproduce the above copyright
12 !notice, this list of conditions and the following disclaimer in the
13 !documentation and/or other materials provided with the distribution.
14 !
15 !THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
16 !"AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
17 !LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
18 !FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
19 !COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
20 !INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
21 !BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
22 !LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
23 !CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
24 !LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
25 !ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
26 !POSSIBILITY OF SUCH DAMAGE.
27 !
28 !AUTHORS
29 !Björn Dahlgren
30 !Ondřej Čertík
31 
33 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 #include "sll_working_precision.h"
35 
36  implicit none
37 
38  public :: &
41 
42  private
43 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
44 
45 contains
46 
54  subroutine sll_s_apply_fd(nin, maxorder, xdata, ydata, xtgt, out)
55  sll_int32, intent(in) :: nin, maxorder
56  sll_real64, intent(in) :: xdata(0:nin - 1), ydata(0:nin - 1), xtgt
57  sll_real64, intent(out) :: out(0:maxorder)
58 
59  sll_int32 :: j
60  sll_real64 :: c(0:nin - 1, 0:maxorder)
61 
62  call sll_s_populate_weights(xtgt, xdata, nin - 1, maxorder, c)
63  forall (j=0:maxorder) out(j) = sum(c(0:, j)*ydata)
64 
65  end subroutine
66 
80  subroutine sll_s_populate_weights(z, x, nd, m, c)
81 
82  sll_real64, intent(in) :: z
83  sll_int32, intent(in) :: nd, m
84  sll_real64, intent(in) :: x(0:nd)
85  sll_real64, intent(out) :: c(0:nd, 0:m)
86 
87  sll_int32 :: i, j, k, mn
88  sll_real64 :: c1, c2, c3, c4, c5
89 
90  c1 = 1.0_f64
91  c4 = x(0) - z
92  c = 0.0_f64
93  c(0, 0) = 1.0_f64
94  do i = 1, nd
95  mn = min(i, m)
96  c2 = 1.0_f64
97  c5 = c4
98  c4 = x(i) - z
99  do j = 0, i - 1
100  c3 = x(i) - x(j)
101  c2 = c2*c3
102  if (j == i - 1) then
103  do k = mn, 1, -1
104  c(i, k) = c1*(real(k, kind=8)*c(i - 1, k - 1) - c5*c(i - 1, k))/c2
105  end do
106  c(i, 0) = -c1*c5*c(i - 1, 0)/c2
107  end if
108  do k = mn, 1, -1
109  c(j, k) = (c4*c(j, k) - real(k, kind=8)*c(j, k - 1))/c3
110  end do
111  c(j, 0) = c4*c(j, 0)/c3
112  end do
113  c1 = c2
114  end do
115  end subroutine sll_s_populate_weights
116 
117 end module
subroutine, public sll_s_apply_fd(nin, maxorder, xdata, ydata, xtgt, out)
Apply finite difference formula to compute derivative.
subroutine, public sll_s_populate_weights(z, x, nd, m, c)
    Report Typos and Errors