Loading [MathJax]/extensions/TeX/AMSsymbols.js
Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
All Classes Namespaces Files Functions Variables Macros Modules Pages
sll_m_advection_1d_periodic.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 ! for the moment mimic of sll_m_periodic_interpolator_1d.F90
19 
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 #include "sll_assert.h"
24 #include "sll_memory.h"
25 #include "sll_working_precision.h"
26 
28 
29  use sll_m_periodic_interp, only: &
33 
34  implicit none
35 
38 
39  private
40 
42 
43  sll_int32 :: num_cells
44  sll_real64 :: xmin
45  sll_real64 :: xmax
46  type(sll_t_periodic_interp_work) :: per_interp
47 
48  contains
49 
50  procedure, pass(adv) :: init => initialize_periodic_1d_advector
51  procedure, pass(adv) :: advect_1d_constant => periodic_advect_1d_constant
52  procedure, pass(adv) :: advect_1d => periodic_advect_1d_fake
53  procedure, pass(adv) :: delete => delete_periodic_1d_advector
55 
56 contains
57 
59  num_cells, &
60  xmin, &
61  xmax, &
62  type, &
63  order) &
64  result(adv)
65  type(sll_t_advector_1d_periodic), pointer :: adv
66  sll_int32, intent(in) :: num_cells
67  sll_real64, intent(in) :: xmin
68  sll_real64, intent(in) :: xmax
69  sll_int32, intent(in) :: type
70  sll_int32, intent(in) :: order
71  sll_int32 :: ierr
72 
73  sll_allocate(adv, ierr)
74  call initialize_periodic_1d_advector(adv, num_cells, xmin, xmax, type, order)
75 
77 
79  adv, &
80  num_cells, &
81  xmin, &
82  xmax, &
83  type, &
84  order)
85 
86  class(sll_t_advector_1d_periodic) :: adv
87  sll_int32, intent(in) :: num_cells
88  sll_real64, intent(in) :: xmin
89  sll_real64, intent(in) :: xmax
90  sll_int32, intent(in) :: type
91  sll_int32, intent(in) :: order
92 
93  call sll_s_periodic_interp_init(adv%per_interp, num_cells, type, order)
94 
95  adv%num_cells = num_cells
96  adv%xmin = xmin
97  adv%xmax = xmax
98 
99  end subroutine initialize_periodic_1d_advector
100 
101  subroutine periodic_advect_1d_constant(adv, A, dt, input, output)
102 
103  class(sll_t_advector_1d_periodic) :: adv
104  sll_real64, intent(in) :: a
105  sll_real64, intent(in) :: dt
106  sll_real64, dimension(:), intent(in) :: input
107  sll_real64, dimension(:), intent(out) :: output
108 
109  sll_real64 :: shift
110  sll_real64 :: xmin
111  sll_real64 :: xmax
112  sll_int32 :: num_cells
113 
114  num_cells = adv%num_cells
115  xmin = adv%xmin
116  xmax = adv%xmax
117  shift = a*dt/(xmax - xmin)*real(num_cells, f64)
118 
119  call sll_s_periodic_interp( &
120  adv%per_interp, &
121  output(1:num_cells), &
122  input(1:num_cells), &
123  shift)
124 
125  ! complete by periodicity
126  if (size(output) > num_cells) then
127  output(num_cells + 1) = output(1)
128  end if
129 
130  end subroutine periodic_advect_1d_constant
131 
132  subroutine periodic_advect_1d_fake(adv, A, dt, input, output)
133 
134  class(sll_t_advector_1d_periodic) :: adv
135  sll_real64, dimension(:), intent(in) :: a
136  sll_real64, intent(in) :: dt
137  sll_real64, dimension(:), intent(in) :: input
138  sll_real64, dimension(:), intent(out) :: output
139 
140  print *, '#periodic_advect_1d_fake'
141  print *, '#not implemented'
142  print *, maxval(a)
143  print *, dt
144  print *, maxval(input)
145  output = 0._f64
146  print *, adv%num_cells
147  stop
148 
149  end subroutine periodic_advect_1d_fake
150 
152  class(sll_t_advector_1d_periodic), intent(inout) :: adv
153  sll_assert(storage_size(adv) > 0)
154  end subroutine delete_periodic_1d_advector
155 
Abstract class for advection.
type(sll_t_advector_1d_periodic) function, pointer, public sll_f_new_periodic_1d_advector(num_cells, xmin, xmax, type, order)
subroutine periodic_advect_1d_fake(adv, A, dt, input, output)
subroutine periodic_advect_1d_constant(adv, A, dt, input, output)
subroutine initialize_periodic_1d_advector(adv, num_cells, xmin, xmax, type, order)
subroutine, public sll_s_periodic_interp_init(this, N, interpolator, order)
subroutine, public sll_s_periodic_interp(this, u_out, u, alpha)
    Report Typos and Errors