Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_singular_mapping_analytic.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_assert.h"
4 
6 
8 
10 
11  use sll_m_hdf5_io_serial, only: &
16 
17  implicit none
18 
20 
21  private
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 
25  integer, parameter :: wp = f64
26 
28  ! (may contain common components/methods)
30 
31  contains
32 
33  ! Deferred procedures
34  procedure(i_fun_jmat_comp), deferred :: jmat_comp
35 
36  ! Non-deferred procedures
37  procedure :: store_data => s_singular_mapping_analytic__store_data
38 
40 
41  ! Interfaces for deferred procedures
42  abstract interface
43 
44  sll_pure function i_fun_jmat_comp(self, eta) result(jmat_comp)
46  class(sll_c_singular_mapping_analytic), intent(in) :: self
47  real(wp), intent(in) :: eta(2)
48  real(wp) :: jmat_comp(2, 2)
49  end function i_fun_jmat_comp
50 
51  end interface
52 
53 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
54 contains
55 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
56 
57  subroutine s_singular_mapping_analytic__store_data(self, n1, n2, file_id)
58  class(sll_c_singular_mapping_analytic), intent(in) :: self
59  integer, intent(in) :: n1
60  integer, intent(in) :: n2
61  type(sll_t_hdf5_ser_handle), intent(in) :: file_id
62 
63  integer :: i1, i2
64  real(wp) :: eta(2), x(2)
65 
66  real(wp), allocatable :: x1(:, :), x2(:, :), jacobian(:, :)
67 
68  ! For hdf5 I/O
69  integer :: error
70 
71  ! Allocate physical mesh
72  allocate (x1(n1, n2 + 1)) ! repeated point along eta2
73  allocate (x2(n1, n2 + 1)) ! repeated point along eta2
74 
75  ! Allocate Jacobian determinant
76  allocate (jacobian(n1, n2 + 1)) ! repeated point along eta2
77 
78  ! Compute physical mesh and Jacobian determinant
79  do i2 = 1, n2 + 1 ! repeated point along eta2
80  do i1 = 1, n1
81  eta(1) = real(i1 - 1, wp)/real(n1 - 1, wp)
82  eta(2) = real(i2 - 1, wp)*sll_p_twopi/real(n2, wp)
83  x(:) = self%eval(eta)
84  x1(i1, i2) = x(1)
85  x2(i1, i2) = x(2)
86  jacobian(i1, i2) = self%jdet(eta)
87  end do
88  end do
89 
90  ! Store physical mesh
91  call sll_o_hdf5_ser_write_array(file_id, x1, "/x1", error)
92  call sll_o_hdf5_ser_write_array(file_id, x2, "/x2", error)
93 
94  ! Store Jacobian determinant
95  call sll_o_hdf5_ser_write_array(file_id, jacobian, "/jacobian", error)
96 
98 
Write nD array of double precision floats or integers into HDF5 file.
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_twopi
Implements the functions to write hdf5 file to store heavy data.
subroutine, public sll_s_hdf5_ser_file_create(filename, handle, error)
Create new HDF5 file.
subroutine, public sll_s_hdf5_ser_file_close(handle, error)
Close existing HDF5 file.
subroutine s_singular_mapping_analytic__store_data(self, n1, n2, file_id)
integer, parameter wp
Working precision.
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
    Report Typos and Errors