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