File: flapack_lec.pyf.src

package info (click to toggle)
python-scipy 0.10.1%2Bdfsg2-1
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 42,232 kB
  • sloc: cpp: 224,773; ansic: 103,496; python: 85,210; fortran: 79,130; makefile: 272; sh: 43
file content (197 lines) | stat: -rw-r--r-- 8,080 bytes parent folder | download | duplicates (6)
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
! -*- f90 -*-
!
! Contains wrappers for the following LAPACK routines:
!
!  Computational routines for linear equations:
!
!   getrf (general, factorize)
!   getrs (general, solve)
!   gecon (general, estimate condition number) - Not Implemented
!   gerfs (general, error bounds) - Not Implemented
!   getri (general, invert)
!   geequ  (general, equilibrate) - Not Implemented
!   gbtrf, gbtrs, gbcon, gbrfs, gbequ (general band) - Not Implemented
!   gttrf, gttrs, gtcon, gtrfs (general tridiagonal) - Not Implemented
!   potrf (symmetric/Hermitian positive, factorize)
!   potrs (symmetric/Hermitian positive, solve)
!   pocon (symmetric/Hermitian positive, estimate condition number) - Not Implemented
!   porfs (symmetric/Hermitian positive, error bounds) - Not Implemented
!   potri (symmetric/Hermitian positive, invert)
!   poequ (symmetric/Hermitian positive, equilibrate) - Not Implemented
!   pptrf, pptrs, ppcon, pprfs, pptri, ppequ (symmetric/Hermitian positive packed storage) - Not Implemented
!   pbtrf, pbtrs, pbcon, pbrfs, pbequ (symmetric/Hermitian positive band) - Not Implemented
!   pttrf, pttrs, ptcon, ptrfs (symmetric/Hermitian positive tridiagonal) - Not Implemented
!   sytrf, sytrs, sycon, syrfs, sytri (symmetric indefinite) - Not Implemented
!   hetrf, hetrs, hecon, herfs, hetri (Hermitian indefinite) - Not Implemented
!   sptrf, sptrs, spcon, sprfs, sptri (symmetric indefinite packed storage) - Not Implemented
!   hptrf, hptrs, hpcon, hprfs, hptri (Hermitian indefinite packed storage) - Not Implemented
!   trtrs (triangular, solve) - Not Implemented
!   trcon (triangular, estimate condition number) - Not Implemented
!   trrfs (triangular, error bounds) - Not Implemented
!   trtri (triangular, invert)
!   tptrs, tpcon, tprfs, tptri (triangular packed storage) - Not Implemented
!   tbtrs, tbcon, tbrfs (triangular band) - Not Implemented

!
! Factorize
! =========

   subroutine <prefix>getrf(m,n,a,piv,info)

   ! lu,piv,info = getrf(a,overwrite_a=0)
   ! Compute an LU factorization of a  general  M-by-N  matrix  A.
   ! A = P * L * U
     threadsafe
     callstatement {int i;(*f2py_func)(&m,&n,a,&m,piv,&info);for(i=0,n=MIN(m,n);i\<n;--piv[i++]);}
     callprotoargument int*,int*,<ctype>*,int*,int*,int*

     integer depend(a),intent(hide):: m = shape(a,0)
     integer depend(a),intent(hide):: n = shape(a,1)
     <ftype> dimension(m,n),intent(in,out,copy,out=lu) :: a
     integer dimension(MIN(m,n)),depend(m,n),intent(out) :: piv
     integer intent(out):: info

   end subroutine <prefix>getrf

   subroutine <prefix>potrf(n,a,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


     ! <_def=,,\,k,\2>
     ! <_init1=*(a+j*n+i)=0.0;,\0,k=j*n+i;(*(a+k)).r=(*(a+k)).i=0.0;,\2>
     ! <_init2=*(a+i*n+j)=0.0;,\0,k=i*n+j;(*(a+k)).r=(*(a+k)).i=0.0;,\2>
     callstatement (*f2py_func)((lower?"L":"U"),&n,a,&n,&info); if(clean){int i,j<_def>;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 char*,int*,<ctype>*,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)
     <ftype> dimension(n,n),intent(in,out,copy,out=c) :: a
     check(shape(a,0)==shape(a,1)) :: a
     integer intent(out) :: info
     
   end subroutine <prefix>potrf


!
! Solve using factorization
! =========================


   subroutine <prefix>getrs(n,nrhs,lu,piv,b,info,trans)

   ! x,info = getrs(lu,piv,b,trans=0,overwrite_b=0)
   ! 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
     threadsafe
     callstatement {int i;for(i=0;i\<n;++piv[i++]);(*f2py_func)((trans?(trans==2?"C":"T"):"N"),&n,&nrhs,lu,&n,piv,b,&n,&info);for(i=0;i\<n;--piv[i++]);}
     callprotoargument char*,int*,int*,<ctype>*,int*,int*,<ctype>*,int*,int*

     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)
     <ftype> dimension(n,n),intent(in) :: lu
     check(shape(lu,0)==shape(lu,1)) :: lu
     integer dimension(n),intent(in),depend(n) :: piv
     <ftype> dimension(n,nrhs),intent(in,out,copy,out=x),depend(n),check(shape(lu,0)==shape(b,0)) :: b
     integer intent(out):: info
   end subroutine <prefix>getrs

   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

!
! Invert using factorization
! ==========================

   subroutine <prefix>getri(n,lu,piv,work,lwork,info)

   ! inv_a,info = getri(lu,piv,lwork=3*n,overwrite_lu=0)
   ! Find A inverse A^-1.
   ! A = P * L * U

     callstatement {int i;for(i=0;i\<n;++piv[i++]);(*f2py_func)(&n,lu,&n,piv,work,&lwork,&info);for(i=0;i\<n;--piv[i++]);}
     callprotoargument int*,<ctype>*,int*,int*,<ctype>*,int*,int*

     integer depend(lu),intent(hide):: n = shape(lu,0)
     <ftype> dimension(n,n),intent(in,out,copy,out=inv_a) :: lu
     check(shape(lu,0)==shape(lu,1)) :: lu
     integer dimension(n),intent(in),depend(n) :: piv
     integer intent(out):: info
     integer optional,intent(in),depend(n),check(lwork>=n) :: lwork=3*n
     <ftype> dimension(lwork),intent(hide,cache),depend(lwork) :: work

   end subroutine <prefix>getri

   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(c,in,out,copy,out=inv_a) :: c
     check(shape(c,0)==shape(c,1)) :: c
     integer intent(out) :: info
     
   end subroutine <prefix>potri

  subroutine <prefix>trtri(n,c,info,lower,unitdiag)
   
     ! inv_c,info = trtri(c,lower=0,unitdiag=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

     callstatement (*f2py_func)((lower?"L":"U"),(unitdiag?"U":"N"),&n,c,&n,&info)
     callprotoargument char*,char*,int*,<ctype>*,int*,int*

     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(in,out,copy,out=inv_c) :: c
     check(shape(c,0)==shape(c,1)) :: c
     integer intent(out) :: info
     
   end subroutine <prefix>trtri