SF_LINALG

Description

SciFortran module for linear algebra

Quick access

Variables:

operator(.x.), operator(.px.)

Routines:

eig(), eigh(), p_eigh(), eigh_jacobi(), eigvals(), eigvalsh(), svdvals(), svd(), inv(), p_inv(), inv_sym(), inv_her(), inv_triang(), inv_gj(), solve(), lstsq(), inv_tridiag(), check_tridiag(), get_tridiag(), build_tridiag(), det(), deye(), zeye(), eye(), diag(), diagonal(), trace(), zeros(), ones(), kron(), kronecker_product(), operator(.kx.)(), outerprod(), cross_product(), s3_product(), mat_product(), p_mat_product(), deye_tridiag(), zeye_tridiag()

Used modules

  • sf_blacs

Variables

  • sf_linalg/operator (px) [public]
  • sf_linalg/operator (x) [public]

Subroutines and functions

interface  sf_linalg/eig(a, eval, evec[, jobvl, jobvr])
Parameters:
  • a (•, •) [real, complex, in]

  • eval (•) [complex, out]

  • evec (•, •) [complex, out]

Options:
  • jobvl [character(len=1)]

  • jobvr [character(len=1)]

interface  sf_linalg/eigh(am, bm, lam, c, a, e, d, u[, method, jobz, uplo, vl, vu, il, iu, tol, ev, irange, vrange])
Parameters:
  • am (•, •) [real, complex, in] – LHS matrix: Am c = lam Bm c

  • bm (•, •) [real, complex, in] – Bmt temporaries overwritten by dsygvd

  • lam (•) [real, out] – eigenvalues: Am c = lam Bm c

  • c (•, •) [real, complex, out] – eigenvectors: Am c = lam Bm c; c(i,j) = ith component of jth vec.

  • a (•, •) [real, complex, inout/in,required] – M v = E v/v(i,j) = ith component of jth vec.

  • e (size(a, 2)) [real, inout] – eigenvalues

  • d (•) [real, in,required]

  • u (max(1, -1 + size(d))) [real]

Options:
  • method [character(len=*)]

  • jobz [character(len=1)]

  • uplo [character(len=1)]

  • vl [real]

  • vu [real]

  • il [integer]

  • iu [integer]

  • tol [real]

  • ev (•, •) [real]

  • irange (2) [integer]

  • vrange (2) [integer]

interface  sf_linalg/p_eigh(a, w, nblock[, method, jobz, uplo, vl, vu, il, iu, tol])
Parameters:
  • a (•, •) [real, complex, inout/in,required] – M v = E v/v(i,j) = ith component of jth vec.

  • w (size(a, 2)) [real, inout] – eigenvalues

  • nblock [integer]

Options:
  • method [character(len=*)]

  • jobz [character(len=1)]

  • uplo [character(len=1)]

  • vl [real]

  • vu [real]

  • il [integer]

  • iu [integer]

  • tol [real]

interface  sf_linalg/eigh_jacobi(a, d, v, nrot, u, sweep)
Parameters:
  • a (•, •) [real, complex, inout]

  • d (size(a, 1)) [real, out]

  • v (size(a, 1), size(a, 2)) [real, out]

  • nrot [integer, out]

  • u (size(a, 1), size(a, 2)) [complex, out]

  • sweep [integer]

interface  sf_linalg/eigvals(a)
Parameters:

a (•, •) [real, complex, in]

interface  sf_linalg/eigvalsh(a)
Parameters:

a (•, •) [real, complex, in]

interface  sf_linalg/svdvals(a)
Parameters:

a (•, •) [real, complex, in]

interface  sf_linalg/svd(a, s, u, vtransp)
Parameters:
  • a (•, •) [real, complex, in] – use a temporary as dgesvd destroys its input

  • s (•) [real, out] – has size min(m, n) –> sigma matrix is (n x m) with sigma_ii = s_i

  • u (•, •) [real, complex, out]

  • vtransp (•, •) [real, complex, out]

interface  sf_linalg/inv(am)
Parameters:

am (•, •) [real, complex, inout] – matrix to be inverted

interface  sf_linalg/p_inv(a, nblock)
Parameters:
  • a (•, •) [real, complex, inout]

  • nblock [integer]

interface  sf_linalg/inv_sym(a[, uplo])
Parameters:

a (•, •) [real, complex]

Options:

uplo [character(len=*)]

interface  sf_linalg/inv_her(a[, uplo])
Parameters:

a (•, •) [complex]

Options:

uplo [character(len=*)]

interface  sf_linalg/inv_triang(a[, uplo, diag])
Parameters:

a (•, •) [real, complex]

Options:
  • uplo [character(len=1)]

  • diag [character(len=1)] – not a unit triangular matrix

interface  sf_linalg/inv_gj(a)
Parameters:

a (•, •) [real, complex, inout] – cmplx

interface  sf_linalg/solve(a, b[, trans])
Parameters:
  • a (•, •) [real, complex, in]

  • b (various shapes) [real, complex, inout]

Options:

trans [character(len=1)]

