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_old.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  private
14 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
15 
17  type, abstract :: sll_scalar_field_2d_base
18  ! consider eliminating this transformation from this base class,
19  ! it is already in the derived classes and is confusing...
20  class(sll_c_coordinate_transformation_2d_base), pointer :: coord_trans
21  contains
22  procedure(function_get_mesh), deferred, pass :: get_cartesian_mesh
23  procedure(function_get_transformation), deferred, pass :: &
24  get_transformation
25  procedure(function_get_jacobian_matrix), deferred, pass :: &
26  get_jacobian_matrix
27  procedure(function_evaluation_real), deferred, pass :: value_at_point
28  procedure(function_evaluation_integer), deferred, pass :: value_at_indices
29  procedure(first_derivative_eta1_evaluation_real), deferred, pass :: &
30  first_deriv_eta1_value_at_point
31  procedure(first_derivative_eta2_evaluation_real), deferred, pass :: &
32  first_deriv_eta2_value_at_point
33  procedure(first_derivative_eta1_evaluation_integer), deferred, pass :: &
34  first_deriv_eta1_value_at_indices
35  procedure(first_derivative_eta2_evaluation_integer), deferred, pass :: &
36  first_deriv_eta2_value_at_indices
37  procedure(set_field_data_subroutine), deferred, pass :: set_field_data
38  procedure(field_2d_message_pass), deferred, pass :: &
39  update_interpolation_coefficients
40  procedure(field_2d_file_output), deferred, pass :: write_to_file
41  procedure(field_2d_subroutine), deferred, pass :: delete
42  ! here we can continue with derivatives or whatever else that might
43  ! be desired.
45 
47  class(sll_scalar_field_2d_base), pointer :: base
49 
51  abstract interface
52  function function_get_mesh(field) result(res)
55  class(sll_scalar_field_2d_base), intent(in) :: field
56  class(sll_t_cartesian_mesh_2d), pointer :: res
57  end function function_get_mesh
58  end interface
59 
60  abstract interface
61  subroutine set_field_data_subroutine(field, values)
64  class(sll_scalar_field_2d_base), intent(inout) :: field
65  sll_real64, dimension(:, :), intent(in) :: values
66  end subroutine set_field_data_subroutine
67  end interface
68 
69  abstract interface
70  subroutine field_2d_message_pass(field)
72  class(sll_scalar_field_2d_base), intent(inout) :: field
73  end subroutine field_2d_message_pass
74  end interface
75 
76  abstract interface
77  function function_get_transformation(field) result(res)
80  class(sll_scalar_field_2d_base), intent(in) :: field
81  class(sll_c_coordinate_transformation_2d_base), pointer :: res
82  end function function_get_transformation
83  end interface
84 
85  abstract interface
86  function function_get_jacobian_matrix(field, eta1, eta2) result(res)
89  class(sll_scalar_field_2d_base), intent(in) :: field
90  sll_real64, intent(in) :: eta1
91  sll_real64, intent(in) :: eta2
92  sll_real64, dimension(2, 2) :: res
93  end function function_get_jacobian_matrix
94  end interface
95 
96  abstract interface
97  function function_evaluation_real(field, eta1, eta2) result(res)
100  class(sll_scalar_field_2d_base), intent(in) :: field
101  sll_real64, intent(in) :: eta1
102  sll_real64, intent(in) :: eta2
103  sll_real64 :: res
104  end function function_evaluation_real
105  end interface
106 
107  abstract interface
108  function function_evaluation_integer(field, i, j) result(res)
111  class(sll_scalar_field_2d_base), intent(in) :: field
112  sll_int32, intent(in) :: i
113  sll_int32, intent(in) :: j
114  sll_real64 :: res
115  end function function_evaluation_integer
116  end interface
117 
118  abstract interface
119  function first_derivative_eta1_evaluation_real(field, eta1, eta2) result(res)
122  class(sll_scalar_field_2d_base), intent(in) :: field
123  sll_real64, intent(in) :: eta1
124  sll_real64, intent(in) :: eta2
125  sll_real64 :: res
127  end interface
128 
129  abstract interface
130  function first_derivative_eta2_evaluation_real(field, eta1, eta2) result(res)
133  class(sll_scalar_field_2d_base), intent(in) :: field
134  sll_real64, intent(in) :: eta1
135  sll_real64, intent(in) :: eta2
136  sll_real64 :: res
138  end interface
139 
140  abstract interface
141  function first_derivative_eta1_evaluation_integer(field, i, j) result(res)
144  class(sll_scalar_field_2d_base), intent(in) :: field
145  sll_int32, intent(in) :: i
146  sll_int32, intent(in) :: j
147  sll_real64 :: res
149  end interface
150 
151  abstract interface
152  function first_derivative_eta2_evaluation_integer(field, i, j) result(res)
155  class(sll_scalar_field_2d_base), intent(in) :: field
156  sll_int32, intent(in) :: i
157  sll_int32, intent(in) :: j
158  sll_real64 :: res
160  end interface
161 
162  abstract interface
163  function return_integer(field) result(res)
166  class(sll_scalar_field_2d_base), intent(in) :: field
167  sll_int32 :: res
168  end function return_integer
169  end interface
170 
171  abstract interface
172  subroutine field_2d_file_output(field, tag)
175  class(sll_scalar_field_2d_base), intent(in) :: field
176  sll_int32, intent(in) :: tag
177  end subroutine field_2d_file_output
178  end interface
179 
180  abstract interface
181  subroutine field_2d_subroutine(field)
183  class(sll_scalar_field_2d_base), intent(inout) :: field
184  end subroutine field_2d_subroutine
185  end interface
186 
Cartesian mesh basic types.
Abstract class for coordinate transformations.
Module to select the kind parameter.
    Report Typos and Errors