TABLE OF CONTENTS


m_abi_linalg/abi_d2zpotrf [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_d2zpotrf

FUNCTION

INPUTS

PARENTS

SOURCE

115 subroutine abi_d2zpotrf(uplo,n,a,lda,info,x_cplx)
116 
117 !This section has been created automatically by the script Abilint (TD).
118 !Do not modify the following lines by hand.
119 #undef ABI_FUNC
120 #define ABI_FUNC 'abi_d2zpotrf'
121 !End of the abilint section
122 
123  implicit none
124 
125 !Arguments ------------------------------------
126  character(len=1), intent(in) :: uplo
127  integer, intent(in) :: n,lda
128  integer, intent(out) :: info
129  integer, intent(in), optional :: x_cplx
130  real(dp),target, intent(inout) :: a(lda,*)  ! FIXME should be x_cplx * lda
131 
132  !Local Variables -----------------------------
133  integer  :: cplx_
134 
135 ! *********************************************************************
136 
137  cplx_=1 ; if(PRESENT(x_cplx)) cplx_ = x_cplx
138 
139 #ifdef HAVE_LINALG_PLASMA
140  if (XPLASMA_ISON) then
141    if(cplx_ == 2) then
142       info = PLASMA_zpotrf_c(uplo_plasma(uplo),n,c_loc(a),lda)
143    else
144       info = PLASMA_dpotrf_c(uplo_plasma(uplo),n,c_loc(a),lda)
145    end if
146    return
147  end if
148 #endif
149 
150  if(cplx_ == 2) then
151     call zpotrf(uplo,n,a,lda,info)
152  else
153     call dpotrf(uplo,n,a,lda,info)
154  end if
155 
156 end subroutine abi_d2zpotrf

m_abi_linalg/abi_dpotrf [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_dpotrf

FUNCTION

INPUTS

PARENTS

SOURCE

37 subroutine abi_dpotrf(uplo,n,a,lda,info)
38 
39 !This section has been created automatically by the script Abilint (TD).
40 !Do not modify the following lines by hand.
41 #undef ABI_FUNC
42 #define ABI_FUNC 'abi_dpotrf'
43 !End of the abilint section
44 
45  implicit none
46  !Arguments ------------------------------------
47  character(len=1), intent(in) :: uplo
48  integer, intent(in) :: n,lda
49  integer, intent(out) :: info
50  real(dp), intent(inout) :: a(*)
51 
52 ! *********************************************************************
53 
54 #ifdef HAVE_LINALG_PLASMA
55  if (XPLASMA_ISON) then
56    ! write(std_out,*) "  abi_dpotrf => PLASMA dpotrf will be called "
57    call PLASMA_dpotrf(uplo_plasma(uplo),n,a,lda,info)
58    return
59  end if
60 #endif
61 
62  call dpotrf(uplo,n,a,lda,info)
63 
64 end subroutine abi_dpotrf

m_abi_linalg/abi_xpotrf [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

  abi_xpotrf

FUNCTION

  abi_xpotrf is the generic function for computing the
  Cholesky factorization of a real symmetric (or hermitian)
    positive definite matrix A.
    The factorization has the form
      A = U**T * U,  if UPLO = 'U', or
      A = L  * L**T,  if UPLO = 'L',
    where U is an upper triangular matrix and L is lower triangular.

COPYRIGHT

  Copyright (C) 2001-2018 ABINIT group (LNguyen,FDahm (CS))
  This file is distributed under the terms of the
  GNU General Public License, see ~ABINIT/Infos/copyright
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE


m_abi_linalg/abi_zpotrf [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zpotrf

FUNCTION

INPUTS

PARENTS

SOURCE

171 subroutine abi_zpotrf(uplo,n,a,lda,info)
172 
173 !This section has been created automatically by the script Abilint (TD).
174 !Do not modify the following lines by hand.
175 #undef ABI_FUNC
176 #define ABI_FUNC 'abi_zpotrf'
177 !End of the abilint section
178 
179  implicit none
180  !Arguments ------------------------------------
181  character(len=1), intent(in) :: uplo
182  integer, intent(in) :: lda,n
183  integer, intent(out) :: info
184  complex(dpc), intent(inout) :: a(*)
185 
186 ! *********************************************************************
187 
188 #ifdef HAVE_LINALG_PLASMA
189  if (XPLASMA_ISON) then
190    ! write(*,*) "  abi_zpotrf => PLASMA zpotrf will be called "
191    call PLASMA_zpotrf(uplo_plasma(uplo),n,a,lda,info)
192    return
193  end if
194 #endif
195 
196  call zpotrf(uplo,n,a,lda,info)
197 
198 end subroutine abi_zpotrf

m_abi_linalg/abi_zpotrf_2d [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zpotrf_2d

FUNCTION

INPUTS

PARENTS

SOURCE

 79 subroutine abi_zpotrf_2d(uplo,n,a,lda,info)
 80 
 81 !This section has been created automatically by the script Abilint (TD).
 82 !Do not modify the following lines by hand.
 83 #undef ABI_FUNC
 84 #define ABI_FUNC 'abi_zpotrf_2d'
 85 
 86 !End of the abilint section
 87 
 88  implicit none
 89 
 90  !Arguments ------------------------------------
 91  character(len=1), intent(in) :: uplo
 92  integer, intent(in) :: lda,n
 93  integer, intent(out) :: info
 94  complex(dpc), intent(inout) :: a(lda,*)
 95 
 96 ! *********************************************************************
 97 
 98  call abi_zpotrf(uplo,n,a(1,1),lda,info)
 99 
100 end subroutine abi_zpotrf_2d