interface  sf_linalg/lstsq(a, b)
Parameters:
  • a (•, •) [real, complex, in]

  • b (•) [real, complex, in] – only one right-hand side

interface  sf_linalg/inv_tridiag(sub_, diag_, over_, inv, ainv, amat, nb, n[, n, nb])
Options:
  • n [integer, optional/default=1 + shape(sub, 0)]

  • nb [integer, optional/default=1 + shape(sub, 0)]

Parameters:
  • sub_ (various shapes) [real, complex]

  • diag_ (various shapes) [real, complex]

  • over_ (various shapes) [real, complex]

  • inv (n) [real, complex]

  • ainv (nb, n, n) [real, complex]

  • amat (•, •) [real, complex, inout]

interface  sf_linalg/check_tridiag(amat, nblock, nsize)
Parameters:
  • amat (•, •) [real, complex]

  • nblock [integer, in,required]

  • nsize [integer, in,required]

interface  sf_linalg/get_tridiag(amat, sub, diag, nblock, nsize[, over])
Parameters:
  • amat (•, •) [real, complex, in,required]

  • sub (various shapes) [real, complex]

  • diag (various shapes) [real, complex]

  • nblock [integer, in,required/default=1 + shape(sub, 0)]

  • nsize [integer, in,required]

Options:

over (various shapes) [real, complex]

interface  sf_linalg/build_tridiag(sub, diag[, over, nblock, nsize])
Parameters:
  • sub (various shapes) [real, complex]

  • diag (various shapes) [real, complex, in,required]

Options:
  • over (various shapes) [real, complex]

  • nblock [integer, optional/default=1 + shape(sub, 0)]

  • nsize [integer]

interface  sf_linalg/det(a)
Parameters:

a (•, •) [real, complex, in]

interface  sf_linalg/deye(n, i, j)
Parameters:
  • n [integer, in]

  • i [integer, in]

  • j [integer, in]

interface  sf_linalg/zeye(n, i, j)
Parameters:
  • n [integer, in]

  • i [integer, in]

  • j [integer, in]

interface  sf_linalg/eye(n, i, j)
Parameters:
  • n [integer, in]

  • i [integer, in]

  • j [integer, in]

interface  sf_linalg/diag(x)
Parameters:

x (•) [real, complex, in]

interface  sf_linalg/diagonal(a)
Parameters:

a (•, •) [real, complex, in]

interface  sf_linalg/trace(a)
Parameters:

a (•, •) [real, complex, in]

interface  sf_linalg/zeros(n, n1, n2, n3, n4, n5, n6, n7)
Parameters:
  • n [integer, in]

  • n1 [integer, in]

  • n2 [integer, in]

  • n3 [integer, in]

  • n4 [integer, in]

  • n5 [integer, in]

  • n6 [integer, in]

  • n7 [integer, in]

interface  sf_linalg/ones(n, n1, n2, n3, n4, n5, n6, n7)
Parameters:
  • n [integer, in]

  • n1 [integer, in]

  • n2 [integer, in]

  • n3 [integer, in]

  • n4 [integer, in]

  • n5 [integer, in]

  • n6 [integer, in]

  • n7 [integer, in]

interface  sf_linalg/kron(a, b)
Parameters:
  • a (•, •) [integer, real, complex, in]

  • b (•, •) [integer, real, complex, in]

interface  sf_linalg/kronecker_product(a, b)
Parameters:
  • a (•, •) [integer, real, complex, in]

  • b (•, •) [integer, real, complex, in]

interface  sf_linalg/operator(.kx.)(a, b)
Parameters:
  • a (•, •) [integer, real, complex, in]

  • b (•, •) [integer, real, complex, in]

interface  sf_linalg/outerprod(a, b)
Parameters:
  • a (•) [real, complex, in]

  • b (•) [real, complex, in]

interface  sf_linalg/cross_product(a, b)
Parameters:
  • a (3) [real, complex]

  • b (3) [real, complex]

interface  sf_linalg/s3_product(a, b, c)
Parameters:
  • a (3) [real, complex, in]

  • b (3) [real, complex, in]

  • c (3) [real, complex, in]

interface  sf_linalg/mat_product(a, b, c[, alfa, beta])
Parameters:
  • a (•, •) [real, complex, inout] – [N,K]

  • b (•, •) [real, complex, inout] – [K,M]

  • c (•, •) [real, complex, inout] – [N,M]

Options:
  • alfa [real, complex]

  • beta [real, complex]

interface  sf_linalg/p_mat_product(a, b, c, nblock[, alfa, beta])
Parameters:
  • a (•, •) [real, complex, inout] – [N,K]

  • b (•, •) [real, complex, inout] – [K,M]

  • c (•, •) [real, complex, inout] – [N,M]

  • nblock [integer]

Options:
  • alfa [real]

  • beta [real]

function  sf_linalg/deye_tridiag(nblock, n)
Parameters:
  • nblock [integer]

  • n [integer]

Return:

eye_block (nblock, n, n) [real]

function  sf_linalg/zeye_tridiag(nblock, n)
Parameters:
  • nblock [integer]

  • n [integer]

Return:

eye_block (nblock, n, n) [complex]