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_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 
14  implicit none
15 
17 
18  private
19 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20 
21  ! Working precision
22  integer, parameter :: wp = f64
23 
25 
26  real(wp), allocatable :: a(:, :)
27  real(wp), allocatable :: at(:, :)
28 
29  ! Logical flag telling whether A was given transposed or not
30  logical :: transposed = .false.
31 
32  integer :: s1
33  integer :: s2
34 
35  contains
36 
42 
44 
45 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
46 contains
47 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
48 
49  ! Initialize linear operator
50  subroutine s_linear_operator_matrix_dense_to_dense__init(self, s1, s2, transposed)
51  class(sll_t_linear_operator_matrix_dense_to_dense), intent(inout) :: self
52  integer, intent(in) :: s1
53  integer, intent(in) :: s2
54  logical, optional, intent(in) :: transposed
55 
56  if (present(transposed)) self%transposed = transposed
57 
58  if (self%transposed) then
59  allocate (self%At(s1, s2))
60  else
61  allocate (self%A(s1, s2))
62  end if
63 
64  self%s1 = s1
65  self%s2 = s2
66 
68 
69  ! Get shape of linear operator
71  class(sll_t_linear_operator_matrix_dense_to_dense), intent(in) :: self
72  integer :: s(2)
73 
74  if (self%transposed) then
75  s = shape(self%At)
76  else
77  s = shape(self%A)
78  end if
79 
81 
82  ! Implement y=Ax, with A dense matrix (2D), x and y dense vectors (1D)
84  class(sll_t_linear_operator_matrix_dense_to_dense), intent(in) :: self
85  class(sll_c_vector_space), intent(in) :: x ! dense
86  class(sll_c_vector_space), intent(inout) :: y ! dense
87 
88  integer :: nx(1), ny(1), n(2)
89 
90  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_dense_to_dense % dot"
91  character(len=64) :: err_msg
92 
93  n = self%get_shape()
94 
95  select type (x)
96 
98 
99  ! Check dimensions
100  nx = shape(x%array)
101  sll_assert(n(2) == nx(1))
102 
103  select type (y)
104 
106 
107  ! Check dimensions
108  ny = shape(y%array)
109  sll_assert(n(1) == ny(1))
110 
111  if (self%transposed) then
112  y%array = matmul(x%array, self%At)
113  else
114  y%array = matmul(self%A, x%array)
115  end if
116 
117  class default
118  err_msg = "y must be of type sll_t_vector_space_real_array_1d"
119  sll_error(this_sub_name, err_msg)
120 
121  end select
122 
123  class default
124  err_msg = "x must be of type sll_t_vector_space_real_array_1d"
125  sll_error(this_sub_name, err_msg)
126 
127  end select
128 
130 
131  ! Implement y=y+Ax, with A dense matrix (2D), x and y dense vectors (1D)
133  class(sll_t_linear_operator_matrix_dense_to_dense), intent(in) :: self
134  class(sll_c_vector_space), intent(in) :: x ! dense
135  class(sll_c_vector_space), intent(inout) :: y ! dense
136 
137  integer :: nx(1), ny(1), n(2)
138 
139  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_dense_to_dense % dot"
140  character(len=64) :: err_msg
141 
142  n = self%get_shape()
143 
144  select type (x)
145 
147 
148  ! Check dimensions
149  nx = shape(x%array)
150  sll_assert(n(2) == nx(1))
151 
152  select type (y)
153 
155 
156  ! Check dimensions
157  ny = shape(y%array)
158  sll_assert(n(1) == ny(1))
159 
160  if (self%transposed) then
161  y%array = y%array + matmul(x%array, self%At)
162  else
163  y%array = y%array + matmul(self%A, x%array)
164  end if
165 
166  class default
167  err_msg = "y must be of type sll_t_vector_space_real_array_1d"
168  sll_error(this_sub_name, err_msg)
169 
170  end select
171 
172  class default
173  err_msg = "x must be of type sll_t_vector_space_real_array_1d"
174  sll_error(this_sub_name, err_msg)
175 
176  end select
177 
179 
180  ! Free objects
182  class(sll_t_linear_operator_matrix_dense_to_dense), intent(inout) :: self
183 
184  if (self%transposed) then
185  deallocate (self%At)
186  else
187  deallocate (self%A)
188  end if
189 
191 
integer function, dimension(2) f_linear_operator_matrix_dense_to_dense__get_shape(self)
subroutine s_linear_operator_matrix_dense_to_dense__init(self, s1, s2, transposed)
Abstract type implementing a generic vector space.
Vector space for wrapping 1D 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