3 #include "sll_assert.h"
4 #include "sll_errors.h"
28 real(
wp),
allocatable :: a(:, :)
51 integer,
intent(in) :: s1
52 integer,
intent(in) :: s2
54 allocate (self%A(s1, s2))
76 integer :: j1, j2, i, j
78 character(len=*),
parameter :: this_sub_name =
"sll_t_linear_operator_matrix_stencil_to_dense % dot"
79 character(len=64) :: err_msg
85 associate(p1 => -lbound(x%array, 1) + 1, &
86 p2 => -lbound(x%array, 2) + 1)
88 associate(nx1 => ubound(x%array, 1) - p1, &
89 nx2 => ubound(x%array, 2) - p2)
92 sll_assert(self%s2 == nx1*nx2)
99 sll_assert(self%s1 ==
size(y%array))
102 do i = 1,
size(y%array)
105 j = (j1 - 1)*nx2 + j2
106 y%array(i) = y%array(i) + self%A(i, j)*x%array(j1, j2)
112 err_msg =
"y must be of type sll_t_vector_space_real_array_1d"
113 sll_error(this_sub_name, err_msg)
122 err_msg =
"x must be of type sll_t_vector_space_real_array_2d"
123 sll_error(this_sub_name, err_msg)
132 class(sll_c_vector_space),
intent(in) :: x
133 class(sll_c_vector_space),
intent(inout) :: y
135 integer :: j1, j2, i, j
137 character(len=*),
parameter :: this_sub_name =
"sll_t_linear_operator_matrix_stencil_to_dense % dot"
138 character(len=64) :: err_msg
142 type is (sll_t_vector_space_real_array_2d)
144 associate(p1 => -lbound(x%array, 1) + 1, &
145 p2 => -lbound(x%array, 2) + 1)
147 associate(nx1 => ubound(x%array, 1) - p1, &
148 nx2 => ubound(x%array, 2) - p2)
151 sll_assert(self%s2 == nx1*nx2)
155 type is (sll_t_vector_space_real_array_1d)
158 sll_assert(self%s1 ==
size(y%array))
160 do i = 1,
size(y%array)
163 j = (j1 - 1)*nx2 + j2
164 y%array(i) = y%array(i) + self%A(i, j)*x%array(j1, j2)
170 err_msg =
"y must be of type sll_t_vector_space_real_array_1d"
171 sll_error(this_sub_name, err_msg)
180 err_msg =
"x must be of type sll_t_vector_space_real_array_2d"
181 sll_error(this_sub_name, err_msg)
190 real(wp),
intent(inout) :: A(:, :)
192 sll_assert(
size(a, 1) == self%s1)
193 sll_assert(
size(a, 2) == self%s2)
subroutine s_linear_operator_matrix_stencil_to_dense__dot_incr(self, x, y)
integer function, dimension(2) f_linear_operator_matrix_stencil_to_dense__get_shape(self)
subroutine s_linear_operator_matrix_stencil_to_dense__free(self)
subroutine s_linear_operator_matrix_stencil_to_dense__to_array(self, A)
subroutine s_linear_operator_matrix_stencil_to_dense__init(self, s1, s2)
subroutine s_linear_operator_matrix_stencil_to_dense__dot(self, x, y)
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.
Vector space for wrapping 1D Fortran real arrays.
Vector space for wrapping 2D Fortran real arrays.