!%f90 -*- f90 -*- python module generic_clapack interface function gesv(n,nrhs,A,piv,B,info,rowmajor) ! LU,piv,X,info = gesv(A,B,rowmajor=1,copy_A=1,copy_B=1) ! Solve A * X = B. ! A * P = L * U ! U is unit upper diagonal triangular, L is lower triangular, ! piv pivots columns. fortranname clapack_gesv integer intent(c,hide) :: gesv callstatement info = (*f2py_func)(102-rowmajor,n,nrhs,A,n,piv,B,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer depend(A),intent(hide):: n = shape(A,0) integer depend(B),intent(hide):: nrhs = shape(B,1) dimension(n,n),check(shape(A,0)==shape(A,1)) :: A integer dimension(n),depend(n),intent(out) :: piv dimension(n,nrhs),check(shape(A,0)==shape(B,0)),depend(n) :: B integer intent(out)::info intent(in,out,copy,out=X) B intent(c,in,out,copy,out=LU) A end function gesv function getrf(m,n,A,piv,info,rowmajor) ! LU,piv,info = getrf(A,rowmajor=1,copy_A=1) ! Compute an LU factorization of a general M-by-N matrix A. ! A * P = L * U fortranname clapack_getrf integer intent(c,hide) :: getrf callstatement info = (*f2py_func)(102-rowmajor,m,n,A,(rowmajor?n:m),piv) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer depend(A),intent(hide):: m = shape(A,0) integer depend(A),intent(hide):: n = shape(A,1) dimension(m,n),intent(c,in,out,copy,out=LU) :: A integer dimension((mgetrf function getrs(n,nrhs,LU,piv,B,info,trans,rowmajor) ! X,info = getrs(LU,piv,B,trans=0,rowmajor=1,copy_B=1) ! Solve A * X = B if trans=0 ! Solve A^T * X = B if trans=1 ! Solve A^H * X = B if trans=2 ! A * P = L * U fortranname clapack_getrs integer intent(c,hide) :: getrs callstatement info = (*f2py_func)(102-rowmajor,110+trans,n,nrhs,LU,n,piv,B,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(trans>=0 && trans <=2) :: trans = 0 integer depend(LU),intent(hide):: n = shape(LU,0) integer depend(B),intent(hide):: nrhs = shape(B,1) dimension(n,n),intent(c,in) :: LU check(shape(LU,0)==shape(LU,1)) :: LU integer dimension(n),intent(in),depend(n) :: piv dimension(n,nrhs),intent(in,out,copy,out=X),depend(n),check(shape(LU,0)==shape(B,0)) :: B integer intent(out):: info end function getrs function getri(n,LU,piv,info,rowmajor) ! invA,info = getri(LU,piv,rowmajor=1,copy_LU=1) ! Find A inverse A^-1. ! A * P = L * U fortranname clapack_getri integer intent(c,hide) :: getri callstatement info = (*f2py_func)(102-rowmajor,n,LU,n,piv) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer depend(LU),intent(hide):: n = shape(LU,0) dimension(n,n),intent(c,in,out,copy,out=invA) :: LU check(shape(LU,0)==shape(LU,1)) :: LU integer dimension(n),intent(in),depend(n) :: piv integer intent(out):: info end function getri function posv(n,nrhs,A,B,info,lower,rowmajor) ! C,X,info = posv(A,B,lower=0,rowmajor=1,copy_A=1,copy_B=1) ! Solve A * X = B. ! A is symmetric positive defined ! A = U^T * U, C = U if lower = 0 ! A = L * L^T, C = L if lower = 1 ! C is triangular matrix of the corresponding Cholesky decomposition. fortranname clapack_posv integer intent(c,hide) :: posv callstatement info = (*f2py_func)(102-rowmajor,121+lower,n,nrhs,A,n,B,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(lower==0||lower==1) :: lower = 0 integer depend(A),intent(hide):: n = shape(A,0) integer depend(B),intent(hide):: nrhs = shape(B,1) dimension(n,n),intent(c,in,out,copy,out=C) :: A check(shape(A,0)==shape(A,1)) :: A dimension(n,nrhs),intent(in,out,copy,out=X),depend(n):: B check(shape(A,0)==shape(B,0)) :: B integer intent(out) :: info end function posv function potrf(n,A,info,lower,rowmajor) ! C,info = potrf(A,lower=0,rowmajor=1,copy_A=1) ! Compute Cholesky decomposition of symmetric positive defined matrix: ! A = U^T * U, C = U if lower = 0 ! A = L * L^T, C = L if lower = 1 ! C is triangular matrix of the corresponding Cholesky decomposition. fortranname clapack_potrf integer intent(c,hide) :: potrf callstatement info = (*f2py_func)(102-rowmajor,121+lower,n,A,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(lower==0||lower==1) :: lower = 0 integer depend(A),intent(hide):: n = shape(A,0) dimension(n,n),intent(c,in,out,copy,out=C) :: A check(shape(A,0)==shape(A,1)) :: A integer intent(out) :: info end function potrf function potrs(n,nrhs,C,B,info,lower,rowmajor) ! X,info = potrs(C,B,lower=0,rowmajor=1,copy_B=1) ! Solve A * X = B. ! A is symmetric positive defined ! A = U^T * U, C = U if lower = 0 ! A = L * L^T, C = L if lower = 1 ! C is triangular matrix of the corresponding Cholesky decomposition. fortranname clapack_potrs integer intent(c,hide) :: potrs callstatement info = (*f2py_func)(102-rowmajor,121+lower,n,nrhs,C,n,B,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(lower==0||lower==1) :: lower = 0 integer depend(C),intent(hide):: n = shape(C,0) integer depend(B),intent(hide):: nrhs = shape(B,1) dimension(n,n),intent(c,in) :: C check(shape(C,0)==shape(C,1)) :: C dimension(n,nrhs),intent(in,out,copy,out=X),depend(n):: B check(shape(C,0)==shape(B,0)) :: B integer intent(out) :: info end function potrs function potri(n,C,info,lower,rowmajor) ! invA,info = potri(C,lower=0,rowmajor=1,copy_C=1) ! Compute A inverse A^-1. ! A = U^T * U, C = U if lower = 0 ! A = L * L^T, C = L if lower = 1 ! C is triangular matrix of the corresponding Cholesky decomposition. fortranname clapack_potri integer intent(c,hide) :: potri callstatement info = (*f2py_func)(102-rowmajor,121+lower,n,C,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(lower==0||lower==1) :: lower = 0 integer depend(C),intent(hide):: n = shape(C,0) dimension(n,n),intent(c,in,out,copy,out=invA) :: C check(shape(C,0)==shape(C,1)) :: C integer intent(out) :: info end function potri function lauum(n,C,info,lower,rowmajor) ! A,info = lauum(C,lower=0,rowmajor=1,copy_C=1) ! Compute product ! U^T * U, C = U if lower = 0 ! L * L^T, C = L if lower = 1 ! C is triangular matrix of the corresponding Cholesky decomposition. fortranname clapack_lauum integer intent(c,hide) :: lauum callstatement info = (*f2py_func)(102-rowmajor,121+lower,n,C,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(lower==0||lower==1) :: lower = 0 integer depend(C),intent(hide):: n = shape(C,0) dimension(n,n),intent(c,in,out,copy,out=A) :: C check(shape(C,0)==shape(C,1)) :: C integer intent(out) :: info end function lauum function trtri(n,C,info,lower,unitdiag,rowmajor) ! invC,info = trtri(C,lower=0,unitdiag=0,rowmajor=1,copy_C=1) ! Compute C inverse C^-1 where ! C = U if lower = 0 ! C = L if lower = 1 ! C is non-unit triangular matrix if unitdiag = 0 ! C is unit triangular matrix if unitdiag = 1 fortranname clapack_trtri integer intent(c,hide) :: trtri callstatement info = (*f2py_func)(102-rowmajor,121+lower,131+unitdiag,n,C,n) integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1 integer optional,intent(in),check(lower==0||lower==1) :: lower = 0 integer optional,intent(in),check(unitdiag==0||unitdiag==1) :: unitdiag = 0 integer depend(C),intent(hide):: n = shape(C,0) dimension(n,n),intent(c,in,out,copy,out=invC) :: C check(shape(C,0)==shape(C,1)) :: C integer intent(out) :: info end function trtri end interface end python module generic_clapack ! This file was auto-generated with f2py (version:2.10.173). ! See http://cens.ioc.ee/projects/f2py2e/