Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_poisson_2d_fem_sps_stencil_new_assembler.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 
5 
7 
8  implicit none
9 
11 
12  private
13 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
14 
15  ! Working precision
16  integer, parameter :: wp = f64
17 
19 
20  integer :: n1
21  integer :: n2
22  integer :: p1
23  integer :: p2
24 
25  class(sll_t_poisson_2d_fem_sps_weak_form), pointer :: weak_form
26 
27  contains
28 
32 
34 
35 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
36 contains
37 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38 
39  ! Initializer
40  subroutine s_poisson_2d_fem_sps_stencil_new_assembler__init(self, n1, n2, p1, p2, weak_form)
41  class(sll_t_poisson_2d_fem_sps_stencil_new_assembler), intent(inout) :: self
42  integer, intent(in) :: n1
43  integer, intent(in) :: n2
44  integer, intent(in) :: p1
45  integer, intent(in) :: p2
46  class(sll_t_poisson_2d_fem_sps_weak_form), target, intent(in) :: weak_form
47 
48  self%n1 = n1
49  self%n2 = n2
50  self%p1 = p1
51  self%p2 = p2
52 
53  self%weak_form => weak_form
54 
55  end subroutine
56 
57  ! Add element in stiffness and mass matrices
59  self, &
60  k1, &
61  k2, &
62  data_1d_eta1, &
63  data_1d_eta2, &
64  int_volume, &
65  inv_metric, &
66  A, &
67  M)
68  class(sll_t_poisson_2d_fem_sps_stencil_new_assembler), intent(in) :: self
69  integer, intent(in) :: k1
70  integer, intent(in) :: k2
71  real(wp), intent(in) :: data_1d_eta1(:, :, :, :)
72  real(wp), intent(in) :: data_1d_eta2(:, :, :, :)
73  real(wp), intent(in) :: int_volume(:, :, :, :)
74  real(wp), intent(in) :: inv_metric(:, :, :, :, :, :)
75  real(wp), intent(inout) :: A(-self%p1:, -self%p2:, :, :)
76  real(wp), intent(inout) :: M(-self%p1:, -self%p2:, :, :)
77 
78  integer :: i1, i2, j1, j2, s1, s2
79  real(wp) :: Ael, Mel
80  integer :: il1, il2, jl1, jl2
81 
82  associate(n1 => self%n1, &
83  n2 => self%n2, &
84  p1 => self%p1, &
85  p2 => self%p2)
86 
87  ! Cycle over basis functions in 2D trial space (same as test space)
88  do jl2 = 1, p2 + 1
89 
90  j2 = modulo(k2 + jl2 - 2, n2) + 1
91 
92  do jl1 = 1, p1 + 1
93 
94  j1 = k1 + jl1 - 1
95 
96  ! Cycle over basis functions in 2D test space
97  do il2 = 1, p2 + 1
98 
99  i2 = modulo(k2 + il2 - 2, n2) + 1
100 
101  do il1 = 1, p1 + 1
102 
103  i1 = k1 + il1 - 1
104 
105  call self%weak_form%element_mat( &
106  test_values_and_derivs_eta1=data_1d_eta1(:, :, il1, k1), &
107  test_values_and_derivs_eta2=data_1d_eta2(:, :, il2, k2), &
108  trial_values_and_derivs_eta1=data_1d_eta1(:, :, jl1, k1), &
109  trial_values_and_derivs_eta2=data_1d_eta2(:, :, jl2, k2), &
110  int_volume=int_volume(:, :, k1, k2), &
111  inv_metric=inv_metric(:, :, k1, k2, :, :), &
112  aij=ael, &
113  mij=mel)
114 
115  s1 = j1 - i1
116  s2 = modulo(j2 - i2 + p2, n2) - p2
117 
118  a(s1, s2, i1, i2) = a(s1, s2, i1, i2) + ael
119  m(s1, s2, i1, i2) = m(s1, s2, i1, i2) + mel
120 
121  end do
122  end do ! End cycle over basis functions in 2D test space
123 
124  end do
125  end do ! End cycle over basis functions in 2D trial space
126 
127  end associate
128 
130 
131  ! Add element in stiffness and mass matrices
133  self, &
134  k1, &
135  k2, &
136  data_1d_eta1, &
137  data_1d_eta2, &
138  data_2d_rhs, &
139  int_volume, &
140  bs)
141  class(sll_t_poisson_2d_fem_sps_stencil_new_assembler), intent(in) :: self
142  integer, intent(in) :: k1
143  integer, intent(in) :: k2
144  real(wp), intent(in) :: data_1d_eta1(:, :, :, :)
145  real(wp), intent(in) :: data_1d_eta2(:, :, :, :)
146  real(wp), intent(in) :: data_2d_rhs(:, :, :, :)
147  real(wp), intent(in) :: int_volume(:, :, :, :)
148  real(wp), intent(inout) :: bs(:, :)
149 
150  integer :: i1, i2, p1, p2
151  real(wp) :: bi
152  integer :: il1, il2
153 
154  ! Extract degree of B-splines
155  p1 = size(data_1d_eta1, 3) - 1
156  p2 = size(data_1d_eta2, 3) - 1
157 
158  associate(n1 => self%n1, n2 => self%n2)
159 
160  ! Cycle over basis functions in 2D test space
161  do il2 = 1, p2 + 1
162 
163  i2 = modulo(il2 - 2 + k2, n2) + 1
164 
165  do il1 = 1, p1 + 1
166 
167  i1 = il1 - 1 + k1
168 
169  call self%weak_form%element_rhs( &
170  test_values_and_derivs_eta1=data_1d_eta1(:, :, il1, k1), &
171  test_values_and_derivs_eta2=data_1d_eta2(:, :, il2, k2), &
172  data_2d_rhs=data_2d_rhs(:, :, k1, k2), &
173  int_volume=int_volume(:, :, k1, k2), &
174  bi=bi)
175 
176  bs(i1, i2) = bs(i1, i2) + bi
177 
178  end do
179  end do ! End cycle over basis functions in 2D test space
180 
181  end associate
182 
184 
subroutine s_poisson_2d_fem_sps_stencil_new_assembler__init(self, n1, n2, p1, p2, weak_form)
subroutine s_poisson_2d_fem_sps_stencil_new_assembler__add_element_mat(self, k1, k2, data_1d_eta1, data_1d_eta2, int_volume, inv_metric, A, M)
subroutine s_poisson_2d_fem_sps_stencil_new_assembler__add_element_rhs(self, k1, k2, data_1d_eta1, data_1d_eta2, data_2d_rhs, int_volume, bs)
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