Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_scalar_field_2d_base.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_working_precision.h"
4 
5  use sll_m_cartesian_meshes, only: &
7 
10 
11  implicit none
12 
13  public :: &
15 
16  private
17 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
18 
20  type, abstract :: sll_c_scalar_field_2d_base
21 
22  contains
23  procedure(function_get_mesh), deferred, pass :: get_cartesian_mesh
24  procedure(function_get_transformation), deferred, pass :: &
25  get_transformation
26  procedure(function_get_jacobian_matrix), deferred, pass :: &
27  get_jacobian_matrix
28  procedure(function_evaluation_real), deferred, pass :: value_at_point
29  procedure(function_evaluation_integer), deferred, pass :: value_at_indices
30  procedure(first_derivative_eta1_evaluation_real), deferred, pass :: &
31  first_deriv_eta1_value_at_point
32  procedure(first_derivative_eta2_evaluation_real), deferred, pass :: &
33  first_deriv_eta2_value_at_point
34  procedure(first_derivative_eta1_evaluation_integer), deferred, pass :: &
35  first_deriv_eta1_value_at_indices
36  procedure(first_derivative_eta2_evaluation_integer), deferred, pass :: &
37  first_deriv_eta2_value_at_indices
38  procedure(set_field_data_subroutine), deferred, pass :: set_field_data
39  procedure(field_2d_message_pass), deferred, pass :: &
40  update_interpolation_coefficients
41  procedure(field_2d_file_output), deferred, pass :: write_to_file
42  procedure(field_2d_subroutine), deferred, pass :: free
43  ! here we can continue with derivatives or whatever else that might
44  ! be desired.
46 
48  class(sll_c_scalar_field_2d_base), pointer :: base
50 
52  abstract interface
53  function function_get_mesh(field) result(res)
56  class(sll_c_scalar_field_2d_base), intent(in) :: field
57  class(sll_t_cartesian_mesh_2d), pointer :: res
58  end function function_get_mesh
59  end interface
60 
61  abstract interface
62  subroutine set_field_data_subroutine(field, values)
65  class(sll_c_scalar_field_2d_base), intent(inout) :: field
66  sll_real64, dimension(:, :), intent(in) :: values
67  end subroutine set_field_data_subroutine
68  end interface
69 
70  abstract interface
71  subroutine field_2d_message_pass(field)
73  class(sll_c_scalar_field_2d_base), intent(inout) :: field
74  end subroutine field_2d_message_pass
75  end interface
76 
77  abstract interface
78  function function_get_transformation(field) result(res)
81  class(sll_c_scalar_field_2d_base), intent(in) :: field
82  class(sll_c_coordinate_transformation_2d_base), pointer :: res
83  end function function_get_transformation
84  end interface
85 
86  abstract interface
87  function function_get_jacobian_matrix(field, eta1, eta2) result(res)
90  class(sll_c_scalar_field_2d_base), intent(in) :: field
91  sll_real64, intent(in) :: eta1
92  sll_real64, intent(in) :: eta2
93  sll_real64, dimension(2, 2) :: res
94  end function function_get_jacobian_matrix
95  end interface
96 
97  abstract interface
98  function function_evaluation_real(field, eta1, eta2) result(res)
101  class(sll_c_scalar_field_2d_base), intent(in) :: field
102  sll_real64, intent(in) :: eta1
103  sll_real64, intent(in) :: eta2
104  sll_real64 :: res
105  end function function_evaluation_real
106  end interface
107 
108  abstract interface
109  function function_evaluation_integer(field, i, j) result(res)
112  class(sll_c_scalar_field_2d_base), intent(in) :: field
113  sll_int32, intent(in) :: i
114  sll_int32, intent(in) :: j
115  sll_real64 :: res
116  end function function_evaluation_integer
117  end interface
118 
119  abstract interface
120  function first_derivative_eta1_evaluation_real(field, eta1, eta2) result(res)
123  class(sll_c_scalar_field_2d_base), intent(in) :: field
124  sll_real64, intent(in) :: eta1
125  sll_real64, intent(in) :: eta2
126  sll_real64 :: res
128  end interface
129 
130  abstract interface
131  function first_derivative_eta2_evaluation_real(field, eta1, eta2) result(res)
134  class(sll_c_scalar_field_2d_base), intent(in) :: field
135  sll_real64, intent(in) :: eta1
136  sll_real64, intent(in) :: eta2
137  sll_real64 :: res
139  end interface
140 
141  abstract interface
142  function first_derivative_eta1_evaluation_integer(field, i, j) result(res)
145  class(sll_c_scalar_field_2d_base), intent(in) :: field
146  sll_int32, intent(in) :: i
147  sll_int32, intent(in) :: j
148  sll_real64 :: res
150  end interface
151 
152  abstract interface
153  function first_derivative_eta2_evaluation_integer(field, i, j) result(res)
156  class(sll_c_scalar_field_2d_base), intent(in) :: field
157  sll_int32, intent(in) :: i
158  sll_int32, intent(in) :: j
159  sll_real64 :: res
161  end interface
162 
163  abstract interface
164  function return_integer(field) result(res)
167  class(sll_c_scalar_field_2d_base), intent(in) :: field
168  sll_int32 :: res
169  end function return_integer
170  end interface
171 
172  abstract interface
173  subroutine field_2d_file_output(field, tag)
176  class(sll_c_scalar_field_2d_base), intent(in) :: field
177  sll_int32, intent(in) :: tag
178  end subroutine field_2d_file_output
179  end interface
180 
181  abstract interface
182  subroutine field_2d_subroutine(field)
184  class(sll_c_scalar_field_2d_base), intent(inout) :: field
185  end subroutine field_2d_subroutine
186  end interface
187 
Cartesian mesh basic types.
Abstract class for coordinate transformations.
Module to select the kind parameter.
    Report Typos and Errors