Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_vector_space_real_array_1d.F90
Go to the documentation of this file.
1 
6 !
8 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
9 #include "sll_errors.h"
10 
11  use sll_m_working_precision, only: f64
12 
14 
15  implicit none
16 
18 
19  private
20 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 
23  integer, parameter :: wp = f64
24 
30 
31  real(wp), allocatable :: array(:)
32 
33  contains
34 
36  procedure :: copy => copy__real
37  procedure :: incr => incr__real
38  procedure :: scal => scal__real
39 
41  procedure :: add => add__real
42  procedure :: mult => mult__real
43  procedure :: mult_add => mult_add__real
44  procedure :: incr_mult => incr_mult__real
45  procedure :: lcmb => lcmb__real
46  procedure :: incr_lcmb => incr_lcmb__real
47 
49  procedure :: norm => norm__real
50  procedure :: inner => inner__real
52 
54 
55  ! Error messages
56  character(len=*), parameter :: wrong_type_x = "x not of type 'sll_t_vector_space_real_array_1d'"
57  character(len=*), parameter :: wrong_type_y = "y not of type 'sll_t_vector_space_real_array_1d'"
58 
59 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 contains
61 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
62 
63  subroutine copy__real(self, x)
64  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
65  class(sll_c_vector_space), intent(in) :: x
66 
67  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % copy"
68 
69  select type (x)
70 
72 
73  ! Automatic allocation
74  self%array(:) = x%array
75 
76  class default
77 
78  sll_error(this_sub_name, wrong_type_x)
79 
80  end select
81 
82  end subroutine copy__real
83 
84  !-----------------------------------------------------------------------------
85  subroutine incr__real(self, x)
86  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
87  class(sll_c_vector_space), intent(in) :: x
88 
89  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % incr"
90 
91  select type (x)
92 
94 
95  self%array(:) = self%array(:) + x%array(:)
96 
97  class default
98 
99  sll_error(this_sub_name, wrong_type_x)
100 
101  end select
102 
103  end subroutine incr__real
104 
105  !-----------------------------------------------------------------------------
106  subroutine scal__real(self, a)
107  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
108  real(wp), intent(in) :: a
109 
110  self%array(:) = self%array(:)*a
111 
112  end subroutine scal__real
113 
114  !-----------------------------------------------------------------------------
115  subroutine add__real(self, x, y)
116  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
117  class(sll_c_vector_space), intent(in) :: x
118  class(sll_c_vector_space), intent(in) :: y
119 
120  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % add"
121 
122  select type (x)
123 
125 
126  select type (y)
127 
129 
130  self%array(:) = x%array(:) + y%array(:)
131 
132  class default
133 
134  sll_error(this_sub_name, wrong_type_y)
135 
136  end select
137 
138  class default
139 
140  sll_error(this_sub_name, wrong_type_x)
141 
142  end select
143 
144  end subroutine add__real
145 
146  !-----------------------------------------------------------------------------
147  subroutine mult__real(self, a, x)
148  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
149  real(wp), intent(in) :: a
150  class(sll_c_vector_space), intent(in) :: x
151 
152  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % mult"
153 
154  select type (x)
155 
157 
158  self%array(:) = a*x%array(:)
159 
160  class default
161 
162  sll_error(this_sub_name, wrong_type_x)
163 
164  end select
165 
166  end subroutine mult__real
167 
168  !-----------------------------------------------------------------------------
169  subroutine mult_add__real(self, a, x, y)
170  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
171  real(wp), intent(in) :: a
172  class(sll_c_vector_space), intent(in) :: x
173  class(sll_c_vector_space), intent(in) :: y
174 
175  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % mult_add"
176 
177  select type (x)
178 
180 
181  select type (y)
182 
184 
185  self%array(:) = a*x%array(:) + y%array(:)
186 
187  class default
188 
189  sll_error(this_sub_name, wrong_type_y)
190 
191  end select
192 
193  class default
194 
195  sll_error(this_sub_name, wrong_type_x)
196 
197  end select
198 
199  end subroutine mult_add__real
200 
201  !-----------------------------------------------------------------------------
202  subroutine incr_mult__real(self, a, x)
203  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
204  real(wp), intent(in) :: a
205  class(sll_c_vector_space), intent(in) :: x
206 
207  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % incr_mult"
208 
209  select type (x)
210 
212 
213  self%array(:) = self%array(:) + a*x%array(:)
214 
215  class default
216 
217  sll_error(this_sub_name, wrong_type_x)
218 
219  end select
220 
221  end subroutine incr_mult__real
222 
223  !----------------------------------------------------------------------------
224  subroutine lcmb__real(self, a, x)
225  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
226  real(wp), intent(in) :: a(:)
227  class(sll_c_vector_space), intent(in) :: x(:)
228 
229  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % lcmb"
230 
231  integer :: i
232 
233  select type (x)
234 
236 
237  ! Construct linear combination incrementally
238  self%array(:) = a(1)*x(1)%array(:)
239  do i = 2, size(a)
240  self%array(:) = self%array(:) + a(i)*x(i)%array(:)
241  end do
242 
243  class default
244 
245  sll_error(this_sub_name, wrong_type_x)
246 
247  end select
248 
249  end subroutine lcmb__real
250 
251  !-----------------------------------------------------------------------------
252  subroutine incr_lcmb__real(self, a, x)
253  class(sll_t_vector_space_real_array_1d), intent(inout) :: self
254  real(wp), intent(in) :: a(:)
255  class(sll_c_vector_space), intent(in) :: x(:)
256 
257  character(len=*), parameter :: this_sub_name = "sll_t_vector_space_real_array_1d % incr_lcmb"
258 
259  integer :: i
260 
261  select type (x)
262 
264 
265  ! Construct linear combination incrementally
266  do i = 1, size(a)
267  self%array(:) = self%array(:) + a(i)*x(i)%array(:)
268  end do
269 
270  class default
271 
272  sll_error(this_sub_name, wrong_type_x)
273 
274  end select
275 
276  end subroutine incr_lcmb__real
277 
278  !-----------------------------------------------------------------------------
279  function norm__real(self) result(res)
280  class(sll_t_vector_space_real_array_1d), intent(in) :: self
281  real(wp) :: res
282 
283  res = norm2(self%array)
284 
285  end function norm__real
286 
287  !-----------------------------------------------------------------------------
288  function inner__real(self, x) result(res)
289  class(sll_t_vector_space_real_array_1d), intent(in) :: self
290  class(sll_c_vector_space), intent(in) :: x
291  real(wp) :: res
292 
293  character(len=*), parameter :: this_fun_name = "sll_t_vector_space_real_array_1d % inner"
294 
295  select type (x)
296 
298 
299  res = dot_product(self%array, x%array)
300 
301  class default
302 
303  sll_error(this_fun_name, wrong_type_x)
304 
305  end select
306 
307  end function inner__real
308 
Abstract type implementing a generic vector space.
Vector space for wrapping 1D Fortran real arrays.
integer, parameter wp
Working precision.
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
function scal(x, y, n)
Definition: sol.f:173
Abstract base class for all vector spaces.
    Report Typos and Errors