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_stencil_to_dense.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_stencil_to_dense), 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_stencil_to_dense), 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 stencil vector (2D) and y dense vector (1D)
72  class(sll_t_linear_operator_matrix_stencil_to_dense), intent(in) :: self
73  class(sll_c_vector_space), intent(in) :: x ! stencil
74  class(sll_c_vector_space), intent(inout) :: y ! dense
75 
76  integer :: j1, j2, i, j
77 
78  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_stencil_to_dense % dot"
79  character(len=64) :: err_msg
80 
81  select type (x)
82 
84 
85  associate(p1 => -lbound(x%array, 1) + 1, &
86  p2 => -lbound(x%array, 2) + 1)
87 
88  associate(nx1 => ubound(x%array, 1) - p1, &
89  nx2 => ubound(x%array, 2) - p2)
90 
91  ! Check dimensions
92  sll_assert(self%s2 == nx1*nx2)
93 
94  select type (y)
95 
97 
98  ! Check dimensions
99  sll_assert(self%s1 == size(y%array))
100 
101  y%array = 0.0_wp
102  do i = 1, size(y%array)
103  do j2 = 1, nx2
104  do j1 = 1, nx1
105  j = (j1 - 1)*nx2 + j2
106  y%array(i) = y%array(i) + self%A(i, j)*x%array(j1, j2)
107  end do
108  end do
109  end do
110 
111  class default
112  err_msg = "y must be of type sll_t_vector_space_real_array_1d"
113  sll_error(this_sub_name, err_msg)
114 
115  end select
116 
117  end associate
118 
119  end associate
120 
121  class default
122  err_msg = "x must be of type sll_t_vector_space_real_array_2d"
123  sll_error(this_sub_name, err_msg)
124 
125  end select
126 
128 
129  ! Implement y=y+Ax, with A dense matrix (2D), x stencil vector (2D) and y dense vector (1D)
131  class(sll_t_linear_operator_matrix_stencil_to_dense), intent(in) :: self
132  class(sll_c_vector_space), intent(in) :: x ! stencil
133  class(sll_c_vector_space), intent(inout) :: y ! dense
134 
135  integer :: j1, j2, i, j
136 
137  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_stencil_to_dense % dot"
138  character(len=64) :: err_msg
139 
140  select type (x)
141 
142  type is (sll_t_vector_space_real_array_2d)
143 
144  associate(p1 => -lbound(x%array, 1) + 1, &
145  p2 => -lbound(x%array, 2) + 1)
146 
147  associate(nx1 => ubound(x%array, 1) - p1, &
148  nx2 => ubound(x%array, 2) - p2)
149 
150  ! Check dimensions
151  sll_assert(self%s2 == nx1*nx2)
152 
153  select type (y)
154 
155  type is (sll_t_vector_space_real_array_1d)
156 
157  ! Check dimensions
158  sll_assert(self%s1 == size(y%array))
159 
160  do i = 1, size(y%array)
161  do j2 = 1, nx2
162  do j1 = 1, nx1
163  j = (j1 - 1)*nx2 + j2
164  y%array(i) = y%array(i) + self%A(i, j)*x%array(j1, j2)
165  end do
166  end do
167  end do
168 
169  class default
170  err_msg = "y must be of type sll_t_vector_space_real_array_1d"
171  sll_error(this_sub_name, err_msg)
172 
173  end select
174 
175  end associate
176 
177  end associate
178 
179  class default
180  err_msg = "x must be of type sll_t_vector_space_real_array_2d"
181  sll_error(this_sub_name, err_msg)
182 
183  end select
184 
186 
187  ! Convert dense matrix to array (trivial)
189  class(sll_t_linear_operator_matrix_stencil_to_dense), intent(in) :: self
190  real(wp), intent(inout) :: A(:, :)
191 
192  sll_assert(size(a, 1) == self%s1)
193  sll_assert(size(a, 2) == self%s2)
194 
195  a = self%A
196 
198 
199  ! Free objects
201  class(sll_t_linear_operator_matrix_stencil_to_dense), intent(inout) :: self
202 
203  deallocate (self%A)
204 
206 
integer function, dimension(2) f_linear_operator_matrix_stencil_to_dense__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