Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
matrix.f90
Go to the documentation of this file.
1! Copyright (c) 2024, 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!
34module matrix
36 use num_types, only : rp, xp
40 use utils, only : neko_error
41 use, intrinsic :: iso_c_binding
42 implicit none
43 private
44
45 type, public :: matrix_t
46 real(kind=rp), allocatable :: x(:,:)
47 type(c_ptr) :: x_d = c_null_ptr
48 integer, private :: nrows = 0
49 integer, private :: ncols = 0
50 integer, private :: n = 0
51 contains
53 procedure, pass(m) :: init => matrix_init
55 procedure, pass(m) :: free => matrix_free
57 procedure, pass(m) :: size => matrix_size
59 procedure, pass(m) :: copy_from => matrix_copy_from
61 procedure, pass(m) :: get_nrows => matrix_nrows
63 procedure, pass(m) :: get_ncols => matrix_ncols
65 procedure, pass(m) :: matrix_assign_matrix
67 procedure, pass(m) :: matrix_assign_scalar
69 procedure, pass(m) :: inverse => matrix_bcknd_inverse
70 procedure, pass(m) :: inverse_on_host => cpu_matrix_inverse
71
72 generic :: assignment(=) => matrix_assign_matrix, &
74
76 procedure, pass(m), private :: alloc => matrix_allocate
77 end type matrix_t
78
79 type, public :: matrix_ptr_t
80 type(matrix_t), pointer :: ptr
81 end type matrix_ptr_t
82
83contains
84
88 subroutine matrix_init(m, nrows, ncols)
89 class(matrix_t), intent(inout) :: m
90 integer, intent(in) :: nrows
91 integer, intent(in) :: ncols
92
93 call m%alloc(nrows, ncols)
94 m%x = 0.0_rp
95 if (neko_bcknd_device .eq. 1) then
96 call device_cfill(m%x_d, 0.0_rp, m%n)
97 end if
98
99 end subroutine matrix_init
100
104 subroutine matrix_allocate(m, nrows, ncols)
105 class(matrix_t), intent(inout) :: m
106 integer, intent(in) :: nrows
107 integer, intent(in) :: ncols
108
109 call m%free()
110
111 allocate(m%x(nrows, ncols))
112 m%nrows = nrows
113 m%ncols = ncols
114 m%n = nrows*ncols
115
116 if (neko_bcknd_device .eq. 1) then
117 call device_map(m%x, m%x_d, m%n)
118 call device_cfill(m%x_d, 0.0_rp, m%n)
119 call device_sync()
120 end if
121
122 end subroutine matrix_allocate
123
125 subroutine matrix_free(m)
126 class(matrix_t), intent(inout) :: m
127
128 if (allocated(m%x)) then
129 deallocate(m%x)
130 end if
131
132 if (c_associated(m%x_d)) then
133 call device_deassociate(m%x)
134 call device_free(m%x_d)
135 end if
136
137 m%nrows = 0
138 m%ncols = 0
139 m%n = 0
140
141 end subroutine matrix_free
142
144 pure function matrix_size(m) result(s)
145 class(matrix_t), intent(in) :: m
146 integer :: s
147 s = m%n
148 end function matrix_size
149
150
155 subroutine matrix_copy_from(m, memdir, sync)
156 class(matrix_t), intent(inout) :: m
157 integer, intent(in) :: memdir
158 logical, intent(in) :: sync
159
160 if (neko_bcknd_device .eq. 1) then
161 call device_memcpy(m%x, m%x_d, m%n, memdir, sync)
162 end if
163
164 end subroutine matrix_copy_from
165
167 pure function matrix_nrows(m) result(nr)
168 class(matrix_t), intent(in) :: m
169 integer :: nr
170 nr = m%nrows
171 end function matrix_nrows
172
174 pure function matrix_ncols(m) result(nc)
175 class(matrix_t), intent(in) :: m
176 integer :: nc
177 nc = m%ncols
178 end function matrix_ncols
179
181 subroutine matrix_assign_matrix(m, w)
182 class(matrix_t), intent(inout) :: m
183 type(matrix_t), intent(in) :: w
184
185 if (allocated(m%x)) then
186 call m%free()
187 end if
188
189 if (.not. allocated(m%x)) then
190
191 m%nrows = w%nrows
192 m%ncols = w%ncols
193 m%n = w%n
194 allocate(m%x(m%nrows, m%ncols))
195
196 if (neko_bcknd_device .eq. 1) then
197 call device_map(m%x, m%x_d, m%n)
198 end if
199
200 end if
201
202 if (neko_bcknd_device .eq. 1) then
203 call device_copy(m%x_d, w%x_d, m%n)
204 else
205 m%x = w%x
206 end if
207
208 end subroutine matrix_assign_matrix
209
211 subroutine matrix_assign_scalar(m, s)
212 class(matrix_t), intent(inout) :: m
213 real(kind=rp), intent(in) :: s
214
215 if (.not. allocated(m%x)) then
216 call neko_error('matrix not allocated')
217 end if
218
219 if (neko_bcknd_device .eq. 1) then
220 call device_cfill(m%x_d, s, m%n)
221 else
222 m%x = s
223 end if
224
225 end subroutine matrix_assign_scalar
226
227 subroutine matrix_bcknd_inverse(m, bcknd)
228 class(matrix_t), intent(inout) :: m
229 integer, optional :: bcknd
230
231 if (neko_bcknd_device .eq. 1 .and. &
232 bcknd .eq. neko_bcknd_device) then
233 call neko_error("matrix_bcknd_inverse not &
234 &implemented on accelarators.")
235 else
236 call cpu_matrix_inverse(m)
237 end if
238
239 end subroutine matrix_bcknd_inverse
240
241 subroutine cpu_matrix_inverse(m)
242 ! Gauss-Jordan matrix inversion with full pivoting
243 ! Num. Rec. p. 30, 2nd Ed., Fortran
244 ! m%x is an sqaure matrix
245 ! rmult is m work array of length nrows = ncols
246 class(matrix_t), intent(inout) :: m
247 integer :: indr(m%nrows), indc(m%ncols), ipiv(m%ncols)
248 real(kind=xp) :: rmult(m%nrows), amx, tmp, piv, eps
249 integer :: i, j, k, ir, jc
250
251 if (.not. (m%ncols .eq. m%nrows)) then
252 call neko_error("Fatal error: trying to invert m matrix that is not &
253 &square")
254 end if
255
256 eps = 1e-9_rp
257 ipiv = 0
258
259 do k = 1, m%nrows
260 amx = 0.0_rp
261 do i = 1, m%nrows ! Pivot search
262 if (ipiv(i) .ne. 1) then
263 do j = 1, m%nrows
264 if (ipiv(j) .eq. 0) then
265 if (abs(m%x(i, j)) .ge. amx) then
266 amx = abs(m%x(i, j))
267 ir = i
268 jc = j
269 end if
270 else if (ipiv(j) .gt. 1) then
271 return
272 end if
273 end do
274 end if
275 end do
276 ipiv(jc) = ipiv(jc) + 1
277
278 ! Swap rows
279 if (ir .ne. jc) then
280 do j = 1, m%ncols
281 tmp = m%x(ir, j)
282 m%x(ir, j) = m%x(jc, j)
283 m%x(jc, j) = tmp
284 end do
285 end if
286 indr(k) = ir
287 indc(k) = jc
288
289 if (abs(m%x(jc, jc)) .lt. eps) then
290 call neko_error("matrix_inverse error: small Gauss Jordan Piv")
291 end if
292 piv = 1.0_xp/m%x(jc, jc)
293 m%x(jc, jc) = 1.0_xp
294 do j = 1, m%ncols
295 m%x(jc, j) = m%x(jc, j)*piv
296 end do
297
298 do j = 1, m%ncols
299 tmp = m%x(jc, j)
300 m%x(jc, j) = m%x(1 , j)
301 m%x(1 , j) = tmp
302 end do
303 do i = 2, m%nrows
304 rmult(i) = m%x(i, jc)
305 m%x(i, jc) = 0.0_rp
306 end do
307
308 do j = 1, m%ncols
309 do i = 2, m%nrows
310 m%x(i, j) = m%x(i, j) - rmult(i)*m%x(1, j)
311 end do
312 end do
313
314 do j = 1, m%ncols
315 tmp = m%x(jc, j)
316 m%x(jc, j) = m%x(1 , j)
317 m%x(1 , j) = tmp
318 end do
319 end do
320
321 ! Unscramble matrix
322 do j = m%nrows, 1, -1
323 if (indr(j) .ne. indc(j)) then
324 do i = 1, m%nrows
325 tmp = m%x(i, indr(j))
326 m%x(i, indr(j)) = m%x(i, indc(j))
327 m%x(i, indc(j)) = tmp
328 end do
329 end if
330 end do
331
332 return
333 end subroutine cpu_matrix_inverse
334
335end module matrix
Deassociate a Fortran array from a device pointer.
Definition device.F90:95
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
Synchronize a device or stream.
Definition device.F90:107
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_cfill(a_d, c, n, strm)
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:219
Defines a matrix.
Definition matrix.f90:34
subroutine matrix_free(m)
Deallocate a matrix.
Definition matrix.f90:126
subroutine matrix_bcknd_inverse(m, bcknd)
Definition matrix.f90:228
pure integer function matrix_ncols(m)
Returns the number of columns in the matrix.
Definition matrix.f90:175
subroutine matrix_assign_matrix(m, w)
Assignment .
Definition matrix.f90:182
subroutine cpu_matrix_inverse(m)
Definition matrix.f90:242
pure integer function matrix_nrows(m)
Returns the number of rows in the matrix.
Definition matrix.f90:168
pure integer function matrix_size(m)
Returns the number of entries in the matrix.
Definition matrix.f90:145
subroutine matrix_copy_from(m, memdir, sync)
Easy way to copy between host and device.
Definition matrix.f90:156
subroutine matrix_assign_scalar(m, s)
Assignment .
Definition matrix.f90:212
subroutine matrix_init(m, nrows, ncols)
Initialise a matrix of size nrows*ncols.
Definition matrix.f90:89
subroutine matrix_allocate(m, nrows, ncols)
Allocate a matrix of size nrows*ncols.
Definition matrix.f90:105
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public xp
Definition num_types.f90:14
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35