Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_coordinate_transformation_2d_base.F90
Go to the documentation of this file.
1 
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 #include "sll_memory.h"
24 #include "sll_working_precision.h"
25 
26  use sll_m_cartesian_meshes, only: &
28 
29  implicit none
30 
31  public :: &
33  sll_p_io_gmsh, &
35  sll_p_io_mtv, &
36  sll_p_io_vtk, &
37  sll_p_io_xdmf, &
38  sll_i_transformation_func_nopass
39 
40  private
41 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
42 
43  sll_int32, parameter :: sll_p_io_xdmf = 0
44  sll_int32, parameter :: sll_p_io_vtk = 1
45  sll_int32, parameter :: sll_p_io_gnuplot = 2
46  sll_int32, parameter :: sll_p_io_mtv = 3
47  sll_int32, parameter :: sll_p_io_gmsh = 4
48 
58  type(sll_t_cartesian_mesh_2d), pointer :: mesh => null()
59 ! type(sll_t_cartesian_mesh_2d), pointer :: mesh2d_minimal => null()
60  !logical to remember when the mesh has already been written to file
62  character(len=64) :: label
64  logical :: written = .false.
65  contains
67  procedure(get_cartesian_mesh_ct), deferred, pass :: get_cartesian_mesh
69  procedure(geometry_function_ct), deferred, pass :: x1
71  procedure(geometry_function_ct), deferred, pass :: x2
73  procedure(geometry_function_ct), deferred, pass :: jacobian
75  procedure(geometry_function_indices_ct), deferred, pass :: x1_at_node
77  procedure(geometry_function_indices_ct), deferred, pass :: x2_at_node
79  procedure(geometry_function_indices_ct), deferred, pass :: jacobian_at_node
81  procedure(matrix_geometry_function_ct), deferred, pass :: jacobian_matrix
83  procedure(matrix_geometry_function_ct), deferred, pass :: &
84  inverse_jacobian_matrix
86  procedure(geometry_function_indices_ct), deferred, pass :: x1_at_cell
88  procedure(geometry_function_indices_ct), deferred, pass :: x2_at_cell
90  procedure(geometry_function_indices_ct), deferred, pass :: jacobian_at_cell
92  procedure(write_transformation_signature), deferred, pass :: write_to_file
94  procedure(read_transformation), deferred, pass :: read_from_file
96  procedure(transformation_subroutine), deferred, pass :: delete
98 
99  !************************************************************************
100  !
101  ! Function signatures
102  !
103  !************************************************************************
104 
105  !************************************************************************
106  ! 2D CASE:
107  !************************************************************************
108 
109 #ifndef DOXYGEN_SHOULD_SKIP_THIS
110 
111  abstract interface
112  function get_cartesian_mesh_ct(transf) result(res)
115  class(sll_c_coordinate_transformation_2d_base), intent(in) :: transf
116  class(sll_t_cartesian_mesh_2d), pointer :: res
117  end function get_cartesian_mesh_ct
118  end interface
119 
120  abstract interface
121  function geometry_function_ct(transf, eta1, eta2) result(res)
125  sll_real64, intent(in) :: eta1
126  sll_real64, intent(in) :: eta2
127  sll_real64 :: res
128  end function geometry_function_ct
129  end interface
130 
131  abstract interface
132  function geometry_function_indices_ct(transf, i, j) result(res)
136  sll_int32, intent(in) :: i
137  sll_int32, intent(in) :: j
138  sll_real64 :: res
139  end function geometry_function_indices_ct
140  end interface
141 
142  abstract interface
143  function matrix_geometry_function_ct(transf, eta1, eta2)
146  class(sll_c_coordinate_transformation_2d_base), intent(inout) :: transf
147  sll_real64, intent(in) :: eta1
148  sll_real64, intent(in) :: eta2
149  sll_real64 :: matrix_geometry_function_ct(2, 2)
150  end function matrix_geometry_function_ct
151  end interface
152 
153  ! It is mandatory to pass the params array since making it optional
154  ! poses problems when storing the parameters array inside the object to
155  ! pass it when the function is called.
156  abstract interface
157  function sll_i_transformation_func_nopass(eta1, eta2, params) result(res)
159  sll_real64, intent(in) :: eta1
160  sll_real64, intent(in) :: eta2
161  sll_real64, dimension(:), intent(in) :: params
162  sll_real64 :: res
163  end function sll_i_transformation_func_nopass
164  end interface
165 
166  abstract interface
167  function matrix_geometry_function_nopass(eta1, eta2) result(res)
169  sll_real64, intent(in) :: eta1
170  sll_real64, intent(in) :: eta2
171  sll_real64 :: res(2, 2)
172  end function matrix_geometry_function_nopass
173  end interface
174 
175  ! WE SHOULD PROBABLY HAVE A SINGLE FILE WITH ALL THE SIGNATURES THAT WE
176  ! GENERALLY USE AND DO NOT DEPEND ON A BASE OR DERIVED TYPE.
177 
178  ! The following interface is meant to specify the signature of the
179  ! (possibly user-defined) functions that should be passed as arguments
180  ! to initialize the analytic transformations.
181  abstract interface
182  function two_arg_scalar_function(eta1, eta2)
184  sll_real64 :: two_arg_scalar_function
185  sll_real64, intent(in) :: eta1
186  sll_real64, intent(in) :: eta2
187  end function two_arg_scalar_function
188  end interface
189 
190  abstract interface
191  function two_arg_message_passing_func(transf, eta1, eta2)
194  sll_real64 :: two_arg_message_passing_func
196  sll_real64, intent(in) :: eta1
197  sll_real64, intent(in) :: eta2
198  end function two_arg_message_passing_func
199  end interface
200 
201  abstract interface
202  subroutine write_transformation_signature(transf, output_format)
205  class(sll_c_coordinate_transformation_2d_base), intent(inout) :: transf
206  sll_int32, intent(in), optional :: output_format
207  end subroutine write_transformation_signature
208  end interface
209 
210  abstract interface
211  subroutine transformation_subroutine(transf)
213  class(sll_c_coordinate_transformation_2d_base), intent(inout) :: transf
214  end subroutine transformation_subroutine
215  end interface
216 
217  abstract interface
218  subroutine read_transformation(transf, filename)
220  class(sll_c_coordinate_transformation_2d_base), intent(inout) :: transf
221  character(len=*), intent(in) :: filename
222  end subroutine read_transformation
223  end interface
224 
225 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
226 
Cartesian mesh basic types.
Abstract class for coordinate transformations.
Module to select the kind parameter.
    Report Typos and Errors