Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2011 Bill Paxton & The MESA Team
4 : !
5 : ! This program is free software: you can redistribute it and/or modify
6 : ! it under the terms of the GNU Lesser General Public License
7 : ! as published by the Free Software Foundation,
8 : ! either version 3 of the License, or (at your option) any later version.
9 : !
10 : ! This program is distributed in the hope that it will be useful,
11 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
12 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
13 : ! See the GNU Lesser General Public License for more details.
14 : !
15 : ! You should have received a copy of the GNU Lesser General Public License
16 : ! along with this program. If not, see <https://www.gnu.org/licenses/>.
17 : !
18 : ! ***********************************************************************
19 :
20 : module test_mtx_support
21 :
22 : use mtx_lib
23 : use mtx_def
24 : use utils_lib, only: mesa_error
25 :
26 : implicit none
27 :
28 : contains
29 :
30 1 : subroutine test_format_conversion
31 : use mtx_def
32 : integer, parameter :: n = 6
33 : integer, parameter :: nzmax = n*n, ndim = n
34 :
35 122 : real(dp) :: a(ndim, n), a2(ndim, n), values(nzmax)
36 : integer, parameter :: ml = 1, mu = 2, ldbb = 2*ml + mu + 1
37 115 : real(dp) :: b(ndim, n), bb(ldbb, n), bb2(ldbb, n)
38 : integer :: ierr, nz, iptr(n + 1), jind(nzmax), i, j, k, kk, hint
39 :
40 1 : write (*, *) 'test_format_conversion'
41 :
42 7 : a(1, 1:n) = [10d0, 0d0, 0d0, 0d0, 0d0, 0d0]
43 7 : a(2, 1:n) = [0d0, 12d0, -3d0, -1d0, 0d0, 0d0]
44 7 : a(3, 1:n) = [0d0, 0d0, 15d0, 0d0, 0d0, 0d0]
45 7 : a(4, 1:n) = [-2d0, 0d0, 0d0, 10d0, -1d0, 0d0]
46 7 : a(5, 1:n) = [-1d0, 0d0, 0d0, -5d0, 1d0, -1d0]
47 7 : a(6, 1:n) = [-1d0, -2d0, 0d0, 0d0, 0d0, 6d0]
48 :
49 7 : b(1, 1:n) = [10d0, 0d0, 0d0, 0d0, 0d0, 0d0]
50 7 : b(2, 1:n) = [-2d0, 12d0, -3d0, -1d0, 0d0, 0d0]
51 7 : b(3, 1:n) = [0d0, 1d0, 15d0, 0d0, 0d0, 0d0]
52 7 : b(4, 1:n) = [0d0, 0d0, 0d0, 10d0, -1d0, 0d0]
53 7 : b(5, 1:n) = [0d0, 0d0, 0d0, -5d0, 1d0, -1d0]
54 7 : b(6, 1:n) = [0d0, 0d0, 0d0, 0d0, 0d0, 6d0]
55 :
56 1 : ierr = 0
57 :
58 1 : write (*, *) 'dense_to_row_sparse'
59 1 : call dense_to_row_sparse(n, ndim, a, nzmax, nz, iptr, jind, values, ierr)
60 1 : if (ierr /= 0) call mesa_error(__FILE__, __LINE__)
61 43 : a2 = -1
62 1 : write (*, *) 'find_loc_in_row_sparse'
63 7 : do i = 1, n
64 6 : hint = 0
65 22 : do k = iptr(i), iptr(i + 1) - 1
66 15 : j = jind(k)
67 15 : call find_loc_in_sparse(compressed_row_sparse, n, nzmax, iptr, jind, i, j, hint, kk, ierr)
68 15 : if (kk /= k .or. ierr /= 0) then
69 0 : write (*, *) 'failure in find_loc_in_row_sparse', i, j, k, kk
70 0 : call mesa_error(__FILE__, __LINE__)
71 : end if
72 21 : hint = k
73 : end do
74 : end do
75 :
76 1 : write (*, *) 'dense_to_band'
77 1 : call dense_to_band(n, ndim, b, ml, mu, bb, ldbb, ierr)
78 1 : if (ierr /= 0) call mesa_error(__FILE__, __LINE__)
79 :
80 1 : write (*, *) 'band_to_dense'
81 43 : a2 = -1
82 1 : call band_to_dense(n, ml, mu, bb, ldbb, ndim, a2, ierr)
83 1 : if (ierr /= 0) call mesa_error(__FILE__, __LINE__)
84 43 : if (any(b /= a2)) call mesa_error(__FILE__, __LINE__)
85 :
86 1 : write (*, *) 'band_to_column_sparse'
87 1 : call band_to_column_sparse(n, ml, mu, bb, ldbb, nzmax, nz, iptr, jind, values, ierr)
88 1 : if (ierr /= 0) call mesa_error(__FILE__, __LINE__)
89 :
90 1 : write (*, *) 'column_sparse_to_band'
91 37 : bb2 = -1
92 1 : call column_sparse_to_band(n, ml, mu, bb2, ldbb, nz, iptr, jind, values, ierr)
93 1 : if (ierr /= 0) call mesa_error(__FILE__, __LINE__)
94 37 : if (any(bb /= bb2)) call mesa_error(__FILE__, __LINE__)
95 :
96 1 : write (*, *) 'band_to_row_sparse'
97 1 : call band_to_row_sparse(n, ml, mu, bb, ldbb, nzmax, nz, iptr, jind, values, ierr)
98 1 : if (ierr /= 0) call mesa_error(__FILE__, __LINE__)
99 :
100 1 : write (*, *) 'okay'
101 1 : write (*, *)
102 :
103 1 : end subroutine test_format_conversion
104 :
105 : end module test_mtx_support
|