[SciPy-dev] New clapack and flapack

Pearu Peterson pearu at cens.ioc.ee
Sat Jan 19 13:07:39 EST 2002


Hi Eric,

I have finished new interfaces to both ATLAS C LAPACK and Fortran
LAPACK. The generic pyf files and the reduced interface_gen.py are
included in attachements (its is much faster now). Here is how to build
the extension modules from scratch (you'll need the latest f2py, though):

1.
$ python interface_gen.py   # this creates clapack.pyf and flapack.pyf

2.
$ f2py -c clapack.pyf -L/usr/local/lib/atlas -llapack -lf77blas -lcblas  -latlas -lg2c

$ f2py -c flapack.pyf  -L/usr/local/lib/atlas -llapack -lf77blas -lcblas -latlas -lg2c

This builds clapack.so and flapack.so.

3.
In python:
>>> import clapack,flapack
>>> print clapack.__doc__ 
This module 'clapack' is auto-generated with f2py (version:2.11.174-1131).
Functions:
  LU,piv,X,info = sgesv(A,B,rowmajor=1,copy_A=1,copy_B=1)
  LU,piv,X,info = dgesv(A,B,rowmajor=1,copy_A=1,copy_B=1)
  LU,piv,X,info = cgesv(A,B,rowmajor=1,copy_A=1,copy_B=1)
  LU,piv,X,info = zgesv(A,B,rowmajor=1,copy_A=1,copy_B=1)
  LU,piv,info = sgetrf(A,rowmajor=1,copy_A=1)
  LU,piv,info = dgetrf(A,rowmajor=1,copy_A=1)
  LU,piv,info = cgetrf(A,rowmajor=1,copy_A=1)
  LU,piv,info = zgetrf(A,rowmajor=1,copy_A=1)
....
>>> print flapack.__doc__ 
This module 'flapack' is auto-generated with f2py (version:2.11.174-1131).
Functions:
  LU,piv,X,info = sgesv(A,B,copy_A=1,copy_B=1)
  LU,piv,X,info = dgesv(A,B,copy_A=1,copy_B=1)
  LU,piv,X,info = cgesv(A,B,copy_A=1,copy_B=1)
  LU,piv,X,info = zgesv(A,B,copy_A=1,copy_B=1)
  LU,piv,info = sgetrf(A,copy_A=1)
  LU,piv,info = dgetrf(A,copy_A=1)
  LU,piv,info = cgetrf(A,copy_A=1)
  LU,piv,info = zgetrf(A,copy_A=1)
....

Few remarks:
1) The sources to these extension modules are quite big (more than 5000
lines in both). Should we factor these modules to
   cslapack.so cclapack.so cdlapack.so czlapack.so
   fslapack.so fclapack.so fdlapack.so fzlapack.so
(any better name convention?). What do you think?
2) How shall I commit all this to the scipy CVS repository? In a
subdirectory of linalg? Eventually it should be in linalg directory, I
think.

Regards,
	Pearu
-------------- next part --------------
#!/usr/bin/env python

import string,re,os

def all_subroutines(interface_in):
    # remove comments
    comment_block_exp = re.compile(r'/\*(?:\s|.)*?\*/')
    subroutine_exp = re.compile(r'subroutine (?:\s|.)*?end subroutine.*')
    function_exp = re.compile(r'function (?:\s|.)*?end function.*')
    
    interface = comment_block_exp.sub('',interface_in)    
    subroutine_list = subroutine_exp.findall(interface)
    function_list = function_exp.findall(interface)
    #function_list = []
    subroutine_list = subroutine_list + function_list 
    subroutine_list = map(lambda x: string.strip(x),subroutine_list)
    return subroutine_list

def real_convert(val_string):
    return val_string

def complex_convert(val_string):
    return '(' + val_string + ',0.)'

def convert_types(interface_in,converter):
    regexp = re.compile(r'<type_convert=(.*?)>')
    interface = interface_in[:]
    while 1:
        sub = regexp.search(interface)
        if sub is None: break
        converted = converter(sub.group(1))
        interface = string.replace(interface,sub.group(),converted)
    return interface

