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
|
! -*- f90 -*-
!
! Iterative Package for SciPy
! Hongze Liu, Travis E. Oliphant,
! Brigham Young University
! 2004
!
python module _iterative ! in
interface ! in :_iterative
subroutine <_c>bicgrevcom(n,b,x,work,ldw,iter,resid,info,ndx1,ndx2,sclr1,sclr2,ijob) ! in :iterative:BiCG.f
integer, intent(hide), depend(b) :: n=len(b)
<_t> dimension(n) :: b
<_t> dimension(n), intent(in,out) :: x
<_t> intent(inout), dimension(ldw*6) :: work
integer, intent(hide), depend(n) :: ldw=MAX(1,n)
integer, intent(in,out) :: iter
<rt=real, double precision, real, double precision>, intent(in,out) :: resid
integer, intent(in, out) :: info
integer, intent(in, out) :: ndx1
integer, intent(in, out) :: ndx2
<_t>, intent(out) :: sclr1
<_t>, intent(out) :: sclr2
integer, intent(in, out) :: ijob
end subroutine <_c>bicgrevcom
subroutine <_c>bicgstabrevcom(n,b,x,work,ldw,iter,resid,info,ndx1,ndx2,sclr1,sclr2,ijob) ! in :iterative:BiCGSTAB.f
integer, intent(hide), depend(b) :: n=len(b)
<_t> dimension(n) :: b
<_t> dimension(n), intent(in,out) :: x
<_t> intent(inout), dimension(ldw*7) :: work
integer, intent(hide), depend(n) :: ldw=MAX(1,n)
integer, intent(in,out) :: iter
<rt=real, double precision, real, double precision>, intent(in,out) :: resid
integer, intent(in, out) :: info
integer, intent(in, out) :: ndx1
integer, intent(in, out) :: ndx2
<_t>, intent(out) :: sclr1
<_t>, intent(out) :: sclr2
integer, intent(in, out) :: ijob
end subroutine <_c>bicgstabrevcom
subroutine <_c>cgrevcom(n,b,x,work,ldw,iter,resid,info,ndx1,ndx2,sclr1,sclr2,ijob) ! in :iterative:CG.f
integer, intent(hide), depend(b) :: n=len(b)
<_t> dimension(n) :: b
<_t> dimension(n), intent(in,out) :: x
<_t> intent(inout), dimension(ldw*4) :: work
integer, intent(hide), depend(n) :: ldw=MAX(1,n)
integer, intent(in,out) :: iter
<rt=real, double precision, real, double precision>, intent(in,out) :: resid
integer, intent(in, out) :: info
integer, intent(in, out) :: ndx1
integer, intent(in, out) :: ndx2
<_t>, intent(out) :: sclr1
<_t>, intent(out) :: sclr2
integer, intent(in, out) :: ijob
end subroutine <_c>cgrevcom
subroutine <_c>cgsrevcom(n,b,x,work,ldw,iter,resid,info,ndx1,ndx2,sclr1,sclr2,ijob) ! in :iterative:CGS.f
integer, intent(hide), depend(b) :: n=len(b)
<_t> dimension(n) :: b
<_t> dimension(n), intent(in,out) :: x
<_t> intent(inout), dimension(ldw*7) :: work
integer, intent(hide), depend(n) :: ldw=MAX(1,n)
integer, intent(in,out) :: iter
<rt=real, double precision, real, double precision>, intent(in,out) :: resid
integer, intent(in, out) :: info
integer, intent(in, out) :: ndx1
integer, intent(in, out) :: ndx2
<_t>, intent(out) :: sclr1
<_t>, intent(out) :: sclr2
integer, intent(in, out) :: ijob
end subroutine <_c>cgsrevcom
subroutine <_c>qmrrevcom(n,b,x,work,ldw,iter,resid,info,ndx1,ndx2,sclr1,sclr2,ijob) ! in :iterative:QMR.f
integer, intent(hide), depend(b) :: n=len(b)
<_t> dimension(n) :: b
<_t> dimension(n), intent(in,out) :: x
<_t> intent(inout), dimension(ldw*11) :: work
integer, intent(hide), depend(n) :: ldw=MAX(1,n)
integer, intent(in,out) :: iter
<rt=real, double precision, real, double precision>, intent(in,out) :: resid
integer, intent(in, out) :: info
integer, intent(in, out) :: ndx1
integer, intent(in, out) :: ndx2
<_t>, intent(out) :: sclr1
<_t>, intent(out) :: sclr2
integer, intent(in, out) :: ijob
end subroutine <_c>qmrrevcom
subroutine <_c>gmresrevcom(n,b,x,restrt,work,ldw,work2,ldw2,iter,resid,info,ndx1,ndx2,sclr1,sclr2,ijob) ! in :iterative:GMRESREVCOM.f
integer, intent(hide), depend(b) :: n=len(b)
<_t> dimension(n) :: b
<_t> dimension(n), intent(in,out) :: x
integer, intent(in), depend(n), check((0<restrt) && (restrt<=n)) :: restrt
<_t> intent(inout), dimension(ldw*(6+restrt)) :: work
integer intent(hide) :: ldw=MAX(1,n)
<_t> intent(inout), depend(restrt,ldw2), dimension(ldw2*(2*restrt+2)) :: work2
integer intent(hide), depend(restrt) :: ldw2=MAX(2,restrt+1)
integer intent(in, out) :: iter
<rt=real, double precision, real, double precision>, intent(in,out) :: resid
integer intent(in, out) :: info
integer intent(in, out) :: ndx1
integer intent(in, out) :: ndx2
<_t> intent(out) :: sclr1
<_t> intent(out) :: sclr2
integer intent(in, out) :: ijob
end subroutine <_c>gmresrevcom
subroutine <_c>stoptest2(n,r,b,bnrm2,resid,tol,info) ! in STOPTEST2.f
integer, intent(hide), depend(b) :: n=len(b)
<_t>, dimension(n), intent(in) :: r
<_t>, dimension(n), intent(in) :: b
<rt=real, double precision, real, double precision>, intent(in, out) :: bnrm2
<rt>, intent(out) :: resid
<rt>, intent(in) :: tol
integer, intent(in, out) :: info
end subroutine <_c>stoptest2
end interface
end python module _iterative
! This file was auto-generated with f2py (version:2.39.235_1703).
! See http://cens.ioc.ee/projects/f2py2e/
|