Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_linear_operator_matrix_dense_to_stencil.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_assert.h"
4 #include "sll_errors.h"
5 
7 
9 
11 
13 
15 
16  implicit none
17 
19 
20  private
21 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
22 
23  ! Working precision
24  integer, parameter :: wp = f64
25 
27 
28  real(wp), allocatable :: a(:, :)
29 
30  integer :: s1
31  integer :: s2
32 
33  contains
34 
41 
43 
44 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
45 contains
46 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 
48  ! Initialize linear operator
50  class(sll_t_linear_operator_matrix_dense_to_stencil), intent(inout) :: self
51  integer, intent(in) :: s1
52  integer, intent(in) :: s2
53 
54  allocate (self%A(s1, s2))
55 
56  self%s1 = s1
57  self%s2 = s2
58 
60 
61  ! Get shape of linear operator
63  class(sll_t_linear_operator_matrix_dense_to_stencil), intent(in) :: self
64  integer :: s(2)
65 
66  s = shape(self%A)
67 
69 
70  ! Implement y=Ax, with A dense matrix (2D), x dense vector (1D) and y stencil vector (2D)
72  class(sll_t_linear_operator_matrix_dense_to_stencil), intent(in) :: self
73  class(sll_c_vector_space), intent(in) :: x ! dense
74  class(sll_c_vector_space), intent(inout) :: y ! stencil
75 
76  integer :: j1, j2, j, k
77 
78  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_dense_to_stencil % dot"
79  character(len=64) :: err_msg
80 
81  select type (x)
82 
84 
85  ! Check dimensions
86  sll_assert(self%s2 == size(x%array))
87 
88  select type (y)
89 
91 
92  associate(p1 => -lbound(y%array, 1) + 1, &
93  p2 => -lbound(y%array, 2) + 1)
94 
95  associate(ny1 => ubound(y%array, 1) - p1, &
96  ny2 => ubound(y%array, 2) - p2)
97 
98  ! Check dimensions
99  sll_assert(self%s1 == ny1*ny2)
100 
101  y%array = 0.0_wp
102  do j2 = 1, ny2
103  do j1 = 1, ny1
104  j = (j1 - 1)*ny2 + j2
105  do k = 1, size(x%array)
106  y%array(j1, j2) = y%array(j1, j2) + self%A(j, k)*x%array(k)
107  end do
108  end do
109  end do
110 
111  ! Update buffer regions
112  y%array(1 - p1:0, :) = y%array(ny1 - p1 + 1:ny1, :)
113  y%array(ny1 + 1:ny1 + p1, :) = y%array(1:p1, :)
114  y%array(:, 1 - p2:0) = y%array(:, ny2 - p2 + 1:ny2)
115  y%array(:, ny2 + 1:ny2 + p2) = y%array(:, 1:p2)
116 
117  end associate
118 
119  end associate
120 
121  class default
122  err_msg = "y must be of type sll_t_vector_space_real_array_2d"
123  sll_error(this_sub_name, err_msg)
124 
125  end select
126 
127  class default
128  err_msg = "x must be of type sll_t_vector_space_real_array_1d"
129  sll_error(this_sub_name, err_msg)
130 
131  end select
132 
134 
135  ! Implement y=y+Ax, with A dense matrix (2D), x dense vector (1D) and y stencil vector (2D)
137  class(sll_t_linear_operator_matrix_dense_to_stencil), intent(in) :: self
138  class(sll_c_vector_space), intent(in) :: x ! dense
139  class(sll_c_vector_space), intent(inout) :: y ! stencil
140 
141  integer :: j1, j2, j, k
142 
143  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_dense_to_stencil % dot_incr"
144  character(len=64) :: err_msg
145 
146  select type (x)
147 
148  type is (sll_t_vector_space_real_array_1d)
149 
150  ! Check dimensions
151  sll_assert(self%s2 == size(x%array))
152 
153  select type (y)
154 
155  type is (sll_t_vector_space_real_array_2d)
156 
157  associate(p1 => -lbound(y%array, 1) + 1, &
158  p2 => -lbound(y%array, 2) + 1)
159 
160  associate(ny1 => ubound(y%array, 1) - p1, &
161  ny2 => ubound(y%array, 2) - p2)
162 
163  ! Check dimensions
164  sll_assert(self%s1 == ny1*ny2)
165 
166  do j2 = 1, ny2
167  do j1 = 1, ny1
168  j = (j1 - 1)*ny2 + j2
169  do k = 1, size(x%array)
170  y%array(j1, j2) = y%array(j1, j2) + self%A(j, k)*x%array(k)
171  end do
172  end do
173  end do
174 
175  ! Update buffer regions
176  y%array(1 - p1:0, :) = y%array(ny1 - p1 + 1:ny1, :)
177  y%array(ny1 + 1:ny1 + p1, :) = y%array(1:p1, :)
178  y%array(:, 1 - p2:0) = y%array(:, ny2 - p2 + 1:ny2)
179  y%array(:, ny2 + 1:ny2 + p2) = y%array(:, 1:p2)
180 
181  end associate
182 
183  end associate
184 
185  class default
186  err_msg = "y must be of type sll_t_vector_space_real_array_2d"
187  sll_error(this_sub_name, err_msg)
188 
189  end select
190 
191  class default
192  err_msg = "x must be of type sll_t_vector_space_real_array_1d"
193  sll_error(this_sub_name, err_msg)
194 
195  end select
196 
198 
199  ! Convert dense matrix to array (trivial)
201  class(sll_t_linear_operator_matrix_dense_to_stencil), intent(in) :: self
202  real(wp), intent(inout) :: A(:, :)
203 
204  sll_assert(size(a, 1) == self%s1)
205  sll_assert(size(a, 2) == self%s2)
206 
207  a = self%A
208 
210 
211  ! Free objects
213  class(sll_t_linear_operator_matrix_dense_to_stencil), intent(inout) :: self
214 
215  deallocate (self%A)
216 
218 
integer function, dimension(2) f_linear_operator_matrix_dense_to_stencil__get_shape(self)
Abstract type implementing a generic vector space.
Vector space for wrapping 1D Fortran real arrays.
Vector space for wrapping 2D Fortran real arrays.
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
Abstract base class for all vector spaces.
    Report Typos and Errors