Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_advection_2d_tensor_product.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 
23 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
24 #include "sll_memory.h"
25 #include "sll_working_precision.h"
26 
27  use sll_m_advection_1d_base, only: &
29 
30  use sll_m_advection_2d_base, only: &
32 
33  implicit none
34 
35  public :: &
38 
39  private
40 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
41 
43 
44  class(sll_c_advector_1d), pointer :: advect_x1
45  class(sll_c_advector_1d), pointer :: advect_x2
46  sll_int32 :: npts1
47  sll_int32 :: npts2
48  sll_real64, dimension(:), pointer :: buf1d
49 
50  contains
51 
52  procedure, pass(adv) :: init => initialize_advector_2d_tensor_product
53  procedure, pass(adv) :: advect_2d => tensor_product_advect_2d
54 
56 
57 contains
58 
59  function sll_f_new_advector_2d_tensor_product(advect_x1, advect_x2, &
60  npts1, npts2) result(adv)
61 
62  type(sll_t_advector_2d_tensor_product), pointer :: adv
63  class(sll_c_advector_1d), pointer :: advect_x1
64  class(sll_c_advector_1d), pointer :: advect_x2
65  sll_int32, intent(in) :: npts1
66  sll_int32, intent(in) :: npts2
67 
68  sll_int32 :: ierr
69 
70  sll_allocate(adv, ierr)
71 
72  call adv%init(advect_x1, advect_x2, npts1, npts2)
73 
75 
77  adv, &
78  advect_x1, &
79  advect_x2, &
80  npts1, &
81  npts2)
82 
83  class(sll_t_advector_2d_tensor_product), intent(inout) :: adv
84  class(sll_c_advector_1d), pointer :: advect_x1
85  class(sll_c_advector_1d), pointer :: advect_x2
86  sll_int32, intent(in) :: npts1
87  sll_int32, intent(in) :: npts2
88 
89  sll_int32 :: ierr
90 
91  adv%advect_x1 => advect_x1
92  adv%advect_x2 => advect_x2
93 
94  adv%npts1 = npts1
95  adv%npts2 = npts2
96 
97  sll_allocate(adv%buf1d(max(npts1, npts2)), ierr)
98 
100 
101  subroutine tensor_product_advect_2d(adv, A1, A2, dt, input, output)
102 
104  sll_real64, dimension(:, :), intent(in) :: a1
105  sll_real64, dimension(:, :), intent(in) :: a2
106  sll_real64, intent(in) :: dt
107  sll_real64, dimension(:, :), intent(in) :: input
108  sll_real64, dimension(:, :), intent(out) :: output
109 
110  sll_int32 :: i1
111  sll_int32 :: i2
112  sll_int32 :: npts1
113  sll_int32 :: npts2
114 
115  npts1 = adv%npts1
116  npts2 = adv%npts2
117 
118  do i2 = 1, npts2
119  adv%buf1d(1:npts1) = input(1:npts1, i2)
120  call adv%advect_x1%advect_1d( &
121  a1(1:npts1, i2), &
122  0.5_f64*dt, &
123  adv%buf1d(1:npts1), &
124  output(1:npts1, i2))
125  end do
126 
127  do i1 = 1, npts1
128  adv%buf1d(1:npts2) = output(i1, 1:npts2)
129  call adv%advect_x2%advect_1d( &
130  a2(i1, 1:npts2), &
131  dt, &
132  adv%buf1d(1:npts2), &
133  output(i1, 1:npts2))
134  end do
135 
136  do i2 = 1, npts2
137  adv%buf1d(1:npts1) = output(1:npts1, i2)
138  call adv%advect_x1%advect_1d( &
139  a1(1:npts1, i2), &
140  0.5_f64*dt, &
141  adv%buf1d(1:npts1), &
142  output(1:npts1, i2))
143  end do
144 
145  end subroutine tensor_product_advect_2d
146 
Abstract class for advection.
subroutine tensor_product_advect_2d(adv, A1, A2, dt, input, output)
type(sll_t_advector_2d_tensor_product) function, pointer, public sll_f_new_advector_2d_tensor_product(advect_x1, advect_x2, npts1, npts2)
subroutine initialize_advector_2d_tensor_product(adv, advect_x1, advect_x2, npts1, npts2)
    Report Typos and Errors