Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_scalar_field_1d_old.F90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------
2 ! SELALIB
3 !------------------------------------------------------------------------------
4 !
5 ! MODULE: sll_scalar_field_1d
6 !
12 !
13 ! DESCRIPTION:
14 !
26 !
27 ! REVISION HISTORY:
28 ! DD Mmm YYYY - Initial Version
29 ! TODO_dd_mmm_yyyy - TODO_describe_appropriate_changes - TODO_name
30 !------------------------------------------------------------------------------
32 
33 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 #include "sll_errors.h"
35 #include "sll_memory.h"
36 #include "sll_working_precision.h"
37 
38  use sll_m_cartesian_meshes, only: &
40 
41  use sll_m_gnuplot, only: &
43 
47 
48  implicit none
49 
50  private
51 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
52 
53  !I removed this line because, it not exists in 2d version
54  !should be in base module
55  !integer, parameter :: sll_p_node_centered_field = 0, sll_p_cell_centered_field = 1
56 
58  class(sll_t_cartesian_mesh_1d), pointer :: mesh
59  sll_real64, dimension(:), pointer :: data
60  sll_int32 :: data_position
61  character(len=64) :: name
62  end type scalar_field_1d
63 
64  abstract interface
65  function scalar_function_1d(eta1)
67  sll_real64 :: scalar_function_1d
68  sll_real64, intent(in) :: eta1
69  end function scalar_function_1d
70  end interface
71 
72 contains ! *i****************************************************************
73  ! this used to be new_scalar_field_1d
75  this, &
76  field_name, &
77  mesh, &
78  data_position, &
79  init_function)
80 
81  class(scalar_field_1d), intent(inout) :: this
82  character(len=*), intent(in) :: field_name
83  class(sll_t_cartesian_mesh_1d), pointer :: mesh
84  sll_int32, intent(in) :: data_position
85  procedure(scalar_function_1d), optional :: init_function
86 
87  character(len=*), parameter :: this_sub_name = 'initialize_scalar_field_1d'
88 
89  sll_int32 :: ierr
90  sll_int32 :: num_cells1
91  sll_int32 :: num_pts1
92  sll_int32 :: i1
93  sll_real64 :: eta1
94  sll_real64 :: delta1
95 
96  sll_warning(this_sub_name, 'This scalar field 1d object is deprecated.')
97 
98  this%mesh => mesh
99  this%name = trim(field_name)
100  num_cells1 = mesh%num_cells
101  num_pts1 = mesh%num_cells + 1
102 
103  this%data_position = data_position
104  if (data_position == sll_p_node_centered_field) then
105  sll_allocate(this%data(num_pts1), ierr)
106  if (present(init_function)) then
107  do i1 = 1, num_pts1
108  this%data(i1) = init_function(mesh%eta_min + (i1 - 1)*mesh%delta_eta)
109  end do
110  else
111  this%data = 0.0_f64 ! initialize to zero
112  end if
113  else if (data_position == sll_p_cell_centered_field) then
114  sll_allocate(this%data(num_cells1), ierr)
115  if (present(init_function)) then
116  delta1 = 1.0_f64/real(num_cells1, f64)
117  eta1 = 0.5_f64*delta1
118  do i1 = 1, num_cells1
119  this%data(i1) = init_function(eta1)
120  eta1 = eta1 + delta1
121  end do
122  else
123  this%data = 0.0_f64 ! initialize to zero
124  end if
125  end if
126  end subroutine initialize_scalar_field_1d
127 
128  ! need to do something about deallocating the field proper, when allocated
129  ! in the heap...
130  subroutine delete_scalar_field_1d(this)
131  type(scalar_field_1d), pointer :: this
132  sll_int32 :: ierr
133  nullify (this%mesh)
134  sll_deallocate(this%data, ierr)
135  end subroutine delete_scalar_field_1d
136 
137  subroutine write_scalar_field_1d( &
138  scalar_field, &
139  multiply_by_jacobian, &
140  output_file_name, &
141  output_format)
142  class(scalar_field_1d) :: scalar_field
143  !sll_real64, dimension(:), pointer :: x1_array
144  logical, optional :: multiply_by_jacobian
145  sll_int32, optional :: output_format
146  character(len=*), optional :: output_file_name
147  character(len=64) :: file_name
148 
149  if (present(multiply_by_jacobian)) then
150  print *, 'multiply_by_jacobian option is not implemented'
151  end if
152 
153  if (present(output_format)) then
154  print *, 'There is just gnuplot format available'
155  end if
156  if (present(output_file_name)) then
157  file_name = output_file_name
158  else
159  file_name = scalar_field%name
160  end if
161 
162  call sll_s_gnuplot_write(scalar_field%data, file_name, 1)
163  end subroutine write_scalar_field_1d
164 
165 end module sll_m_scalar_field_1d_old
Cartesian mesh basic types.
Implements the functions to write data file plotable by GNUplot.
subroutine, public sll_s_gnuplot_write(array, array_name, iplot)
Write an array to display with gnuplot.
Implements the geometry and mesh descriptor types.
subroutine initialize_scalar_field_1d(this, field_name, mesh, data_position, init_function)
subroutine write_scalar_field_1d(scalar_field, multiply_by_jacobian, output_file_name, output_format)
Module to select the kind parameter.
    Report Typos and Errors