Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_advection_2d_oblic.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 
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 #include "sll_memory.h"
24 #include "sll_working_precision.h"
25 
26  use sll_m_advection_1d_base, only: &
28 
29  use sll_m_lagrange_interpolation, only: &
31 
32  implicit none
33 
34  public :: &
38 
39  private
40 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
41 
43  sll_int32 :: nc_x1
44  class(sll_c_advector_1d), pointer :: adv_x1
45  sll_int32 :: nc_x2
46  sll_real64 :: x2_min
47  sll_real64 :: x2_max
48  sll_int32 :: stencil_r
49  sll_int32 :: stencil_s
50  sll_real64, dimension(:), pointer :: xx
51  sll_real64, dimension(:, :), pointer :: buf
52 
54 
55 contains
57  Nc_x1, &
58  adv_x1, &
59  Nc_x2, &
60  x2_min, &
61  x2_max, &
62  stencil_r, &
63  stencil_s) &
64  result(adv)
65  type(sll_t_oblic_2d_advector), pointer :: adv
66  sll_int32, intent(in) :: nc_x1
67  class(sll_c_advector_1d), pointer :: adv_x1
68  sll_int32, intent(in) :: nc_x2
69  sll_real64, intent(in) :: x2_min
70  sll_real64, intent(in) :: x2_max
71  sll_int32, intent(in) :: stencil_r
72  sll_int32, intent(in) :: stencil_s
73  !local variables
74  sll_int32 :: ierr
75 
76  sll_allocate(adv, ierr)
77 
79  adv, &
80  nc_x1, &
81  adv_x1, &
82  nc_x2, &
83  x2_min, &
84  x2_max, &
85  stencil_r, &
86  stencil_s)
87  end function sll_f_new_oblic_2d_advector
88 
90  adv, &
91  Nc_x1, &
92  adv_x1, &
93  Nc_x2, &
94  x2_min, &
95  x2_max, &
96  stencil_r, &
97  stencil_s)
98  type(sll_t_oblic_2d_advector), intent(inout) :: adv
99  sll_int32, intent(in) :: nc_x1
100  class(sll_c_advector_1d), pointer :: adv_x1
101  sll_int32, intent(in) :: nc_x2
102  sll_real64, intent(in) :: x2_min
103  sll_real64, intent(in) :: x2_max
104  sll_int32, intent(in) :: stencil_r
105  sll_int32, intent(in) :: stencil_s
106  !local variables
107  sll_int32 :: ierr
108  sll_int32 :: r
109  sll_int32 :: s
110  sll_int32 :: i
111 
112  adv%Nc_x1 = nc_x1
113  adv%adv_x1 => adv_x1
114  adv%Nc_x2 = nc_x2
115  adv%x2_min = x2_min
116  adv%x2_max = x2_max
117  adv%stencil_r = stencil_r
118  adv%stencil_s = stencil_s
119 
120  r = adv%stencil_r
121  s = adv%stencil_s
122 
123  sll_allocate(adv%xx(r:s), ierr)
124  sll_allocate(adv%buf(r:s, nc_x1 + 1), ierr)
125 
126  do i = r, s
127  adv%xx(i) = real(i, f64)
128  end do
129 
130  end subroutine initialize_oblic_2d_advector
131 
150 
152  adv, &
153  !iota, &
154  A1, &
155  A2, &
156  dt, &
157  input, &
158  output)
159  type(sll_t_oblic_2d_advector) :: adv
160  !sll_real64, intent(in) :: iota
161  sll_real64, intent(in) :: a1
162  sll_real64, intent(in) :: a2
163  sll_real64, intent(in) :: dt
164  sll_real64, dimension(:, :), intent(in) :: input
165  sll_real64, dimension(:, :), intent(out) :: output
166  !local variables
167  sll_int32 :: i1
168  sll_int32 :: i2
169  sll_int32 :: i2_loc
170  sll_int32 :: i0
171  sll_real64 :: alpha
172  sll_int32 :: d
173  sll_int32 :: r
174  sll_int32 :: s
175  sll_real64 :: dt_loc
176  sll_int32 :: ell
177  sll_real64, dimension(:, :), pointer :: buf
178  sll_real64, dimension(:), pointer :: xx
179  class(sll_c_advector_1d), pointer :: adv_x1
180  sll_real64 :: delta_x2
181  sll_int32 :: nc_x1
182  sll_int32 :: nc_x2
183 
184  buf => adv%buf
185  xx => adv%xx
186  adv_x1 => adv%adv_x1
187  nc_x1 = adv%Nc_x1
188  nc_x2 = adv%Nc_x2
189  delta_x2 = (adv%x2_max - adv%x2_min)/real(adv%Nc_x2, f64)
190  r = adv%stencil_r
191  s = adv%stencil_s
192  alpha = a2*dt/delta_x2
193  i0 = floor(alpha)
194  alpha = alpha - i0
195  d = s - r
196 
197  do i2 = 1, nc_x2 + 1
198  !choose several dt_loc so that advection in x2 is exact
199  do ell = r, s
200  dt_loc = real(ell + i0, f64)*delta_x2/a2
201  i2_loc = modulo(i2 - ell - i0 - 1, nc_x2) + 1
202  call adv_x1%advect_1d_constant( &
203  a1, &
204  dt_loc, &
205  input(1:nc_x1 + 1, i2_loc), &
206  buf(ell, 1:nc_x1 + 1))
207  end do
208  !interpolate between these values
209  do i1 = 1, nc_x1 + 1
210  output(i1, i2) = sll_f_lagrange_interpolate(alpha, d, xx, buf(r:s, i1))
211  end do
212  end do
213 
214  end subroutine sll_s_oblic_advect_2d_constant
215 
216  subroutine oblic_advect_2d( &
217  adv, &
218  !iota, &
219  A1, &
220  A2, &
221  dt, &
222  input, &
223  output)
224  type(sll_t_oblic_2d_advector) :: adv
225  !sll_real64, intent(in) :: iota
226  sll_real64, dimension(:, :), intent(in) :: a1
227  sll_real64, dimension(:, :), intent(in) :: a2
228  sll_real64, intent(in) :: dt
229  sll_real64, dimension(:, :), intent(in) :: input
230  sll_real64, dimension(:, :), intent(out) :: output
231 
232  print *, '#oblic_advect_2d not implemented for the moment'
233  print *, size(input, 1), size(input, 2), storage_size(adv)
234  output = 0.0_f64 + a1*a2*dt
235 
236  end subroutine oblic_advect_2d
237 
238 end module sll_m_advection_2d_oblic
Abstract class for advection.
subroutine, public sll_s_oblic_advect_2d_constant(adv, A1, A2, dt, input, output)
solves \partial_t f + \nabla A \cdot f = 0, \ A = (A1,A2) for time step dt interpolation in aligned a...
type(sll_t_oblic_2d_advector) function, pointer, public sll_f_new_oblic_2d_advector(Nc_x1, adv_x1, Nc_x2, x2_min, x2_max, stencil_r, stencil_s)
subroutine oblic_advect_2d(adv, A1, A2, dt, input, output)
subroutine initialize_oblic_2d_advector(adv, Nc_x1, adv_x1, Nc_x2, x2_min, x2_max, stencil_r, stencil_s)
real(kind=f64) function, public sll_f_lagrange_interpolate(x, degree, xi, yi)
    Report Typos and Errors