1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207
|
!%f90 -*- f90 -*-
!
! Signatures for f2py wrappers of ATLAS LAPACK functions.
!
! gesv
! getrf
! getrs
! getri
! posv
! potrf
! potrs
! potri
! lauum
! trtri
!
python module clapack
interface
function <prefix>gesv(n,nrhs,a,piv,b,info,rowmajor)
! lu,piv,x,info = gesv(a,b,rowmajor=1,overwrite_a=0,overwrite_b=0)
! Solve A * X = B.
! A * P = L * U
! U is unit upper diagonal triangular, L is lower triangular,
! piv pivots columns.
fortranname clapack_<prefix>gesv
integer intent(c,hide) :: <prefix>gesv
callstatement <prefix>gesv_return_value = info = (*f2py_func)(102-rowmajor,n,nrhs,a,n,piv,b,n)
callprotoargument const int,const int,const int,<ctype>*,const int,int*,<ctype>*,const int
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)
<ftype> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
integer dimension(n),depend(n),intent(out) :: piv
<ftype> 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 <prefix>gesv
function <prefix>posv(n,nrhs,a,b,info,lower,rowmajor)
! c,x,info = posv(a,b,lower=0,rowmajor=1,overwrite_a=0,overwrite_b=0)
! 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_<prefix>posv
integer intent(c,hide) :: <prefix>posv
callstatement <prefix>posv_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,nrhs,a,n,b,n)
callprotoargument const int,const int,const int,const int,<ctype>*,const int,<ctype>*,const int
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)
<ftype> dimension(n,n),intent(c,in,out,copy,out=c) :: a
check(shape(a,0)==shape(a,1)) :: a
<ftype> 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 <prefix>posv
function <prefix>potrf(n,a,info,lower,clean,rowmajor)
! c,info = potrf(a,lower=0,clean=1,rowmajor=1,overwrite_a=0)
! 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.
! clean==1 zeros strictly lower or upper parts of U or L, respectively
! c,info = potrf(a,lower=0,clean=1,rowmajor=1,overwrite_a=0)
! Compute Cholesky decomposition of symmetric positive defined matrix:
! A = U^H * U, C = U if lower = 0
! A = L * L^H, C = L if lower = 1
! C is triangular matrix of the corresponding Cholesky decomposition.
! clean==1 zeros strictly lower or upper parts of U or L, respectively
fortranname clapack_<prefix>potrf
integer intent(c,hide) :: <prefix>potrf
! <_init1=*(a+i*n+j)=0.0;,\0,k=i*n+j;(a+k)-\>r=(a+k)-\>i=0.0;,\2>
! <_init2=*(a+j*n+i)=0.0;,\0,k=j*n+i;(a+k)-\>r=(a+k)-\>i=0.0;,\2>
callstatement <prefix>potrf_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,a,n); if(clean){int i,j<,,\,k,\2>;if(lower){for(i=0;i\<n;++i) for(j=i+1;j\<n;++j) {<_init1>}} else {for(i=0;i\<n;++i) for(j=i+1;j\<n;++j) {<_init2>}}}
callprotoargument const int,const int,const int,<ctype>*,const int
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(clean==0||clean==1) :: clean = 1
integer depend(a),intent(hide):: n = shape(a,0)
<ftype> 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 <prefix>potrf
function <prefix>potrs(n,nrhs,c,b,info,lower,rowmajor)
! x,info = potrs(c,b,lower=0,rowmajor=1,overwrite_b=0)
! 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_<prefix>potrs
integer intent(c,hide) :: <prefix>potrs
callstatement <prefix>potrs_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,nrhs,c,n,b,n)
callprotoargument const int,const int,const int,const int,<ctype>*,const int,<ctype>*,const int
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)
<ftype> dimension(n,n),intent(c,in) :: c
check(shape(c,0)==shape(c,1)) :: c
<ftype> 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 <prefix>potrs
function <prefix>potri(n,c,info,lower,rowmajor)
! inv_a,info = potri(c,lower=0,rowmajor=1,overwrite_c=0)
! 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_<prefix>potri
integer intent(c,hide) :: <prefix>potri
callstatement <prefix>potri_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,c,n)
callprotoargument const int,const int,const int,<ctype>*,const int
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)
<ftype> dimension(n,n),intent(c,in,out,copy,out=inv_a) :: c
check(shape(c,0)==shape(c,1)) :: c
integer intent(out) :: info
end function <prefix>potri
function <prefix>lauum(n,c,info,lower,rowmajor)
! a,info = lauum(c,lower=0,rowmajor=1,overwrite_c=0)
! 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_<prefix>lauum
integer intent(c,hide) :: <prefix>lauum
callstatement <prefix>lauum_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,c,n)
callprotoargument const int,const int,const int,<ctype>*,const int
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)
<ftype> 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 <prefix>lauum
function <prefix>trtri(n,c,info,lower,unitdiag,rowmajor)
! inv_c,info = trtri(c,lower=0,unitdiag=0,rowmajor=1,overwrite_c=0)
! 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_<prefix>trtri
integer intent(c,hide) :: <prefix>trtri
callstatement <prefix>trtri_return_value = info = (*f2py_func)(102-rowmajor,121+lower,131+unitdiag,n,c,n)
callprotoargument const int,const int,const int,const int,<ctype>*,const int
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)
<ftype> dimension(n,n),intent(c,in,out,copy,out=inv_c) :: c
check(shape(c,0)==shape(c,1)) :: c
integer intent(out) :: info
end function <prefix>trtri
end interface
end python module clapack
|