TABLE OF CONTENTS
- m_abi_linalg/abi_d2ztrsm
- m_abi_linalg/abi_d2ztrsm_3d
- m_abi_linalg/abi_dtrsm
- m_abi_linalg/abi_xtrsm
- m_abi_linalg/abi_ztrsm
m_abi_linalg/abi_d2ztrsm [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_d2ztrsm
FUNCTION
INPUTS
PARENTS
SOURCE
171 subroutine abi_d2ztrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb,& 172 & x_cplx) 173 174 !This section has been created automatically by the script Abilint (TD). 175 !Do not modify the following lines by hand. 176 #undef ABI_FUNC 177 #define ABI_FUNC 'abi_d2ztrsm' 178 !End of the abilint section 179 180 implicit none 181 !Arguments------------------------------------- 182 character(len=1), intent(in) :: side,uplo,transa,diag 183 integer, intent(in) :: m,n,lda,ldb 184 complex(dpc), intent(in) :: alpha 185 real(dp),target, intent(in) :: a(lda,*) ! FIXME should be lda * x_cplx 186 real(dp),target, intent(inout) :: b(ldb,*) 187 !Only for lobpcgwf 188 integer, intent(in), optional :: x_cplx 189 190 !Local variables------------------------------- 191 integer :: cplx_ 192 #ifdef HAVE_LINALG_PLASMA 193 integer :: info 194 #endif 195 196 #ifdef DEV_LINALG_TIMING 197 real(dp) :: tsec(2) 198 call timab(TIMAB_XTRSM,1,tsec) 199 #endif 200 201 cplx_=1 ; if(PRESENT(x_cplx)) cplx_ = x_cplx 202 203 if (XPLASMA_ISON) then 204 #ifdef HAVE_LINALG_PLASMA 205 if(cplx_ == 2) then 206 info = PLASMA_ztrsm_c(side_plasma(side),uplo_plasma(uplo),trans_plasma(TRANSA),diag_plasma(diag),& 207 & m,n,alpha,c_loc(a),lda,c_loc(b),ldb) 208 else 209 info = PLASMA_dtrsm_c(side_plasma(side),uplo_plasma(uplo),trans_plasma(TRANSA),diag_plasma(diag),& 210 & m,n,real(alpha,dp),c_loc(a),lda,c_loc(b),ldb) 211 end if 212 #endif 213 else 214 if(cplx_ == 2) then 215 call ztrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb) 216 else 217 call dtrsm(side,uplo,transa,diag,m,n,real(alpha,dp),a,lda,b,ldb) 218 end if 219 end if 220 221 #ifdef DEV_LINALG_TIMING 222 call timab(TIMAB_XTRSM,2,tsec) 223 #endif 224 225 end subroutine abi_d2ztrsm
m_abi_linalg/abi_d2ztrsm_3d [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_d2ztrsm_3d
FUNCTION
INPUTS
PARENTS
SOURCE
243 subroutine abi_d2ztrsm_3d(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb) 244 245 !This section has been created automatically by the script Abilint (TD). 246 !Do not modify the following lines by hand. 247 #undef ABI_FUNC 248 #define ABI_FUNC 'abi_d2ztrsm_3d' 249 !End of the abilint section 250 251 implicit none 252 253 !Arguments------------------------------------- 254 character(len=1), intent(in) :: side,uplo,transa,diag 255 integer, intent(in) :: m,n,lda,ldb 256 complex(dpc), intent(in) :: alpha 257 real(dp), target,intent(in) :: a(2,lda,*) 258 real(dp), target,intent(inout) :: b(2,ldb,*) 259 260 !Local variables------------------------------- 261 #ifdef HAVE_LINALG_PLASMA 262 integer :: info 263 #endif 264 265 #ifdef DEV_LINALG_TIMING 266 real(dp) :: tsec(2) 267 call timab(TIMAB_XTRSM,1,tsec) 268 #endif 269 270 if (XPLASMA_ISON) then 271 #ifdef HAVE_LINALG_PLASMA 272 info = PLASMA_ztrsm_c(side_plasma(side),uplo_plasma(uplo),trans_plasma(TRANSA),diag_plasma(diag),& 273 & m,n,alpha,c_loc(a),lda,c_loc(b),ldb) 274 #endif 275 else 276 call ztrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb) 277 end if 278 279 #ifdef DEV_LINALG_TIMING 280 call timab(TIMAB_XTRSM,2,tsec) 281 #endif 282 283 end subroutine abi_d2ztrsm_3d
m_abi_linalg/abi_dtrsm [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_dtrsm
FUNCTION
INPUTS
PARENTS
SOURCE
99 subroutine abi_dtrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb,& 100 & x_cplx) 101 102 !This section has been created automatically by the script Abilint (TD). 103 !Do not modify the following lines by hand. 104 #undef ABI_FUNC 105 #define ABI_FUNC 'abi_dtrsm' 106 !End of the abilint section 107 108 implicit none 109 110 !Arguments------------------------------------- 111 character(len=1), intent(in) :: side,uplo,transa,diag 112 integer, intent(in) :: m,n,lda,ldb 113 real(dp), intent(in) :: alpha 114 real(dp),target, intent(in) :: a(lda,*) ! FIXME should be lda * x_cplx 115 real(dp),target, intent(inout) :: b(ldb,*) 116 !Only for lobpcgwf 117 integer, intent(in), optional :: x_cplx 118 119 !Local variables------------------------------- 120 integer :: cplx_ 121 #ifdef HAVE_LINALG_PLASMA 122 integer :: info 123 #endif 124 125 #ifdef DEV_LINALG_TIMING 126 real(dp) :: tsec(2) 127 call timab(TIMAB_XTRSM,1,tsec) 128 #endif 129 130 cplx_=1 ; if(PRESENT(x_cplx)) cplx_ = x_cplx 131 132 if (XPLASMA_ISON) then 133 #ifdef HAVE_LINALG_PLASMA 134 if(cplx_ == 2) then 135 info = PLASMA_ztrsm_c(side_plasma(side),uplo_plasma(uplo),trans_plasma(TRANSA),diag_plasma(diag),& 136 & m,n,cmplx(alpha,0.d0,dpc),c_loc(a),lda,c_loc(b),ldb) 137 else 138 info = PLASMA_dtrsm_c(side_plasma(side),uplo_plasma(uplo),trans_plasma(TRANSA),diag_plasma(diag),& 139 & m,n,alpha,c_loc(a),lda,c_loc(b),ldb) 140 end if 141 #endif 142 else 143 if(cplx_ == 2) then 144 call ztrsm(side,uplo,transa,diag,m,n,cmplx(alpha,0.d0,dpc),a,lda,b,ldb) 145 else 146 call dtrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb) 147 end if 148 end if 149 150 #ifdef DEV_LINALG_TIMING 151 call timab(TIMAB_XTRSM,2,tsec) 152 #endif 153 154 end subroutine abi_dtrsm
m_abi_linalg/abi_xtrsm [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_xtrsm
FUNCTION
abi_xtrsm is the generic function that solve : * op( A )*X = alpha*B, or X*op( A ) = alpha*B, * * where alpha is a scalar, X and B are m by n matrices, A is a unit, or * non-unit, upper or lower triangular matrix and op( A ) is one of * * op( A ) = A or op( A ) = A**T. * * The matrix X is overwritten on B.
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_ztrsm [ Functions ]
[ Top ] [ m_abi_linalg ] [ Functions ]
NAME
abi_ztrsm
FUNCTION
INPUTS
PARENTS
SOURCE
39 subroutine abi_ztrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb) 40 41 !This section has been created automatically by the script Abilint (TD). 42 !Do not modify the following lines by hand. 43 #undef ABI_FUNC 44 #define ABI_FUNC 'abi_ztrsm' 45 !End of the abilint section 46 47 implicit none 48 49 !Arguments------------------------------------- 50 character(len=1), intent(in) :: side 51 character(len=1), intent(in) :: uplo 52 character(len=1), intent(in) :: transa 53 character(len=1), intent(in) :: diag 54 integer, intent(in) :: m,n,ldb,lda 55 complex(dpc), intent(in) :: alpha 56 complex(dpc),target,intent(in) :: a(lda,*) 57 complex(dpc),target,intent(inout) :: b(ldb,*) 58 59 !Local variables------------------------------- 60 #ifdef HAVE_LINALG_PLASMA 61 integer :: info 62 #endif 63 64 #ifdef DEV_LINALG_TIMING 65 real(dp) :: tsec(2) 66 call timab(TIMAB_XTRSM,1,tsec) 67 #endif 68 69 if (XPLASMA_ISON) then 70 #ifdef HAVE_LINALG_PLASMA 71 info = PLASMA_ztrsm_c(side_plasma(side),uplo_plasma(uplo),trans_plasma(transa),diag_plasma(diag),& 72 & m,n,alpha,c_loc(a),lda,c_loc(b),ldb) 73 #endif 74 else 75 call ztrsm(side,uplo,transa,diag,m,n,alpha,a,lda,b,ldb) 76 end if 77 78 #ifdef DEV_LINALG_TIMING 79 call timab(TIMAB_XTRSM,2,tsec) 80 #endif 81 82 end subroutine abi_ztrsm