TABLE OF CONTENTS


m_abi_linalg/abi_d2zgemm [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_d2zgemm

FUNCTION

INPUTS

SOURCE

157 subroutine abi_d2zgemm(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC,&
158                        x_cplx)
159 
160 !Arguments ------------------------------------
161  character(len=1), intent(in) :: transa
162  character(len=1), intent(in) :: transb
163  integer, intent(in) :: lda
164  integer, intent(in) :: ldb
165  integer, intent(in) :: ldc
166  integer, intent(in) :: m
167  integer, intent(in) :: n
168  integer, intent(in) :: k
169  complex(dpc),intent(in) :: alpha
170  complex(dpc),intent(in) :: beta
171  real(dp),target, intent(in)    :: a(lda,*)    ! FIXME should be lda * x_cplx
172  real(dp),target, intent(in)    :: b(ldb,*)
173  real(dp),target, intent(inout) :: c(ldc,*)
174 !Only for lobpcgwf
175  integer, intent(in), optional :: x_cplx
176 
177 !Local variables-------------------------------
178  integer  :: cplx_,info
179 #ifdef DEV_LINALG_TIMING
180  real(dp) :: tsec(2)
181  call timab(TIMAB_XGEMM,1,tsec)
182 #endif
183 
184  cplx_=1 ; if(PRESENT(x_cplx)) cplx_ = x_cplx
185 
186  if (ABI_LINALG_PLASMA_ISON) then
187    info = -1
188 #ifdef HAVE_LINALG_PLASMA
189    !write(std_out,*) "Will call plasma_[zd]gemm"
190    if (cplx_ == 2) then
191       info = PLASMA_zgemm_c(trans_plasma(TRANSA),trans_plasma(TRANSB),M,N,K,&
192 &      ALPHA,c_loc(A),LDA,c_loc(B),LDB,BETA,c_loc(C),LDC)
193    else
194       info = PLASMA_dgemm_c(trans_plasma(TRANSA),trans_plasma(TRANSB),M,N,K,&
195 &      real(ALPHA,dp),c_loc(A),LDA,c_loc(B),LDB,real(BETA,dp),c_loc(C),LDC)
196    end if
197 #endif
198    ABI_CHECK(info==0,"PLASMA_[z,d]gemm_c returned info !=0")
199  else
200    if (cplx_ == 2) then
201      if (use_zgemm3m(m,n,k)) then
202        call _ZGEMM3M(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
203      else
204        call zgemm(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
205      end if
206    else
207      call dgemm(TRANSA,TRANSB,M,N,K,real(ALPHA,dp),A,LDA,B,LDB,real(BETA,dp),C,LDC)
208    end if
209  end if
210 
211 #ifdef DEV_LINALG_TIMING
212  call timab(TIMAB_XGEMM,2,tsec)
213 #endif
214 
215 end subroutine abi_d2zgemm

m_abi_linalg/abi_xgemm [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

  abi_xgemm

FUNCTION

  abi_xgemm is the generic function that solves:

    C := alpha*op( A )*op( B ) + beta*C,

  where  op( X ) is one of

    op( X ) = X   or   op( X ) = X**T,

 alpha and beta are scalars, and A, B and C are matrices, with op( A )
 an m by k matrix,  op( B )  a  k by n matrix and  C an m by n matrix.

COPYRIGHT

  Copyright (C) 2001-2024 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_zgemm_2d [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zgemm_2d

FUNCTION

INPUTS

SOURCE

40  subroutine abi_zgemm_2d(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
41 
42  !Arguments------------------------------------
43  character(len=1),intent(in) :: TRANSA
44  character(len=1),intent(in) :: TRANSB
45  integer,intent(in) :: K
46  integer,intent(in) :: LDA
47  integer,intent(in) :: LDB
48  integer,intent(in) :: LDC
49  integer,intent(in) :: M
50  integer,intent(in) :: N
51  complex(dpc),intent(in) :: ALPHA
52  complex(dpc),intent(in) :: BETA
53  complex(dpc),target,intent(in) :: A(lda,*)
54  complex(dpc),target,intent(in) :: B(ldb,*)
55  complex(dpc),target,intent(inout) :: C(ldc,*)
56 
57  integer :: info
58 #ifdef DEV_LINALG_TIMING
59  real(dp) :: tsec(2)
60  call timab(TIMAB_XGEMM,1,tsec)
61 #endif
62 
63  if (ABI_LINALG_PLASMA_ISON) then
64    info = -1
65 #ifdef HAVE_LINALG_PLASMA
66    !write(std_out,*)"Will call PLASMA_zgemm_c"
67    info = PLASMA_zgemm_c(trans_plasma(TRANSA),trans_plasma(TRANSB),M,N,K,ALPHA,&
68 &    c_loc(A),LDA,c_loc(B),LDB,BETA,c_loc(C),LDC)
69 #endif
70    ABI_CHECK(info==0,"PLASMA_zgemm_c returned info !=0")
71  else
72    if (use_zgemm3m(m,n,k)) then
73      call _ZGEMM3M(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
74    else
75      call zgemm(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
76    end if
77  end if
78 
79 #ifdef DEV_LINALG_TIMING
80  call timab(TIMAB_XGEMM,2,tsec)
81 #endif
82 
83 end subroutine abi_zgemm_2d

m_abi_linalg/abi_zgemm_2r [ Functions ]

[ Top ] [ m_abi_linalg ] [ Functions ]

NAME

 abi_zgemm_2r

FUNCTION

INPUTS

SOURCE

 99  subroutine abi_zgemm_2r(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
100 
101  !Arguments------------------------------------
102  character(len=1),intent(in) :: TRANSA
103  character(len=1),intent(in) :: TRANSB
104  integer,intent(in) :: K
105  integer,intent(in) :: LDA
106  integer,intent(in) :: LDB
107  integer,intent(in) :: LDC
108  integer,intent(in) :: M
109  integer,intent(in) :: N
110  complex(dpc),intent(in) :: ALPHA
111  complex(dpc),intent(in) :: BETA
112  real(dp),target,intent(in) :: A(*)
113  real(dp),target,intent(in) :: B(*)
114  real(dp),target,intent(inout) :: C(*)
115 
116  integer :: info
117 #ifdef DEV_LINALG_TIMING
118  real(dp) :: tsec(2)
119  call timab(TIMAB_XGEMM,1,tsec)
120 #endif
121 
122  if (ABI_LINALG_PLASMA_ISON) then
123    info = -1
124 #ifdef HAVE_LINALG_PLASMA
125    !write(std_out,*)"Will call PLASMA_zgemm_c"
126    info = PLASMA_zgemm_c(trans_plasma(TRANSA),trans_plasma(TRANSB),M,N,K,ALPHA,&
127 &    c_loc(A),LDA,c_loc(B),LDB,BETA,c_loc(C),LDC)
128 #endif
129    ABI_CHECK(info==0,"PLASMA_zgemm_c returned info !=0")
130  else
131    if (use_zgemm3m(m,n,k)) then
132      call _ZGEMM3M(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
133    else
134      call zgemm(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
135    end if
136  end if
137 
138 #ifdef DEV_LINALG_TIMING
139  call timab(TIMAB_XGEMM,2,tsec)
140 #endif
141 
142 end subroutine abi_zgemm_2r