linalg 1.8.2
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
sparskit.f90
1
3module sparskit
4 implicit none
5
6 ! BLASSM.F
7 interface
8
31 subroutine amub(nrow, ncol, job, a, ja, ia, b, jb, ib, c, jc, ic, &
32 nzmax, iw, ierr)
33 use iso_fortran_env, only : int32, real64
34 integer(int32), intent(in) :: nrow, ncol, job, nzmax
35 integer(int32), intent(in) :: ja(*), ia(nrow+1), jb(*), ib(*)
36 integer(int32), intent(out) :: jc(*), ic(*), iw(ncol), ierr
37 real(real64), intent(in) :: a(*), b(*)
38 real(real64), intent(out) :: c(*)
39 end subroutine
40
63 subroutine aplb(nrow, ncol, job, a, ja, ia, b, jb, ib, c, jc, ic, &
64 nzmax, iw, ierr)
65 use iso_fortran_env, only : int32, real64
66 integer(int32), intent(in) :: nrow, ncol, job, nzmax
67 integer(int32), intent(in) :: ja(*), ia(nrow+1), jb(*), ib(nrow+1)
68 real(real64), intent(in) :: a(*), b(*)
69 real(real64), intent(out) :: c(*)
70 integer(int32), intent(out) :: jc(*), ic(nrow+1), iw(ncol), ierr
71 end subroutine
72
96 subroutine aplsb(nrow, ncol, a, ja, ia, s, b, jb, ib, c, jc, ic, &
97 nzmax, iw, ierr)
98 use iso_fortran_env, only : int32, real64
99 integer(int32), intent(in) :: nrow, ncol, nzmax
100 integer(int32), intent(in) :: ja(*), ia(nrow+1), jb(*), ib(nrow+1)
101 integer(int32), intent(out) :: ierr
102 integer(int32), intent(out) :: jc(*), ic(nrow+1), iw(ncol)
103 real(real64), intent(in) :: s
104 real(real64), intent(in) :: a(*), b(*)
105 real(real64), intent(out) :: c(*)
106 end subroutine
107 end interface
108
109 ! FORMATS.F
110 interface
111
124 subroutine csrcsc2(n, n2, job, ipos, a, ja, ia, a0, ja0, ia0)
125 use iso_fortran_env, only : int32, real64
126 integer(int32), intent(in) :: n, n2, job, ipos, ja(*), ia(n+1)
127 integer(int32), intent(out) :: ja0(*), ia0(n2+1)
128 real(real64), intent(in) :: a(*)
129 real(real64), intent(out) :: a0(*)
130 end subroutine
131
150 subroutine bndcsr(n, abd, nabd, lowd, ml, mu, a, ja, ia, len, ierr)
151 use iso_fortran_env, only : int32, real64
152 integer(int32), intent(in) :: n, nabd, lowd, ml, mu, len
153 real(real64), intent(in) :: abd(nabd,*)
154 real(real64), intent(out) :: a(*)
155 integer(int32), intent(out) :: ia(n+1), ja(*), ierr
156 end subroutine
157
170 subroutine csrmsr(n, a, ja, ia, ao, jao, wk, iwk)
171 use iso_fortran_env, only : int32, real64
172 integer(int32), intent(in) :: n
173 integer(int32), intent(inout) :: ja(*), ia(n+1)
174 integer(int32), intent(out) :: jao(*), iwk(n+1)
175 real(real64), intent(inout) :: a(*)
176 real(real64), intent(out) :: ao(*), wk(n)
177 end subroutine
178
191 subroutine msrcsr(n, a, ja, ao, jao, iao, wk, iwk)
192 use iso_fortran_env, only : int32, real64
193 integer(int32), intent(in) :: n, ja(*)
194 integer(int32), intent(out) :: jao(*), iao(n+1), iwk(n+1)
195 real(real64), intent(in) :: a(*)
196 real(real64), intent(out) :: ao(*), wk(n)
197 end subroutine
198
212 subroutine coocsr(nrow, nnz, a, ir, jc, ao, jao, iao)
213 use iso_fortran_env, only : int32, real64
214 integer(int32), intent(in) :: nrow, nnz, jc(*)
215 integer(int32), intent(inout) :: ir(*)
216 real(real64), intent(in) :: a(*)
217 integer(int32), intent(out) :: jao(*), iao(*)
218 real(real64) :: ao(*)
219 end subroutine
220 end interface
221
222 ! UNARY.F
223 interface
224
236 function getelm(i, j, a, ja, ia, iadd, sorted) result(rst)
237 use iso_fortran_env, only : int32, real64
238 integer(int32), intent(in) :: i, j, ia(*), ja(*)
239 real(real64), intent(in) :: a(*)
240 integer(int32), intent(out) :: iadd
241 logical, intent(in) :: sorted
242 real(real64) :: rst
243 end function
244
256 subroutine getdia(nrow, ncol, job, a, ja, ia, len, diag, idiag, ioff)
257 use iso_fortran_env, only : int32, real64
258 integer(int32), intent(in) :: nrow, ncol, job, ja(*), ia(*), ioff
259 integer(int32), intent(out) :: len, idiag(*)
260 real(real64), intent(in) :: a(*)
261 real(real64), intent(out) :: diag(*)
262 end subroutine
263
273 subroutine csort(n, a, ja, ia, values)
274 use iso_fortran_env, only : int32, real64
275 integer(int32), intent(in) :: n
276 real(real64), intent(inout) :: a(*)
277 integer(int32), intent(inout) :: ja(*)
278 integer(int32), intent(in) :: ia(*)
279 logical, intent(in) :: values
280 end subroutine
281
299
300 subroutine clncsr(job, value2, nrow, a, ja, ia, indu, iwk)
301 use iso_fortran_env, only : int32, real64
302 integer(int32), intent(in) :: job, value2, nrow
303 real(real64), intent(inout) :: a(*)
304 integer(int32), intent(inout) :: ja(*), ia(*)
305 integer(int32), intent(inout) :: indu(*), iwk(*)
306 end subroutine
307 end interface
308
309 ! ILUT.F
310 interface
311
343 subroutine ilut(n, a, ja, ia, lfil, droptol, alu, jlu, ju, iwk, w, &
344 jw, ierr)
345 use iso_fortran_env, only : int32, real64
346 integer(int32), intent(in) :: n, ja(*), ia(n+1), lfil, iwk
347 integer(int32), intent(out) :: jlu(*), ju(n), jw(2*n), ierr
348 real(real64), intent(in) :: a(*), droptol
349 real(real64), intent(out) :: alu(*), w(n+1)
350 end subroutine
351
409 subroutine ilutp(n, a, ja, ia, lfil, droptol, permtol, mbloc, alu, &
410 jlu, ju, iwk, w, jw, iperm, ierr)
411 use iso_fortran_env, only : int32, real64
412 integer(int32), intent(in) :: n, lfil, iwk, mbloc
413 integer(int32), intent(inout) :: ja(*), ia(n+1)
414 integer(int32), intent(out) :: jlu(*), ju(n), jw(2*n), &
415 iperm(2*n), ierr
416 real(real64), intent(in) :: droptol, permtol
417 real(real64), intent(inout) :: a(*)
418 real(real64), intent(out) :: alu(*), w(n+1)
419 end subroutine
420
452 subroutine ilud(n, a, ja, ia, alph, tol, alu, jlu, ju, iwk, w, jw, ierr)
453 use iso_fortran_env, only : int32, real64
454 integer(int32), intent(in) :: n, iwk, ja(*), ia(n+1)
455 integer(int32), intent(out) :: jlu(*), ju(n), jw(2*n), ierr
456 real(real64), intent(in) :: a(*), alph, tol
457 real(real64), intent(out) :: alu(*), w(2*n)
458 end subroutine
459
499 subroutine iludp(n, a, ja, ia, alph, droptol, permtol, mbloc, alu, &
500 jlu, ju, iwk, w, jw, iperm, ierr)
501 use iso_fortran_env, only : int32, real64
502 integer(int32), intent(in) :: n, iwk, mbloc
503 integer(int32), intent(inout) :: ja(*), ia(n+1)
504 integer(int32), intent(out) :: jlu(*), ju(n), jw(2*n), iperm(2*n), &
505 ierr
506 real(real64), intent(in) :: alph, droptol, permtol
507 real(real64), intent(inout) :: a(*)
508 real(real64), intent(out) :: alu(*), w(2*n)
509 end subroutine
510
540 subroutine pgmres(n, im, rhs, sol, vv, eps, maxits, iout, aa, ja, ia, &
541 alu, jlu, ju, ierr)
542 use iso_fortran_env, only : int32, real64
543 integer(int32), intent(in) :: n, im, maxits, iout, ja(*), ia(n+1), &
544 jlu(*), ju(n)
545 integer(int32), intent(out) :: ierr
546 real(real64), intent(in) :: aa(*), eps, alu(*)
547 real(real64), intent(inout) :: rhs(n), sol(n)
548 real(real64), intent(out) :: vv(n,*)
549 end subroutine
550
559 subroutine lusol(n, y, x, alu, jlu, ju)
560 use iso_fortran_env, only : int32, real64
561 integer(int32), intent(in) :: n, jlu(*), ju(*)
562 real(real64), intent(in) :: y(n), alu(*)
563 real(real64), intent(out) :: x(n)
564 end subroutine
565 end interface
566end module
Computes the matrix product: C = A * B.
Definition sparskit.f90:31
Computes the matrix sum: C = A + B, where the matrices are given in CSR format.
Definition sparskit.f90:63
Computes the matrix sum: C = A + s * B, where the matrices are given in CSR format.
Definition sparskit.f90:96
Converts the LINPACK, BLAS, LAPACK banded matrix format into a CSR format.
Definition sparskit.f90:150
@breif Cleans up a CSR matrix.
Definition sparskit.f90:300
Converte a matrix stored in coordinate format to CSR format.
Definition sparskit.f90:212
Sorces the elements of a CSR matrix in increasing order of their column indices within each row.
Definition sparskit.f90:273
Converts a CSR matrix into a CSC matrix (transposition).
Definition sparskit.f90:124
Converts a CSR matrix to an MSR matrix.
Definition sparskit.f90:170
Extracts the diagonal from a matrix.
Definition sparskit.f90:256
Gets element A(i,j) of matrix A for any pair (i,j).
Definition sparskit.f90:236
Computes the incomplete LU factorization of a sparse matrix in CSR format with standard dropping stra...
Definition sparskit.f90:452
Computes the incomplete LU factorization of a sparse matrix in CSR format with standard dropping stra...
Definition sparskit.f90:499
Computes the incomplete LU factorization of a sparse matrix in CSR format using a dual truncation mec...
Definition sparskit.f90:343
Computes the incomplete LU factorization of a sparse matrix in CSR format using a dual truncation mec...
Definition sparskit.f90:409
Solves the LU-factored system (LU) x = y.
Definition sparskit.f90:559
Converts and MSR matrix to a CSR matrix.
Definition sparskit.f90:191
An ILUT preconditioned GMRES algorithm. This routine utilizes the L and U matrices generated by the I...
Definition sparskit.f90:540
An interface to the SPARSKIT library available at https://www-users.cse.umn.edu/~saad/software/SPARSK...
Definition sparskit.f90:3