File: generic_clapack.pyf

package info (click to toggle)
python-scipy 0.6.0-12
  • links: PTS, VCS
  • area: main
  • in suites: lenny
  • size: 32,016 kB
  • ctags: 46,675
  • sloc: cpp: 124,854; ansic: 110,614; python: 108,664; fortran: 76,260; objc: 424; makefile: 384; sh: 10
file content (291 lines) | stat: -rw-r--r-- 12,626 bytes parent folder | download
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
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
!%f90 -*- f90 -*-
! Signatures for f2py-wrappers of ATLAS C LAPACK functions.
!
! Author: Pearu Peterson
! Created: Jan-Feb 2002
! $Revision: 905 $ $Date: 2003-11-27 07:28:01 -0800 (Thu, 27 Nov 2003) $
!
! Usage:
!   f2py -c clapack.pyf -L/usr/local/lib/atlas -llapack -lf77blas -lcblas -latlas -lg2c

python module generic_clapack

interface

   function <tchar=s,d,c,z>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_<tchar=s,d,c,z>gesv
     integer intent(c,hide) ::  <tchar=s,d,c,z>gesv
     callstatement <tchar=s,d,c,z>gesv_return_value = info = (*f2py_func)(102-rowmajor,n,nrhs,a,n,piv,b,n)
     callprotoargument const int,const int,const int,<type_in_c>*,const int,int*,<type_in_c>*,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)
     <type_in> dimension(n,n),check(shape(a,0)==shape(a,1)) :: a
     integer dimension(n),depend(n),intent(out) :: piv
     <type_in> 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 <tchar=s,d,c,z>gesv

   function <tchar=s,d,c,z>getrf(m,n,a,piv,info,rowmajor)

   ! lu,piv,info = getrf(a,rowmajor=1,overwrite_a=0)
   ! Compute an LU factorization of a  general  M-by-N  matrix  A.
   ! A * P = L * U
     threadsafe
     fortranname  clapack_<tchar=s,d,c,z>getrf
     integer intent(c,hide) ::  <tchar=s,d,c,z>getrf
     callstatement <tchar=s,d,c,z>getrf_return_value = info = (*f2py_func)(102-rowmajor,m,n,a,(rowmajor?n:m),piv)
     callprotoargument const int,const int,const int,<type_in_c>*,const int,int*

     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)
     <type_in> dimension(m,n),intent(c,in,out,copy,out=lu) :: a
     integer dimension((m<n?m:n)),depend(m,n),intent(out) :: piv
     integer intent(out):: info

   end function <tchar=s,d,c,z>getrf

   function <tchar=s,d,c,z>getrs(n,nrhs,lu,piv,b,info,trans,rowmajor)

   ! x,info = getrs(lu,piv,b,trans=0,rowmajor=1,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

     fortranname  clapack_<tchar=s,d,c,z>getrs
     integer intent(c,hide) ::  <tchar=s,d,c,z>getrs
     callstatement <tchar=s,d,c,z>getrs_return_value = info = (*f2py_func)(102-rowmajor,111+trans,n,nrhs,lu,n,piv,b,n)
     callprotoargument const int,const int,const int,const int,<type_in_c>*,const int,int*,<type_in_c>*,const int

     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)
     <type_in> dimension(n,n),intent(c,in) :: lu
     check(shape(lu,0)==shape(lu,1)) :: lu
     integer dimension(n),intent(in),depend(n) :: piv
     <type_in> 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 <tchar=s,d,c,z>getrs


   function <tchar=s,d,c,z>getri(n,lu,piv,info,rowmajor)

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

     fortranname  clapack_<tchar=s,d,c,z>getri
     integer intent(c,hide) ::  <tchar=s,d,c,z>getri
     callstatement <tchar=s,d,c,z>getri_return_value = info = (*f2py_func)(102-rowmajor,n,lu,n,piv)
     callprotoargument const int,const int,<type_in_c>*,const int,const int*

     integer optional,intent(in),check(rowmajor==1||rowmajor==0) :: rowmajor = 1

     integer depend(lu),intent(hide):: n = shape(lu,0)
     <type_in> dimension(n,n),intent(c,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

   end function <tchar=s,d,c,z>getri

   function <tchar=s,d,c,z>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_<tchar=s,d,c,z>posv
     integer intent(c,hide) ::  <tchar=s,d,c,z>posv
     callstatement <tchar=s,d,c,z>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,<type_in_c>*,const int,<type_in_c>*,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)
     <type_in> dimension(n,n),intent(c,in,out,copy,out=c) :: a
     check(shape(a,0)==shape(a,1)) :: a
     <type_in> 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 <tchar=s,d,c,z>posv
        
   function <tchar=s,d>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

     fortranname  clapack_<tchar=s,d>potrf
     integer intent(c,hide) ::  <tchar=s,d>potrf
     callstatement <tchar=s,d>potrf_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,a,n); if(clean){int i,j;if(lower){for(i=0;i<n;++i) for(j=i+1;j<n;++j) *(a+i*n+j)=0.0;} else {for(i=0;i<n;++i) for(j=i+1;j<n;++j) *(a+j*n+i)=0.0;}}
     callprotoargument const int,const int,const int,<type_in_c>*,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)
     <type_in> 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 <tchar=s,d>potrf

   function <tchar=c,z>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^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_<tchar=c,z>potrf
     integer intent(c,hide) ::  <tchar=c,z>potrf
     callstatement <tchar=c,z>potrf_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,a,n); if(clean){int i,j,k;if(lower){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;}} else {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;}}}
     callprotoargument const int,const int,const int,<type_in_c>*,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)
     <type_in> 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 <tchar=c,z>potrf
   
   
   function <tchar=s,d,c,z>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_<tchar=s,d,c,z>potrs
     integer intent(c,hide) ::  <tchar=s,d,c,z>potrs
     callstatement <tchar=s,d,c,z>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,<type_in_c>*,const int,<type_in_c>*,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)
     <type_in> dimension(n,n),intent(c,in) :: c
     check(shape(c,0)==shape(c,1)) :: c
     <type_in> 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 <tchar=s,d,c,z>potrs

   function <tchar=s,d,c,z>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_<tchar=s,d,c,z>potri
     integer intent(c,hide) ::  <tchar=s,d,c,z>potri
     callstatement <tchar=s,d,c,z>potri_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,c,n)
     callprotoargument const int,const int,const int,<type_in_c>*,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)
     <type_in> 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 <tchar=s,d,c,z>potri


   function <tchar=s,d,c,z>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_<tchar=s,d,c,z>lauum
     integer intent(c,hide) ::  <tchar=s,d,c,z>lauum
     callstatement <tchar=s,d,c,z>lauum_return_value = info = (*f2py_func)(102-rowmajor,121+lower,n,c,n)
     callprotoargument const int,const int,const int,<type_in_c>*,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)
     <type_in> 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 <tchar=s,d,c,z>lauum

   function <tchar=s,d,c,z>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_<tchar=s,d,c,z>trtri
     integer intent(c,hide) ::  <tchar=s,d,c,z>trtri
     callstatement <tchar=s,d,c,z>trtri_return_value = info = (*f2py_func)(102-rowmajor,121+lower,131+unitdiag,n,c,n)
     callprotoargument const int,const int,const int,const int,<type_in_c>*,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)
     <type_in> 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 <tchar=s,d,c,z>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/