Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_array_plotting.F90
Go to the documentation of this file.
1 
3 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
4 #include "sll_memory.h"
5 #include "sll_working_precision.h"
6 
7  use sll_m_cartesian_meshes, only: &
9 
10  use sll_m_gnuplot, only: &
12 
13  implicit none
14 
15  public :: &
16  sll_p_x1x2, &
17  sll_p_x1x3, &
18  sll_p_x1x4, &
19  sll_p_x2x3, &
20  sll_p_x2x4, &
21  sll_p_x3x4, &
23 
24  private
25 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
26 
27  sll_int32, parameter :: sll_p_x1x2 = 0
28  sll_int32, parameter :: sll_p_x1x3 = 1
29  sll_int32, parameter :: sll_p_x1x4 = 2
30  sll_int32, parameter :: sll_p_x2x3 = 3
31  sll_int32, parameter :: sll_p_x2x4 = 4
32  sll_int32, parameter :: sll_p_x3x4 = 5
33 
34 contains
35 
47  subroutine sll_s_write_projection_2d(mesh4d, &
48  array, &
49  label, &
50  projection, &
51  slice, &
52  iplot)
53 
54  type(sll_t_cartesian_mesh_4d), intent(in) :: mesh4d
55  character(len=*), intent(in) :: label
56  sll_real64, intent(in) :: array(:, :, :, :)
57  sll_int32, intent(in) :: projection
58  sll_int32, intent(in) :: slice(2)
59  sll_int32, intent(in) :: iplot
60 
61  sll_real64, pointer :: f_split(:, :)
62 
63  sll_real64 :: eta1_min, eta2_min, eta3_min, eta4_min
64  sll_real64 :: eta1_max, eta2_max, eta3_max, eta4_max
65  sll_int32 :: n1, n2, n3, n4
66  sll_int32 :: error
67 
68  n1 = mesh4d%num_cells1 + 1
69  n2 = mesh4d%num_cells2 + 1
70  n3 = mesh4d%num_cells3 + 1
71  n4 = mesh4d%num_cells4 + 1
72 
73  eta1_min = mesh4d%eta1_min; eta1_max = mesh4d%eta1_max
74  eta2_min = mesh4d%eta2_min; eta2_max = mesh4d%eta2_max
75  eta3_min = mesh4d%eta3_min; eta3_max = mesh4d%eta3_max
76  eta4_min = mesh4d%eta4_min; eta4_max = mesh4d%eta4_max
77 
78  select case (projection)
79 
80  case (sll_p_x1x2)
81 
82  sll_allocate(f_split(n1, n2), error)
83  f_split = array(:, :, slice(1), slice(2))
84  call sll_o_gnuplot_2d(eta1_min, eta1_max, n1, &
85  eta2_min, eta2_max, n2, &
86  f_split, &
87  label, iplot, error)
88  case (sll_p_x1x3)
89 
90  sll_allocate(f_split(n1, n3), error)
91  f_split = array(:, slice(1), :, slice(2))
92  call sll_o_gnuplot_2d(eta1_min, eta1_max, n1, &
93  eta3_min, eta3_max, n3, &
94  f_split, &
95  label, iplot, error)
96  case (sll_p_x1x4)
97 
98  sll_allocate(f_split(n1, n4), error)
99  f_split = array(:, slice(1), slice(2), :)
100  call sll_o_gnuplot_2d(eta1_min, eta1_max, n1, &
101  eta4_min, eta4_max, n4, &
102  f_split, &
103  label, iplot, error)
104  case (sll_p_x2x3)
105 
106  sll_allocate(f_split(n2, n3), error)
107  f_split = array(slice(1), :, :, slice(2))
108  call sll_o_gnuplot_2d(eta2_min, eta2_max, n2, &
109  eta3_min, eta3_max, n3, &
110  f_split, &
111  label, iplot, error)
112  case (sll_p_x2x4)
113 
114  sll_allocate(f_split(n2, n4), error)
115  f_split = array(slice(1), :, slice(2), :)
116  call sll_o_gnuplot_2d(eta2_min, eta2_max, n2, &
117  eta4_min, eta4_max, n4, &
118  f_split, &
119  label, iplot, error)
120  case (sll_p_x3x4)
121 
122  sll_allocate(f_split(n3, n4), error)
123  f_split = array(slice(1), slice(2), :, :)
124  call sll_o_gnuplot_2d(eta3_min, eta3_max, n3, &
125  eta4_min, eta4_max, n4, &
126  f_split, &
127  label, iplot, error)
128 
129  end select
130 
131  end subroutine sll_s_write_projection_2d
132 
133 end module sll_m_array_plotting
Write file for gnuplot to display 2d field.
integer(kind=i32), parameter, public sll_p_x1x3
subroutine, public sll_s_write_projection_2d(mesh4d, array, label, projection, slice, iplot)
Write a gnuplot file to display 2d projection of 4d array.
integer(kind=i32), parameter, public sll_p_x1x2
integer(kind=i32), parameter, public sll_p_x2x3
integer(kind=i32), parameter, public sll_p_x3x4
integer(kind=i32), parameter, public sll_p_x2x4
integer(kind=i32), parameter, public sll_p_x1x4
Cartesian mesh basic types.
Implements the functions to write data file plotable by GNUplot.
    Report Typos and Errors