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_new.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  integer :: s3
33 
34  contains
35 
41 
43 
44 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
45 contains
46 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 
48  ! Initialize linear operator
50  class(sll_t_linear_operator_matrix_stencil_to_dense_new), intent(inout) :: self
51  integer, intent(in) :: s1
52  integer, intent(in) :: s2
53  integer, intent(in) :: s3
54 
55  allocate (self%A(s1, s2, s3))
56 
57  self%s1 = s1
58  self%s2 = s2
59  self%s3 = s3
60 
62 
63  ! Get shape of linear operator
65  class(sll_t_linear_operator_matrix_stencil_to_dense_new), intent(in) :: self
66  integer :: s(2)
67 
68  s(1) = size(self%A, 1)
69  s(2) = size(self%A, 2)*size(self%A, 3)
70 
72 
73  ! Implement y=Ax, with A dense matrix (2D), x stencil vector (2D) and y dense vector (1D)
75  class(sll_t_linear_operator_matrix_stencil_to_dense_new), intent(in) :: self
76  class(sll_c_vector_space), intent(in) :: x ! stencil
77  class(sll_c_vector_space), intent(inout) :: y ! dense
78 
79  integer :: j1, j2, i
80  real(wp) :: temp
81 
82  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_stencil_to_dense_new % dot"
83  character(len=64) :: err_msg
84 
85  select type (x)
86 
88 
89  associate(p1 => -lbound(x%array, 1) + 1, &
90  p2 => -lbound(x%array, 2) + 1)
91 
92  associate(nx1 => ubound(x%array, 1) - p1, &
93  nx2 => ubound(x%array, 2) - p2)
94 
95  select type (y)
96 
98 
99  ! Check dimensions
100  sll_assert(self%s1 == size(y%array))
101 
102  y%array = 0.0_wp
103 
104 !$OMP PARALLEL DO PRIVATE(temp)
105  do i = 1, size(y%array)
106 
107  temp = 0.0_wp
108 
109  do j2 = 1, nx2
110  do j1 = 1, p1
111  temp = temp + self%A(i, j1, j2)*x%array(j1, j2)
112  end do
113  end do
114 
115  y%array(i) = y%array(i) + temp
116 
117  end do
118 !$OMP END PARALLEL DO
119 
120  class default
121  err_msg = "y must be of type sll_t_vector_space_real_array_1d"
122  sll_error(this_sub_name, err_msg)
123 
124  end select
125 
126  end associate
127 
128  end associate
129 
130  class default
131  err_msg = "x must be of type sll_t_vector_space_real_array_2d"
132  sll_error(this_sub_name, err_msg)
133 
134  end select
135 
137 
138  ! Implement y=y+Ax, with A dense matrix (2D), x stencil vector (2D) and y dense vector (1D)
140  class(sll_t_linear_operator_matrix_stencil_to_dense_new), intent(in) :: self
141  class(sll_c_vector_space), intent(in) :: x ! stencil
142  class(sll_c_vector_space), intent(inout) :: y ! dense
143 
144  integer :: j1, j2, i
145  real(wp) :: temp
146 
147  character(len=*), parameter :: this_sub_name = "sll_t_linear_operator_matrix_stencil_to_dense_new % dot"
148  character(len=64) :: err_msg
149 
150  select type (x)
151 
152  type is (sll_t_vector_space_real_array_2d)
153 
154  associate(p1 => -lbound(x%array, 1) + 1, &
155  p2 => -lbound(x%array, 2) + 1)
156 
157  associate(nx1 => ubound(x%array, 1) - p1, &
158  nx2 => ubound(x%array, 2) - p2)
159 
160  select type (y)
161 
162  type is (sll_t_vector_space_real_array_1d)
163 
164  ! Check dimensions
165  sll_assert(self%s1 == size(y%array))
166 
167 !$OMP PARALLEL DO PRIVATE(temp)
168  do i = 1, size(y%array)
169 
170  temp = 0.0_wp
171 
172  do j2 = 1, nx2
173  do j1 = 1, p1
174  temp = temp + self%A(i, j1, j2)*x%array(j1, j2)
175  end do
176  end do
177 
178  y%array(i) = y%array(i) + temp
179 
180  end do
181 !$OMP END PARALLEL DO
182 
183  class default
184  err_msg = "y must be of type sll_t_vector_space_real_array_1d"
185  sll_error(this_sub_name, err_msg)
186 
187  end select
188 
189  end associate
190 
191  end associate
192 
193  class default
194  err_msg = "x must be of type sll_t_vector_space_real_array_2d"
195  sll_error(this_sub_name, err_msg)
196 
197  end select
198 
200 
201  ! Free objects
203  class(sll_t_linear_operator_matrix_stencil_to_dense_new), intent(inout) :: self
204 
205  deallocate (self%A)
206 
208 
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