Neko  0.8.99
A portable framework for high-order spectral element flow simulations
vector.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module vector
36  use math, only: sub3, chsign, add3, cmult2, cadd2, cfill, copy
37  use num_types, only: rp
38  use device, only: device_map, device_free, c_ptr, c_null_ptr
41  use utils, only: neko_error
42  use, intrinsic :: iso_c_binding
43  implicit none
44  private
45 
46  type, public :: vector_t
48  real(kind=rp), allocatable :: x(:)
50  type(c_ptr) :: x_d = c_null_ptr
52  integer :: n = 0
53  contains
55  procedure, pass(v) :: init => vector_init
57  procedure, pass(v) :: free => vector_free
59  procedure, pass(v) :: size => vector_size
61  procedure, pass(v) :: vector_assign_vector
63  procedure, pass(v) :: vector_assign_scalar
65  procedure, pass(a) :: vector_add_vector
67  procedure, pass(a) :: vector_add_scalar_left
69  procedure, pass(a) :: vector_add_scalar_right
71  procedure, pass(a) :: vector_sub_vector
73  procedure, pass(a) :: vector_sub_scalar_left
75  procedure, pass(a) :: vector_sub_scalar_right
77  procedure, pass(a) :: vector_cmult_left
79  procedure, pass(a) :: vector_cmult_right
80 
81  generic :: assignment(=) => vector_assign_vector, &
83  generic :: operator(+) => vector_add_vector, &
85  generic :: operator(-) => vector_sub_vector, &
87  generic :: operator(*) => vector_cmult_left, vector_cmult_right
88  end type vector_t
89 
90  type, public :: vector_ptr_t
91  type(vector_t), pointer :: ptr
92  end type vector_ptr_t
93 
94 contains
95 
97  subroutine vector_init(v, n)
98  class(vector_t), intent(inout) :: v
99  integer, intent(in) :: n
100 
101  call v%free()
102 
103  allocate(v%x(n))
104  v%x = 0.0_rp
105 
106  if (neko_bcknd_device .eq. 1) then
107  call device_map(v%x, v%x_d, n)
108  call device_cfill(v%x_d, 0.0_rp, n)
109  end if
110 
111  v%n = n
112 
113  end subroutine vector_init
114 
116  subroutine vector_free(v)
117  class(vector_t), intent(inout) :: v
118 
119  if (allocated(v%x)) then
120  deallocate(v%x)
121  end if
122 
123  if (c_associated(v%x_d)) then
124  call device_free(v%x_d)
125  end if
126 
127  v%n = 0
128 
129  end subroutine vector_free
130 
132  function vector_size(v) result(s)
133  class(vector_t), intent(inout) :: v
134  integer :: s
135  s = v%n
136  end function vector_size
137 
139  subroutine vector_assign_vector(v, w)
140  class(vector_t), intent(inout) :: v
141  type(vector_t), intent(in) :: w
142 
143  if (allocated(v%x)) then
144  call v%free()
145  end if
146 
147  if (.not. allocated(v%x)) then
148 
149  v%n = w%n
150  allocate(v%x(v%n))
151 
152  if (neko_bcknd_device .eq. 1) then
153  call device_map(v%x, v%x_d, v%n)
154  end if
155 
156  end if
157 
158  if (neko_bcknd_device .eq. 1) then
159  call device_copy(v%x_d, w%x_d, v%n)
160  else
161  v%x = w%x
162  end if
163 
164  end subroutine vector_assign_vector
165 
167  subroutine vector_assign_scalar(v, s)
168  class(vector_t), intent(inout) :: v
169  real(kind=rp), intent(in) :: s
170 
171  if (.not. allocated(v%x)) then
172  call neko_error('Vector not allocated')
173  end if
174 
175  if (neko_bcknd_device .eq. 1) then
176  call device_cfill(v%x_d, s, v%n)
177  else
178  call cfill(v%x, s, v%n)
179  end if
180 
181  end subroutine vector_assign_scalar
182 
184  function vector_add_vector(a, b) result(v)
185  class(vector_t), intent(in) :: a, b
186  type(vector_t) :: v
187 
188  if (a%n .ne. b%n) call neko_error("Vectors must be the same length!")
189 
190  v%n = a%n
191  allocate(v%x(v%n))
192 
193  if (neko_bcknd_device .eq. 1) then
194  call device_map(v%x, v%x_d, v%n)
195  end if
196 
197  if (neko_bcknd_device .eq. 1) then
198  call device_add3(v%x_d, a%x_d, b%x_d, v%n)
199  else
200  call add3(v%x, a%x, b%x, v%n)
201  end if
202 
203  end function vector_add_vector
204 
206  function vector_add_scalar_left(a, c) result(v)
207  class(vector_t), intent(in) :: a
208  real(kind=rp), intent(in) :: c
209  type(vector_t) :: v
210 
211  v%n = a%n
212  allocate(v%x(v%n))
213 
214  if (neko_bcknd_device .eq. 1) then
215  call device_map(v%x, v%x_d, v%n)
216  end if
217 
218  if (neko_bcknd_device .eq. 1) then
219  call device_cadd2(v%x_d, a%x_d, c, v%n)
220  else
221  call cadd2(v%x, a%x, c, v%n)
222  end if
223 
224  end function vector_add_scalar_left
225 
227  function vector_add_scalar_right(c, a) result(v)
228  real(kind=rp), intent(in) :: c
229  class(vector_t), intent(in) :: a
230  type(vector_t) :: v
231 
232  v = vector_add_scalar_left(a, c)
233 
234  end function vector_add_scalar_right
235 
237  function vector_sub_vector(a, b) result(v)
238  class(vector_t), intent(in) :: a, b
239  type(vector_t) :: v
240 
241  if (a%n .ne. b%n) call neko_error("Vectors must be the same length!")
242 
243  v%n = a%n
244  allocate(v%x(v%n))
245 
246  if (neko_bcknd_device .eq. 1) then
247  call device_map(v%x, v%x_d, v%n)
248  end if
249 
250  if (neko_bcknd_device .eq. 1) then
251  call device_sub3(v%x_d, a%x_d, b%x_d, v%n)
252  else
253  call sub3(v%x, a%x, b%x, v%n)
254  end if
255 
256  end function vector_sub_vector
257 
259  function vector_sub_scalar_left(a, c) result(v)
260  class(vector_t), intent(in) :: a
261  real(kind=rp), intent(in) :: c
262  type(vector_t) :: v
263 
264  v%n = a%n
265  allocate(v%x(v%n))
266 
267  if (neko_bcknd_device .eq. 1) then
268  call device_map(v%x, v%x_d, v%n)
269  end if
270 
271  if (neko_bcknd_device .eq. 1) then
272  call device_cadd2(v%x_d, a%x_d, -1.0_rp*c, v%n)
273  else
274  call cadd2(v%x, a%x, -1.0_rp*c, a%n)
275  end if
276 
277  end function vector_sub_scalar_left
278 
280  function vector_sub_scalar_right(c, a) result(v)
281  real(kind=rp), intent(in) :: c
282  class(vector_t), intent(in) :: a
283  type(vector_t) :: v
284 
285  v = vector_sub_scalar_left(a, c)
286 
287  if (neko_bcknd_device .eq. 1) then
288  call device_cmult(v%x_d, -1.0_rp, v%n)
289  else
290  v%x = -v%x
291  !call chsign(v%x, v%n)
292  end if
293 
294  end function vector_sub_scalar_right
295 
297  function vector_cmult_left(a, c) result(v)
298  class(vector_t), intent(in) :: a
299  real(kind=rp), intent(in) :: c
300  type(vector_t) :: v
301 
302  v%n = a%n
303  allocate(v%x(v%n))
304 
305  if (neko_bcknd_device .eq. 1) then
306  call device_map(v%x, v%x_d, v%n)
307  end if
308 
309  if (neko_bcknd_device .eq. 1) then
310  call device_cmult2(v%x_d, a%x_d, c, v%n)
311  else
312  call cmult2(v%x, a%x, c, v%n)
313  end if
314 
315  end function vector_cmult_left
316 
318  function vector_cmult_right(c, a) result(v)
319  real(kind=rp), intent(in) :: c
320  class(vector_t), intent(in) :: a
321  type(vector_t) :: v
322 
323  v = vector_cmult_left(a, c)
324 
325  end function vector_cmult_right
326 
327 end module vector
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
subroutine, public device_cmult2(a_d, b_d, c, n)
Multiplication by constant c .
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
subroutine, public device_sub3(a_d, b_d, c_d, n)
Vector subtraction .
subroutine, public device_add3(a_d, b_d, c_d, n)
Vector addition .
subroutine, public device_cadd2(a_d, b_d, c, n)
Add a scalar to vector .
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
subroutine, public device_cfill(a_d, c, n)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition: device.F90:185
Definition: math.f90:60
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition: math.f90:661
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition: math.f90:301
subroutine, public add3(a, b, c, n)
Vector addition .
Definition: math.f90:560
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition: math.f90:602
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition: math.f90:314
subroutine, public chsign(a, n)
Change sign of vector .
Definition: math.f90:400
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:228
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35
Defines a vector.
Definition: vector.f90:34
subroutine vector_init(v, n)
Initialise a vector of size n.
Definition: vector.f90:98
subroutine vector_assign_scalar(v, s)
Assignment .
Definition: vector.f90:168
type(vector_t) function vector_add_scalar_right(c, a)
Scalar-vector addition .
Definition: vector.f90:228
integer function vector_size(v)
Return the number of entries in the vector.
Definition: vector.f90:133
type(vector_t) function vector_cmult_left(a, c)
Vector-scalar multiplication .
Definition: vector.f90:298
subroutine vector_free(v)
Deallocate a vector.
Definition: vector.f90:117
type(vector_t) function vector_add_scalar_left(a, c)
Vector-scalar addition .
Definition: vector.f90:207
type(vector_t) function vector_sub_scalar_right(c, a)
Scalar-vector subtraction .
Definition: vector.f90:281
type(vector_t) function vector_add_vector(a, b)
Vector-vector addition .
Definition: vector.f90:185
type(vector_t) function vector_sub_scalar_left(a, c)
Vector-scalar subtraction .
Definition: vector.f90:260
type(vector_t) function vector_cmult_right(c, a)
Scalar-vector multiplication .
Definition: vector.f90:319
type(vector_t) function vector_sub_vector(a, b)
Vector-vector subtraction .
Definition: vector.f90:238
subroutine vector_assign_vector(v, w)
Assignment .
Definition: vector.f90:140