def generic_expand(generic_interface):
    generic_types ={'s' :('real',            'real', real_convert,
                          'real'), 
                    'd' :('double precision','double precision',real_convert,
                          'double precision'),
                    'c' :('complex',         'complex',complex_convert,
                          'real'), 
                    'z' :('double complex',  'double complex',complex_convert,
                          'double precision'),                   
                    'sc':('complex',         'real',real_convert,
                          'real'),
                    'dz':('double complex',  'double precision',real_convert,
                          'double precision'),
                    'cs':('real',            'complex',complex_convert,
                          'real'),
                    'zd':('double precision','double complex', complex_convert,
                          'double precision')}
    #2. get all subroutines    
    subs = all_subroutines(generic_interface)
    print len(subs)
    #loop through the subs
    type_exp = re.compile(r'<tchar=(.*?)>')    
    interface = ''
    for sub in subs:
        #3. Find the typecodes to use:
        m = type_exp.search(sub)
        if m is None:
            interface = interface + '\n\n' + sub
            continue
        type_chars = m.group(1)
        # get rid of spaces
        type_chars = string.replace(type_chars,' ','')
        # get a list of the characters (or character pairs)
        type_chars = string.split(type_chars,',')
        #print type_chars
        # Now get rid of the special tag that contained the types
        sub = re.sub(type_exp,'<tchar>',sub)
        sub_generic = string.strip(sub)
        for char in type_chars:            
            type_in,type_out,converter, rtype_in = generic_types[char]
            sub = convert_types(sub_generic,converter)
            function_def = string.replace(sub,'<tchar>',char)    
            function_def = string.replace(function_def,'<type_in>',type_in)
            #function_def = string.replace(function_def,'<rtype_in>',rtype_in) 
            #function_def = string.replace(function_def,'<type_out>',type_out)
            interface = interface + '\n\n' + function_def

    return interface

#def interface_to_module(interface_in,module_name,include_list,sdir='.'):
def interface_to_module(interface_in,module_name):
    pre_prefix = "!%f90 -*- f90 -*-\n"
    # includes = ''
#     for file in include_list:
#         f = open(os.path.join(sdir,file))
#         includes += f.read(-1)
#         f.close() 
    # heading and tail of the module definition.
    file_prefix = "\npython module " + module_name +" ! in\n" \
             "    interface  \n"
    file_suffix = "\n    end interface\n" \
             "end module %s" % module_name
    #return  pre_prefix + includes + file_prefix + interface_in + file_suffix
    return  pre_prefix + file_prefix + interface_in + file_suffix


def generate_clapack(sdir,output_path='.'):
    print "generating clapack interface"
    f = open(os.path.join(sdir,'generic_clapack.pyf'))
    module_name = 'clapack'
    generic_interface = f.read(-1)
    f.close()
    generic_interface = generic_expand(generic_interface)
    module_def = interface_to_module(generic_interface,module_name)
    f = open(os.path.join(output_path,module_name+'.pyf'),'w')
    f.write(module_def)
    f.close()
def generate_flapack(sdir,output_path='.'):
    print "generating flapack interface"
    f = open(os.path.join(sdir,'generic_flapack.pyf'))
    module_name = 'flapack'
    generic_interface = f.read(-1)
    f.close()
    generic_interface = generic_expand(generic_interface)
    module_def = interface_to_module(generic_interface,module_name)
    f = open(os.path.join(output_path,module_name+'.pyf'),'w')
    f.write(module_def)
    f.close()

def process_all():
    # process the standard files.
    generate_clapack('.')
    generate_flapack('.')

if __name__ == "__main__":
    process_all()
