Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_primitives.F90
Go to the documentation of this file.
1 
4 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
5 #include "sll_working_precision.h"
6 
7  implicit none
8 
9  public :: &
12 
13  private
14 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
15 
16 contains
17 
23  subroutine sll_s_function_to_primitive(f, x, n, xm)
24 
25  sll_real64, dimension(:), intent(inout) :: f
26  sll_real64, dimension(:), intent(in) :: x
27  sll_int32, intent(in) :: n
28  sll_real64, intent(out) :: xm
29 
30  sll_int32 :: i
31  sll_real64 :: dx
32  sll_real64 :: tmp
33  sll_real64 :: tmp2
34 
35  dx = 1._f64/real(n, f64)
36 
37 !from f compute the mean
38  xm = 0.0_f64
39  do i = 1, n
40  xm = xm + f(i)*(x(i + 1) - x(i))
41  end do
42 
43  f(1) = (f(1) - xm)*(x(2) - x(1))
44  tmp = f(1)
45  f(1) = 0.0_f64
46  do i = 2, n
47  f(i) = (f(i) - xm)*(x(i + 1) - x(i))
48  tmp2 = f(i)
49  f(i) = f(i - 1) + tmp
50  tmp = tmp2
51  end do
52  f(n + 1) = f(n) + tmp
53 
54  end subroutine sll_s_function_to_primitive
55 
61  subroutine sll_s_primitive_to_function(f, x, n, xm)
62 
63  sll_real64, dimension(:), intent(inout) :: f
64  sll_real64, dimension(:), intent(in) :: x
65  sll_int32, intent(in) :: n
66  sll_real64, intent(in) :: xm
67 
68  sll_int32 :: i
69  sll_real64 :: tmp
70  sll_real64 :: dx
71 
72  dx = 1._f64/real(n, f64)
73 
74  tmp = f(1)
75  do i = 1, n - 1
76  f(i) = f(i + 1) - f(i) + xm*(x(i + 1) - x(i))
77  end do
78  f(n) = tmp - f(n) + xm*(x(1) + 1._f64 - x(n))
79 
80 !from mean compute f
81  do i = 1, n
82  f(i) = f(i)/(x(i + 1) - x(i))
83  end do
84 
85  f(n + 1) = f(1)
86 
87  end subroutine sll_s_primitive_to_function
88 
89 end module sll_m_primitives
Functions to compute primitive of 1d function.
subroutine, public sll_s_primitive_to_function(f, x, n, xm)
Compute function value from primitive.
subroutine, public sll_s_function_to_primitive(f, x, n, xm)
Compute primitive of f.
    Report Typos and Errors