Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_initial_distribution_functions.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_working_precision.h"
4 
5  use sll_m_constants, only: &
6  sll_kx, &
8 
9  implicit none
10 
11  private
12 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
13 
14 contains
15  elemental function landau(x, v)
16  sll_real64, intent(in) :: x, v
17  sll_real64 :: landau
18  sll_real64, parameter :: eps = 0.01_f64
19 
20  ! sll_kx is defined in the sll_m_constants module.
21  ! It is set at the mesh initialization
22  landau = (1.0_f64 + eps*cos(sll_kx*x)) &
23  /sqrt(2.0_f64*sll_p_pi)*exp(-0.5_f64*v*v)
24  end function landau
25 
26  elemental function two_stream(x, v) result(fval)
27  sll_real64, intent(in) :: x, v
28  sll_real64 :: fval
29  ! local variables
30  sll_real64, parameter :: eps = 0.01_f64
31  sll_real64, parameter :: xi = 0.90_f64
32  sll_real64, parameter :: v0 = 2.4_f64
33  sll_real64 :: vv
34  ! sll_kx is defined in the sll_m_constants module.
35  ! It is set at the mesh initialization
36  vv = v*v
37  fval = (1.0_f64 + eps*((cos(2*sll_kx*x) + cos(3*sll_kx*x))/1.2_f64 + cos(sll_kx*x)))* &
38  (1.0_f64/sqrt(2*sll_p_pi))*((2 - 2*xi)/(3 - 2*xi))* &
39  (1.0_f64 + .5_f64*vv/(1 - xi))*exp(-.5_f64*vv)
40  !fval=(1+eps*eps*cos(sll_kx*x)+eps*cos(2*sll_kx*x))*0.5_f64/sqrt(2*pi) &
41  ! *(exp(-.5_f64*(vx-v0)**2)+ exp(-.5_f64*(vx+v0)**2))
42  !fval=(1+eps*cos(sll_kx*x))*0.5_f64/sqrt(2*pi)*(exp(-.5_f64*(vx-v0)**2) &
43  ! + exp(-.5_f64*(vx+v0)**2))
44  !fval= 1/sqrt (2 * sll_p_pi) * ( 0.9_f64 * exp (-.5_f64 * vx * vx) &
45  ! + 0.2_f64 * exp(-0.5_f64 * (vx - 4.5_f64)*(vx - 4.5_f64) &
46  ! /(0.5_f64 * 0.5_f64))) * (1. + 0.03_f64 * cos (0.3_f64 * x))
47  ! fval=(1+eps*cos(sll_kx*x))*1/sqrt(2*pi)*exp(-.5_f64*vv)
48  ! fval=exp(-.5_f64*(xx+vv))
49  end function two_stream
50 
51  ! Samuel : I've change elemental by pure because there is an error with
52  ! gcc 4.7.0 in the unit_test because he does'nt like pointer to elemental
53  ! function.
54 
55  pure function sll_m_gaussian(x, v) result(fval)
56  sll_real64, intent(in) :: x, v
57  sll_real64 :: fval
58  ! local variables
59  sll_real64, parameter :: xoffset = 0.5_f64
60  sll_real64, parameter :: voffset = 0.5_f64
61  sll_real64 :: xx, vv
62 
63  xx = (x - xoffset)**2
64  vv = (v - voffset)**2
65  fval = exp(-0.5_f64*(xx + vv)*40._f64)
66  end function sll_m_gaussian
67 
68 end module
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_pi
sll_m_gaussian random generator
elemental real(kind=f64) function two_stream(x, v)
elemental real(kind=f64) function landau(x, v)
    Report Typos and Errors