File: flapack_le.pyf.src

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 (101 lines) | stat: -rw-r--r-- 4,177 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
! -*- f90 -*-
!
! Contains wrappers for the following LAPACK routines:
!
!  Simple Driver Routines for Linear Equations:
!   gesv (general)
!   gbsv (general band)
!   gtsv (general tridiagonla) - Not Implemented
!   posv (symmetric/hermitian positive definite)
!   ppsv (symmetric/hermitian positive definite packed storage) - Not Implemented
!   pbsv (symmetric/hermitian positive definite band) - Not Implemeted
!   ptsv (symmetric/hermitian positive definite tridiagonal) - Not Implemented
!   sysv, hesv (symmetric/hermitian indefinite) - Not Implemented
!   spsv, hpsv (symmetric/hermitian indefinite packed storage) - Not Implemented
!
!
!  Expert Driver Routines for Linear Equations:
!   gesvx (general) - Not Implemented
!   gbsvx (general band) - Not Implemented
!   gtsvx (general tridiagonla) - Not Implemented
!   posvx (symmetric/hermitian positive definite) - Not Implemented
!   ppsvx (symmetric/hermitian positive definite packed storage) - Not Implemented
!   pbsvx (symmetric/hermitian positive definite band) - Not Implemeted
!   ptsvx (symmetric/hermitian positive definite tridiagonal) - Not Implemented
!   sysvx, hesvx (symmetric/hermitian indefinite) - Not Implemented
!   spsvx, hpsvx (symmetric/hermitian indefinite packed storage) - Not Implemented
!

!
! Simple Driver Routines for Linear Equations
! ===========================================

   subroutine <prefix>gesv(n,nrhs,a,piv,b,info)

   ! lu,piv,x,info = gesv(a,b,overwrite_a=0,overwrite_b=0)
   ! Solve A * X = B.
   ! A = P * L * U
   ! U is upper diagonal triangular, L is unit lower triangular,
   ! piv pivots columns.

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

     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(in,out,copy,out=lu) a

   end subroutine <prefix>gesv

   subroutine <prefix>gbsv(n,kl,ku,nrhs,ab,piv,b,info)
   ! 
   ! lub,piv,x,info = gbsv(kl,ku,ab,b,overwrite_ab=0,overwrite_b=0)
   ! Solve A * X = B
   ! A = P * L * U
   ! A is a band matrix of order n with kl subdiagonals and ku superdiagonals
   ! starting at kl-th row.
   ! X, B are n-by-nrhs matrices
   !
     callstatement {int i=2*kl+ku+1;(*f2py_func)(&n,&kl,&ku,&nrhs,ab,&i,piv,b,&n,&info);for(i=0;i\<n;--piv[i++]);}
     callprotoargument int*,int*,int*,int*,<ctype>*,int*,int*,<ctype>*,int*,int*
     integer depend(ab),intent(hide):: n = shape(ab,1)
     integer intent(in) :: kl
     integer intent(in) :: ku
     integer depend(b),intent(hide) :: nrhs = shape(b,1)
     <ftype> dimension(2*kl+ku+1,n),depend(kl,ku), check(2*kl+ku+1==shape(ab,0)) :: ab
     integer dimension(n),depend(n),intent(out) :: piv
     <ftype> dimension(n,nrhs),depend(n),check(shape(ab,1)==shape(b,0)) :: b
     integer intent(out) :: info
     intent(in,out,copy,out=x) b
     intent(in,out,copy,out=lub) ab
   end subroutine <prefix>gbsv

   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