forked from scipy/scipy
-
Notifications
You must be signed in to change notification settings - Fork 0
/
flapack_pos_def.pyf.src
157 lines (127 loc) · 7.34 KB
/
flapack_pos_def.pyf.src
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
! Signatures for f2py-wrappers of FORTRAN LAPACK Positive Definite Matrix functions.
!
subroutine <prefix>posv(n,nrhs,a,b,info,lower)
! c,x,info = posv(a,b,lower=0,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.
callstatement (*f2py_func)((lower?"L":"U"),&n,&nrhs,a,&n,b,&n,&info)
callprotoargument char*,int*,int*,<ctype>*,int*,<ctype>*,int*,int*
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(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 subroutine <prefix>posv
subroutine <prefix>posvx(fact,n,nrhs,a,lda,af,ldaf,equed,s,b,ldb,x,ldx,rcond,ferr,berr,work,irwork,info,lower)
! Solve A * X = B for Symmetric/Hermitian A
! "expert" version of the ?POSV routines
threadsafe
callstatement (*f2py_func)(fact,(lower?"L":"U"),&n,&nrhs,a,&lda,af,&ldaf,equed,s,b,&ldb,x,&ldx,&rcond,ferr,berr,work,irwork,&info)
callprotoargument char*,char*,int*,int*,<ctype>*,int*,<ctype>*,int*,char*,<ctypereal>*,<ctype>*,int*,<ctype>*,int*,<ctypereal>*,<ctypereal>*,<ctypereal>*,<ctype>*,<int,int,float,double>*,int*
character optional,intent(in):: fact = "E"
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),check(shape(a,0)==shape(a,1)),intent(in,copy,out,out=a_s):: a
integer depend(a),intent(hide):: lda = shape(a,0)
<ftype> optional,intent(in,out,out=lu),dimension(n,n),depend(n):: af
integer depend(af),intent(hide):: ldaf = shape(af,0)
character optional,intent(in,out):: equed = "Y"
<ftypereal> optional,dimension(n),depend(n),intent(in,out):: s
<ftype> dimension(n,nrhs),check(shape(b,0)==n),depend(n),intent(in,copy,out,out=b_s):: b
integer depend(b),intent(hide):: ldb = shape(b,0)
<ftype> dimension(n,nrhs),depend(n,nrhs),intent(out):: x
integer depend(x),intent(hide):: ldx = shape(x,0)
<ftypereal> intent(out):: rcond
<ftypereal> intent(out),dimension(nrhs),depend(nrhs):: ferr
<ftypereal> intent(out),dimension(nrhs),depend(nrhs):: berr
<ftype> intent(hide),dimension(<3*n,3*n,2*n,2*n>),depend(n):: work
<integer,integer,real,double precision> intent(hide),dimension(n),depend(n):: irwork
integer intent(out):: info
end subroutine <prefix>posvx
subroutine <prefix>pocon(uplo,n,a,lda,anorm,rcond,work,irwork,info)
! Computes the 1- or inf- norm reciprocal condition number estimate
! for a positive definite symmetric/hermitian matrix.
threadsafe
callstatement (*f2py_func)(uplo,&n,a,&lda,&anorm,&rcond,work,irwork,&info)
callprotoargument char*,int*,<ctype>*,int*,<ctypereal>*,<ctypereal>*,<ctype>*,<int,int,float,double>*,int*
character optional,intent(in):: uplo = 'U'
integer depend(a),intent(hide):: n = shape(a,0)
<ftype> dimension(n,n),check(shape(a,0)==shape(a,1)),intent(in):: a
integer depend(a),intent(hide):: lda = shape(a,0)
<ftypereal> intent(in):: anorm
<ftypereal> intent(out):: rcond
<ftype> depend(n),dimension(<3*n,3*n,2*n,2*n>),intent(hide,cache):: work
<integer,integer,real, double precision> depend(n),dimension(n),intent(hide,cache):: irwork
integer intent(out):: info
end subroutine <prefix>pocon
subroutine <prefix2>potrf(n,a,lda,info,lower,clean)
! c,info = potrf(a,lower=0,clean=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
callstatement (*f2py_func)((lower?"L":"U"),&n,a,&lda,&info); if(clean){int i,j;if(lower){for(i=0;i\<n;++i) for(j=i+1;j\<n;++j) *(a+j*n+i)=0.0;} else {for(i=0;i\<n;++i) for(j=i+1;j<n;++j) *(a+i*n+j)=0.0;}}
callprotoargument char*,int*,<ctype2>*,int*,int*
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)
<ftype2> dimension(n,n),check(shape(a,0)==shape(a,1)),intent(in,out,copy,out=c) :: a
integer depend(n),intent(hide) :: lda = MAX(1,n)
integer intent(out) :: info
end subroutine <prefix2>potrf
subroutine <prefix2c>potrf(n,a,lda,info,lower,clean)
! c,info = potrf(a,lower=0,clean=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
callstatement (*f2py_func)((lower?"L":"U"),&n,a,&lda,&info); if(clean){int i,j,k;if(lower){for(i=0;i\<n;++i) for(j=i+1;j\<n;++j) {k=j*n+i;(a+k)->r=(a+k)->i=0.0;}} else {for(i=0;i\<n;++i) for(j=i+1;j\<n;++j) {k=i*n+j;(a+k)->r=(a+k)->i=0.0;}}}
callprotoargument char*,int*,<ctype2c>*,int*,int*
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)
<ftype2c> dimension(n,n),check(shape(a,0)==shape(a,1)),intent(in,out,copy,out=c) :: a
integer depend(n),intent(hide) :: lda = MAX(1,n)
integer intent(out) :: info
end subroutine <prefix2c>potrf
subroutine <prefix>potrs(n,nrhs,c,b,info,lower)
! x,info = potrs(c,b,lower=0=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.
callstatement (*f2py_func)((lower?"L":"U"),&n,&nrhs,c,&n,b,&n,&info)
callprotoargument char*,int*,int*,<ctype>*,int*,<ctype>*,int*,int*
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(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 subroutine <prefix>potrs
subroutine <prefix>potri(n,c,info,lower)
! inv_a,info = potri(c,lower=0,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.
callstatement (*f2py_func)((lower?"L":"U"),&n,c,&n,&info)
callprotoargument char*,int*,<ctype>*,int*,int*
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(in,out,copy,out=inv_a) :: c
check(shape(c,0)==shape(c,1)) :: c
integer intent(out) :: info
end subroutine <prefix>potri