-------------- next part --------------
!%f90 -*- f90 -*-
python module generic_flapack
interface

   subroutine <tchar=s,d,c,z>gesv(n,nrhs,A,piv,B,info)

   ! LU,piv,X,info = gesv(A,B,copy_A=1,copy_B=1)
   ! Solve A * X = B.
   ! A = P * L * U
   ! U is upper diagonal triangular, L is unit lower triangular,
   ! piv pivots columns.

     callstatement (*f2py_func)(&n,&nrhs,A,&n,piv,B,&n,&info)

     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(in,out,copy,out=LU) A

   end subroutine <tchar=s,d,c,z>gesv

   subroutine <tchar=s,d,c,z>getrf(m,n,A,piv,info)

   ! LU,piv,info = getrf(A,copy_A=1)
   ! Compute an LU factorization of a  general  M-by-N  matrix  A.
   ! A = P * L * U

     callstatement (*f2py_func)(&m,&n,A,&m,piv,&info)

     integer depend(A),intent(hide):: m = shape(A,0)
     integer depend(A),intent(hide):: n = shape(A,1)
     <type_in> dimension(m,n),intent(in,out,copy,out=LU) :: A
     integer dimension((m<n?m:n)),depend(m,n),intent(out) :: piv
     integer intent(out):: info

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

   subroutine <tchar=s,d,c,z>getrs(n,nrhs,LU,piv,B,info,trans)

   ! X,info = getrs(LU,piv,B,trans=0,copy_B=1)
   ! 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

     callstatement (*f2py_func)((trans?(trans==2?'C':'T'):'N'),&n,&nrhs,LU,&n,piv,B,&n,&info)

     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(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 subroutine <tchar=s,d,c,z>getrs


   subroutine <tchar=s,d,c,z>getri(n,LU,piv,work,lwork,info)

   ! invA,info = getri(LU,piv,copy_LU=1)
   ! Find A inverse A^-1.
   ! A = P * L * U

     callstatement (*f2py_func)(&n,LU,&n,piv,work,&lwork,&info)

     integer depend(LU),intent(hide):: n = shape(LU,0)
     <type_in> dimension(n,n),intent(in,out,copy,out=invA) :: 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
     <type_in> dimension(lwork),intent(hide),depend(lwork) :: work

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

   subroutine <tchar=s,d,c,z>posv(n,nrhs,A,B,info,lower)

   ! C,X,info = posv(A,B,lower=0,copy_A=1,copy_B=1)
   ! 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)

     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(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 subroutine <tchar=s,d,c,z>posv
        
   subroutine <tchar=s,d,c,z>potrf(n,A,info,lower)
   
     ! C,info = potrf(A,lower=0=1,copy_A=1)
     ! 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.

     callstatement (*f2py_func)((lower?'L':'U'),&n,A,&n,&info)
     
     integer optional,intent(in),check(lower==0||lower==1) :: lower = 0

     integer depend(A),intent(hide):: n = shape(A,0)
     <type_in> dimension(n,n),intent(in,out,copy,out=C) :: A
     check(shape(A,0)==shape(A,1)) :: A
     integer intent(out) :: info
     
   end subroutine <tchar=s,d,c,z>potrf
   
   
   subroutine <tchar=s,d,c,z>potrs(n,nrhs,C,B,info,lower)

   ! X,info = potrs(C,B,lower=0=1,copy_B=1)
   ! 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)

     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(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 subroutine <tchar=s,d,c,z>potrs

   subroutine <tchar=s,d,c,z>potri(n,C,info,lower)
   
     ! invA,info = potri(C,lower=0,copy_C=1)
     ! 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)

     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=invA) :: C
     check(shape(C,0)==shape(C,1)) :: C
     integer intent(out) :: info
     
   end subroutine <tchar=s,d,c,z>potri


   subroutine <tchar=s,d,c,z>lauum(n,C,info,lower)
   
     ! A,info = lauum(C,lower=0,copy_C=1)
     ! 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.

     callstatement (*f2py_func)((lower?'L':'U'),&n,C,&n,&info)
     
     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(in,out,copy,out=A) :: C
     check(shape(C,0)==shape(C,1)) :: C
     integer intent(out) :: info
     
   end subroutine <tchar=s,d,c,z>lauum

   subroutine <tchar=s,d,c,z>trtri(n,C,info,lower,unitdiag)
   
     ! invC,info = trtri(C,lower=0,unitdiag=1,copy_C=1)
     ! 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)
     
     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(in,out,copy,out=invC) :: C
     check(shape(C,0)==shape(C,1)) :: C
     integer intent(out) :: info
     
   end subroutine <tchar=s,d,c,z>trtri

end interface

end python module generic_flapack

! This file was auto-generated with f2py (version:2.10.173).
! See http://cens.ioc.ee/projects/f2py2e/
-------------- next part --------------
!%f90 -*- f90 -*-
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,copy_A=1,copy_B=1)
   ! 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 info = (*f2py_func)(102-rowmajor,n,nrhs,A,n,piv,B,n)

     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,copy_A=1)
   ! Compute an LU factorization of a  general  M-by-N  matrix  A.
   ! A * P = L * U

     fortranname  clapack_<tchar=s,d,c,z>getrf
     integer intent(c,hide) ::  <tchar=s,d,c,z>getrf
     callstatement info = (*f2py_func)(102-rowmajor,m,n,A,(rowmajor?n:m),piv)

     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,copy_B=1)
   ! 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 info = (*f2py_func)(102-rowmajor,110+trans,n,nrhs,LU,n,piv,B,n)

     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)

   ! invA,info = getri(LU,piv,rowmajor=1,copy_LU=1)
   ! 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 info = (*f2py_func)(102-rowmajor,n,LU,n,piv)

     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=invA) :: 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,copy_A=1,copy_B=1)
   ! 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 info = (*f2py_func)(102-rowmajor,121+lower,n,nrhs,A,n,B,n)

     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,c,z>potrf(n,A,info,lower,rowmajor)
   
     ! C,info = potrf(A,lower=0,rowmajor=1,copy_A=1)
     ! 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.

     fortranname  clapack_<tchar=s,d,c,z>potrf
     integer intent(c,hide) ::  <tchar=s,d,c,z>potrf
     callstatement info = (*f2py_func)(102-rowmajor,121+lower,n,A,n)
     
     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)
     <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,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,copy_B=1)
   ! 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 info = (*f2py_func)(102-rowmajor,121+lower,n,nrhs,C,n,B,n)

     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)
   
     ! invA,info = potri(C,lower=0,rowmajor=1,copy_C=1)
     ! 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 info = (*f2py_func)(102-rowmajor,121+lower,n,C,n)

     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=invA) :: 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,copy_C=1)
     ! 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 info = (*f2py_func)(102-rowmajor,121+lower,n,C,n)
     
     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)
   
     ! invC,info = trtri(C,lower=0,unitdiag=0,rowmajor=1,copy_C=1)
     ! 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 info = (*f2py_func)(102-rowmajor,121+lower,131+unitdiag,n,C,n)
     
     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=invC) :: 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/


More information about the SciPy-Dev mailing list