Neko  0.8.1
A portable framework for high-order spectral element flow simulations
matrix.f90
Go to the documentation of this file.
1 ! Copyright (c) 2023, 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 matrix
35  use neko_config
36  use num_types
37  use device
38  use device_math
39  use utils
40  use, intrinsic :: iso_c_binding
41  implicit none
42  private
43 
44  type, public :: matrix_t
45  real(kind=rp), allocatable :: x(:,:)
46  type(c_ptr) :: x_d = c_null_ptr
47  integer :: nrows = 0
48  integer :: ncols = 0
49  integer :: n = 0
50  contains
52  procedure, pass(m) :: init => matrix_init
54  procedure, pass(m) :: free => matrix_free
56  procedure, pass(m) :: size => matrix_size
58  procedure, pass(m) :: matrix_assign_matrix
60  procedure, pass(m) :: matrix_assign_scalar
61  generic :: assignment(=) => matrix_assign_matrix, &
63  end type matrix_t
64 
65 contains
66 
70  subroutine matrix_init(m, nrows, ncols)
71  class(matrix_t), intent(inout) :: m
72  integer, intent(in) :: nrows
73  integer, intent(in) :: ncols
74 
75  call m%free()
76 
77  allocate(m%x(nrows, ncols))
78  m%x = 0.0_rp
79  m%nrows = nrows
80  m%ncols = ncols
81  m%n = nrows*ncols
82 
83  if (neko_bcknd_device .eq. 1) then
84  call device_map(m%x, m%x_d, m%n)
85  call device_cfill(m%x_d, 0.0_rp, m%n)
86  end if
87 
88  end subroutine matrix_init
89 
91  subroutine matrix_free(m)
92  class(matrix_t), intent(inout) :: m
93 
94  if (allocated(m%x)) then
95  deallocate(m%x)
96  end if
97 
98  if (c_associated(m%x_d)) then
99  call device_free(m%x_d)
100  end if
101 
102  m%nrows = 0
103  m%ncols = 0
104  m%n = 0
105 
106  end subroutine matrix_free
107 
109  function matrix_size(m) result(s)
110  class(matrix_t), intent(inout) :: m
111  integer :: s
112  s = m%n
113  end function matrix_size
114 
116  subroutine matrix_assign_matrix(m, w)
117  class(matrix_t), intent(inout) :: m
118  type(matrix_t), intent(in) :: w
119 
120  if (allocated(m%x)) then
121  call m%free()
122  end if
123 
124  if (.not. allocated(m%x)) then
125 
126  m%nrows = w%nrows
127  m%ncols = w%ncols
128  m%n = w%n
129  allocate(m%x(m%nrows, m%ncols))
130 
131  if (neko_bcknd_device .eq. 1) then
132  call device_map(m%x, m%x_d, m%n)
133  end if
134 
135  end if
136 
137  if (neko_bcknd_device .eq. 1) then
138  call device_copy(m%x_d, w%x_d, m%n)
139  else
140  m%x = w%x
141  end if
142 
143  end subroutine matrix_assign_matrix
144 
146  subroutine matrix_assign_scalar(m, s)
147  class(matrix_t), intent(inout) :: m
148  real(kind=rp), intent(in) :: s
149 
150  if (.not. allocated(m%x)) then
151  call neko_error('matrix not allocated')
152  end if
153 
154  if (neko_bcknd_device .eq. 1) then
155  call device_cfill(m%x_d, s, m%n)
156  else
157  m%x = s
158  end if
159 
160  end subroutine matrix_assign_scalar
161 
162 
163 end module matrix
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
subroutine, public device_copy(a_d, b_d, n)
subroutine, public device_cfill(a_d, c, n)
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:172
Defines a matrix.
Definition: matrix.f90:34
subroutine matrix_free(m)
Deallocate a matrix.
Definition: matrix.f90:92
subroutine matrix_assign_matrix(m, w)
Assignment .
Definition: matrix.f90:117
integer function matrix_size(m)
Returns the number of entries in the matrix.
Definition: matrix.f90:110
subroutine matrix_assign_scalar(m, s)
Assignment .
Definition: matrix.f90:147
subroutine matrix_init(m, nrows, ncols)
Initialise a matrix of size nrows*ncols.
Definition: matrix.f90:71
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