Neko 1.99.1
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 math, only: sub3, chsign, add3, cmult2, cadd2
37 use num_types, only: rp
38 use device, only: device_map, device_free
41 use utils, only: neko_error
42 use, intrinsic :: iso_c_binding
43 implicit none
44 private
45
46 type, public :: matrix_t
47 real(kind=rp), allocatable :: x(:,:)
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) :: get_nrows => matrix_nrows
62 procedure, pass(m) :: get_ncols => matrix_ncols
64 procedure, pass(m) :: matrix_assign_matrix
66 procedure, pass(m) :: matrix_assign_scalar
68 procedure, pass(m) :: inverse => matrix_bcknd_inverse
69
70 generic :: assignment(=) => matrix_assign_matrix, &
72
74 procedure, pass(m), private :: alloc => matrix_allocate
75 end type matrix_t
76
77contains
78
82 subroutine matrix_init(m, nrows, ncols)
83 class(matrix_t), intent(inout) :: m
84 integer, intent(in) :: nrows
85 integer, intent(in) :: ncols
86
87 call m%alloc(nrows, ncols)
88 m%x = 0.0_rp
89 if (neko_bcknd_device .eq. 1) then
90 call device_cfill(m%x_d, 0.0_rp, m%n)
91 end if
92
93 end subroutine matrix_init
94
98 subroutine matrix_allocate(m, nrows, ncols)
99 class(matrix_t), intent(inout) :: m
100 integer, intent(in) :: nrows
101 integer, intent(in) :: ncols
102
103 call m%free()
104
105 allocate(m%x(nrows, ncols))
106 m%nrows = nrows
107 m%ncols = ncols
108 m%n = nrows*ncols
109
110 if (neko_bcknd_device .eq. 1) then
111 call device_map(m%x, m%x_d, m%n)
112 end if
113
114 end subroutine matrix_allocate
115
117 subroutine matrix_free(m)
118 class(matrix_t), intent(inout) :: m
119
120 if (allocated(m%x)) then
121 deallocate(m%x)
122 end if
123
124 if (c_associated(m%x_d)) then
125 call device_free(m%x_d)
126 end if
127
128 m%nrows = 0
129 m%ncols = 0
130 m%n = 0
131
132 end subroutine matrix_free
133
135 pure function matrix_size(m) result(s)
136 class(matrix_t), intent(in) :: m
137 integer :: s
138 s = m%n
139 end function matrix_size
140
142 pure function matrix_nrows(m) result(nr)
143 class(matrix_t), intent(in) :: m
144 integer :: nr
145 nr = m%nrows
146 end function matrix_nrows
147
149 pure function matrix_ncols(m) result(nc)
150 class(matrix_t), intent(in) :: m
151 integer :: nc
152 nc = m%ncols
153 end function matrix_ncols
154
156 subroutine matrix_assign_matrix(m, w)
157 class(matrix_t), intent(inout) :: m
158 type(matrix_t), intent(in) :: w
159
160 if (allocated(m%x)) then
161 call m%free()
162 end if
163
164 if (.not. allocated(m%x)) then
165
166 m%nrows = w%nrows
167 m%ncols = w%ncols
168 m%n = w%n
169 allocate(m%x(m%nrows, m%ncols))
170
171 if (neko_bcknd_device .eq. 1) then
172 call device_map(m%x, m%x_d, m%n)
173 end if
174
175 end if
176
177 if (neko_bcknd_device .eq. 1) then
178 call device_copy(m%x_d, w%x_d, m%n)
179 else
180 m%x = w%x
181 end if
182
183 end subroutine matrix_assign_matrix
184
186 subroutine matrix_assign_scalar(m, s)
187 class(matrix_t), intent(inout) :: m
188 real(kind=rp), intent(in) :: s
189
190 if (.not. allocated(m%x)) then
191 call neko_error('matrix not allocated')
192 end if
193
194 if (neko_bcknd_device .eq. 1) then
195 call device_cfill(m%x_d, s, m%n)
196 else
197 m%x = s
198 end if
199
200 end subroutine matrix_assign_scalar
201
202 subroutine matrix_bcknd_inverse(m, bcknd)
203 class(matrix_t), intent(inout) :: m
204 integer, optional :: bcknd
205
206 if (neko_bcknd_device .eq. 1 .and. &
207 bcknd .eq. neko_bcknd_device) then
208 call neko_error("matrix_bcknd_inverse not &
209 &implemented on accelarators.")
210 else
211 call cpu_matrix_inverse(m)
212 end if
213
214 end subroutine matrix_bcknd_inverse
215
216 subroutine cpu_matrix_inverse(m)
217 ! Gauss-Jordan matrix inversion with full pivoting
218 ! Num. Rec. p. 30, 2nd Ed., Fortran
219 ! m%x is an sqaure matrix
220 ! rmult is m work array of length nrows = ncols
221 class(matrix_t), intent(inout) :: m
222 integer :: indr(m%nrows), indc(m%ncols), ipiv(m%ncols)
223 real(kind=rp) :: rmult(m%nrows), amx, tmp, piv, eps
224 integer :: i, j, k, ir, jc
225
226 if (.not. (m%ncols .eq. m%nrows)) then
227 call neko_error("Fatal error: trying to invert m matrix that is not &
228 &square")
229 end if
230
231 eps = 1e-9_rp
232 ipiv = 0
233
234 do k = 1, m%nrows
235 amx = 0.0_rp
236 do i = 1, m%nrows ! Pivot search
237 if (ipiv(i) .ne. 1) then
238 do j = 1, m%nrows
239 if (ipiv(j) .eq. 0) then
240 if (abs(m%x(i, j)) .ge. amx) then
241 amx = abs(m%x(i, j))
242 ir = i
243 jc = j
244 end if
245 else if (ipiv(j) .gt. 1) then
246 return
247 end if
248 end do
249 end if
250 end do
251 ipiv(jc) = ipiv(jc) + 1
252
253 ! Swap rows
254 if (ir .ne. jc) then
255 do j = 1, m%ncols
256 tmp = m%x(ir, j)
257 m%x(ir, j) = m%x(jc, j)
258 m%x(jc, j) = tmp
259 end do
260 end if
261 indr(k) = ir
262 indc(k) = jc
263
264 if (abs(m%x(jc, jc)) .lt. eps) then
265 call neko_error("matrix_inverse error: small Gauss Jordan Piv")
266 end if
267 piv = 1.0_rp/m%x(jc, jc)
268 m%x(jc, jc) = 1.0_rp
269 do j = 1, m%ncols
270 m%x(jc, j) = m%x(jc, j)*piv
271 end do
272
273 do j = 1, m%ncols
274 tmp = m%x(jc, j)
275 m%x(jc, j) = m%x(1 , j)
276 m%x(1 , j) = tmp
277 end do
278 do i = 2, m%nrows
279 rmult(i) = m%x(i, jc)
280 m%x(i, jc) = 0.0_rp
281 end do
282
283 do j = 1, m%ncols
284 do i = 2, m%nrows
285 m%x(i, j) = m%x(i, j) - rmult(i)*m%x(1, j)
286 end do
287 end do
288
289 do j = 1, m%ncols
290 tmp = m%x(jc, j)
291 m%x(jc, j) = m%x(1 , j)
292 m%x(1 , j) = tmp
293 end do
294 end do
295
296 ! Unscramble matrix
297 do j= m%nrows, 1, -1
298 if (indr(j) .ne. indc(j)) then
299 do i = 1, m%nrows
300 tmp = m%x(i, indr(j))
301 m%x(i, indr(j)) = m%x(i, indc(j))
302 m%x(i, indc(j)) = tmp
303 end do
304 end if
305 end do
306
307 return
308 end subroutine cpu_matrix_inverse
309
310end module matrix
Map a Fortran array to a device (allocate and associate)
Definition device.F90:72
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_cadd2(a_d, b_d, c, n, strm)
Add a scalar to vector .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_cmult2(a_d, b_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_add3(a_d, b_d, c_d, n, strm)
Vector addition .
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:214
Definition math.f90:60
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:429
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition math.f90:480
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:745
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:787
subroutine, public chsign(a, n)
Change sign of vector .
Definition math.f90:585
Defines a matrix.
Definition matrix.f90:34
subroutine matrix_free(m)
Deallocate a matrix.
Definition matrix.f90:118
subroutine matrix_bcknd_inverse(m, bcknd)
Definition matrix.f90:203
pure integer function matrix_ncols(m)
Returns the number of columns in the matrix.
Definition matrix.f90:150
subroutine matrix_assign_matrix(m, w)
Assignment .
Definition matrix.f90:157
subroutine cpu_matrix_inverse(m)
Definition matrix.f90:217
pure integer function matrix_nrows(m)
Returns the number of rows in the matrix.
Definition matrix.f90:143
pure integer function matrix_size(m)
Returns the number of entries in the matrix.
Definition matrix.f90:136
subroutine matrix_assign_scalar(m, s)
Assignment .
Definition matrix.f90:187
subroutine matrix_init(m, nrows, ncols)
Initialise a matrix of size nrows*ncols.
Definition matrix.f90:83
subroutine matrix_allocate(m, nrows, ncols)
Allocate a matrix of size nrows*ncols.
Definition matrix.f90:99
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35