Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_landau_2d_initializer.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_assert.h"
4 #include "sll_working_precision.h"
5 
6  use sll_m_cartesian_meshes, only: &
8 
9  use sll_m_constants, only: &
10  sll_p_pi
11 
14 
19 
20  implicit none
21 
22  public :: &
24 
25  private
26 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 
29  class(sll_c_coordinate_transformation_2d_base), pointer :: transf
30  sll_real64 :: eps
31  sll_real64 :: kx
32  sll_int32 :: is_delta_f
33  contains
34  procedure, pass(init_obj) :: init => initialize_landau_2d
35  procedure, pass(init_obj) :: f_of_x1x2 => f_x1x2_landau_2d
36  end type sll_t_init_landau_2d
37 
38 contains
39 
40  subroutine initialize_landau_2d(init_obj, transf, data_position, eps_val, &
41  kx_val, is_delta_f)
42  class(sll_t_init_landau_2d), intent(inout) :: init_obj
43  class(sll_c_coordinate_transformation_2d_base), pointer :: transf
44  sll_int32 :: data_position
45  sll_real64, intent(in), optional :: eps_val
46  sll_real64, intent(in), optional :: kx_val
47  sll_int32, intent(in), optional :: is_delta_f
48 
49  init_obj%data_position = data_position
50  if (present(eps_val)) then
51  init_obj%eps = eps_val
52  else
53  init_obj%eps = 0.01_f64 ! just some default value
54  end if
55  if (present(kx_val)) then
56  init_obj%kx = kx_val
57  else
58  init_obj%kx = 0.5_f64 ! just some default value
59  end if
60  if (present(is_delta_f)) then
61  init_obj%is_delta_f = is_delta_f
62  else
63  init_obj%is_delta_f = 1 ! default value is false
64  end if
65  init_obj%transf => transf
66  ! kx remains uninitialized because we need mesh information
67  end subroutine
68 
69  subroutine f_x1x2_landau_2d(init_obj, data_out)
70  class(sll_t_init_landau_2d), intent(inout) :: init_obj
71  class(sll_c_coordinate_transformation_2d_base), pointer :: transf
72  sll_real64, dimension(:, :), intent(out) :: data_out
73  sll_int32 :: i
74  sll_int32 :: j
75  sll_int32 :: num_pts1
76  sll_int32 :: num_pts2
77  sll_real64 :: eps
78  sll_real64 :: kx
79  sll_real64 :: x
80  sll_real64 :: v
81 
82  eps = init_obj%eps
83  transf => init_obj%transf
84 
85  if (init_obj%data_position == sll_p_node_centered_field) then
86  num_pts1 = transf%mesh%num_cells1 + 1
87  num_pts2 = transf%mesh%num_cells2 + 1
88  else if (init_obj%data_position == sll_p_cell_centered_field) then
89  num_pts1 = transf%mesh%num_cells1
90  num_pts2 = transf%mesh%num_cells2
91  end if
92  kx = init_obj%kx
93  sll_assert(size(data_out, 1) .ge. num_pts1)
94  sll_assert(size(data_out, 2) .ge. num_pts2)
95  do j = 1, num_pts2
96  do i = 1, num_pts1
97  if (init_obj%data_position == sll_p_node_centered_field) then
98  v = transf%x2_at_node(i, j)
99  x = transf%x1_at_node(i, j)
100  else if (init_obj%data_position == sll_p_cell_centered_field) then
101  v = transf%x2_at_cell(i, j)
102  x = transf%x1_at_cell(i, j)
103  else
104  print *, 'f_x1x2_landau_2d:', init_obj%data_position, 'not defined'
105  end if
106  if (init_obj%is_delta_f == 0) then ! delta_f code
107  data_out(i, j) = eps*cos(kx*x)/sqrt(2*sll_p_pi)*exp(-0.5_f64*v*v)
108  else ! full f
109  data_out(i, j) = &
110  (1.0_f64 + eps*cos(kx*x))/sqrt(2*sll_p_pi)*exp(-0.5_f64*v*v)
111  end if
112  end do
113  end do
114  end subroutine
115 
Cartesian mesh basic types.
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_pi
Abstract class for coordinate transformations.
subroutine initialize_landau_2d(init_obj, transf, data_position, eps_val, kx_val, is_delta_f)
subroutine f_x1x2_landau_2d(init_obj, data_out)
    Report Typos and Errors