Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_advection_1d_BSL.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
18 ! in development; should be at least cubic splines
19 ! attached with computation of characteristics
20 
24 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
25 #include "sll_memory.h"
26 #include "sll_working_precision.h"
27 #include "sll_errors.h"
28 
32 
33  implicit none
34 
36 
37  private
38 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
39 
41 
42  class(sll_c_interpolator_1d), pointer :: interp
43  class(sll_c_characteristics_1d_base), pointer :: charac
44  sll_real64, dimension(:), pointer :: eta_coords
45  sll_real64, dimension(:), pointer :: charac_feet
46  sll_int32 :: npts
47 
48  contains
49 
50  procedure, pass(adv) :: init => initialize_advector_1d_bsl
51  procedure, pass(adv) :: advect_1d => bsl_advect_1d
52  procedure, pass(adv) :: advect_1d_constant => bsl_advect_1d_constant
53  procedure, pass(adv) :: delete => delete_bsl_1d_adv
54 
55  end type sll_t_advector_1d_bsl
56 
57 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
58 contains
59 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 
61  function sll_f_new_advector_1d_bsl(interp, charac, npts, eta_min, eta_max, &
62  eta_coords) result(adv)
63 
64  type(sll_t_advector_1d_bsl), pointer :: adv
65  class(sll_c_interpolator_1d), pointer :: interp
66  class(sll_c_characteristics_1d_base), pointer :: charac
67 
68  sll_int32, intent(in) :: npts
69  sll_real64, intent(in), optional :: eta_min
70  sll_real64, intent(in), optional :: eta_max
71  sll_real64, pointer, optional :: eta_coords(:)
72 
73  sll_int32 :: ierr
74 
75  sll_allocate(adv, ierr)
76 
77  call adv%init(interp, charac, npts, eta_min, eta_max, eta_coords)
78 
79  end function sll_f_new_advector_1d_bsl
80 
81  subroutine initialize_advector_1d_bsl(adv, interp, charac, npts, &
82  eta_min, eta_max, eta_coords)
83 
84  class(sll_t_advector_1d_bsl), intent(inout) :: adv
85  class(sll_c_interpolator_1d), target :: interp
86  class(sll_c_characteristics_1d_base), target :: charac
87  sll_int32, intent(in) :: npts
88 
89  sll_real64, intent(in), optional :: eta_min
90  sll_real64, intent(in), optional :: eta_max
91  sll_real64, dimension(:), pointer, optional :: eta_coords
92 
93  sll_int32 :: ierr
94  sll_int32 :: i
95  sll_real64 :: delta_eta
96 
97  adv%npts = npts
98  adv%interp => interp
99  adv%charac => charac
100 
101  sll_allocate(adv%eta_coords(npts), ierr)
102  sll_allocate(adv%charac_feet(npts), ierr)
103 
104  if (present(eta_min) .and. present(eta_max)) then
105  if (present(eta_coords)) then
106  sll_error('initialize_advector_1d_bsl', 'provide either eta_coords or eta_min and eta_max and not both')
107  else
108  delta_eta = (eta_max - eta_min)/real(npts - 1, f64)
109  do i = 1, npts
110  adv%eta_coords(i) = eta_min + real(i - 1, f64)*delta_eta
111  end do
112  end if
113  else if (present(eta_coords)) then
114  if (size(eta_coords) < npts) then
115  sll_error('initialize_advector_1d_bsl', 'bad size for eta_coords in initialize_sll_t_advector_1d_bsl')
116  else
117  adv%eta_coords(1:npts) = eta_coords(1:npts)
118  end if
119  else
120  sll_warning('initialize_advector_1d_bsl', 'we assume eta_min = 0._f64 eta_max = 1._f64')
121  delta_eta = 1._f64/real(npts - 1, f64)
122  do i = 1, npts
123  adv%eta_coords(i) = real(i - 1, f64)*delta_eta
124  end do
125  end if
126 
127  end subroutine initialize_advector_1d_bsl
128 
129  subroutine bsl_advect_1d(adv, A, dt, input, output)
130 
131  class(sll_t_advector_1d_bsl) :: adv
132  sll_real64, dimension(:), intent(in) :: a
133  sll_real64, intent(in) :: dt
134  sll_real64, dimension(:), intent(in) :: input
135  sll_real64, dimension(:), intent(out) :: output
136 
137  call adv%charac%compute_characteristics(a, dt, adv%eta_coords, adv%charac_feet)
138 
139  call adv%interp%interpolate_array(adv%npts, input, adv%charac_feet, output)
140 
141  end subroutine bsl_advect_1d
142 
143  subroutine bsl_advect_1d_constant(adv, A, dt, input, output)
144 
145  class(sll_t_advector_1d_bsl) :: adv
146  sll_real64, intent(in) :: a
147  sll_real64, intent(in) :: dt
148  sll_real64, dimension(:), intent(in) :: input
149  sll_real64, dimension(:), intent(out) :: output
150  sll_real64, dimension(:), allocatable :: a1
151  sll_int32 :: ierr
152 
153  sll_allocate(a1(adv%npts), ierr)
154 
155  a1 = a
156 
157  call adv%charac%compute_characteristics(a1, dt, adv%eta_coords, adv%charac_feet)
158 
159  call adv%interp%interpolate_array(adv%npts, input, adv%charac_feet, output)
160 
161  sll_deallocate_array(a1, ierr)
162 
163  end subroutine bsl_advect_1d_constant
164 
165  subroutine delete_bsl_1d_adv(adv)
166 
167  class(sll_t_advector_1d_bsl), intent(inout) :: adv
168 
169  sll_int32 :: ierr
170 
171  sll_deallocate(adv%eta_coords, ierr)
172  sll_deallocate(adv%charac_feet, ierr)
173 
174  end subroutine delete_bsl_1d_adv
175 
176 end module sll_m_advection_1d_bsl
Abstract class for advection.
backward semi-lagrangian 1d advection
subroutine initialize_advector_1d_bsl(adv, interp, charac, npts, eta_min, eta_max, eta_coords)
subroutine bsl_advect_1d(adv, A, dt, input, output)
subroutine bsl_advect_1d_constant(adv, A, dt, input, output)
type(sll_t_advector_1d_bsl) function, pointer, public sll_f_new_advector_1d_bsl(interp, charac, npts, eta_min, eta_max, eta_coords)
Abstract class for characteristic derived type.
Module for 1D interpolation and reconstruction.
Abstract class for 1D interpolation and reconstruction.
    Report Typos and Errors