TABLE OF CONTENTS
ABINIT/getghcnd [ Functions ]
NAME
getghcnd
FUNCTION
Compute <G|H_ND|C> for input vector |C> expressed in reciprocal space Result is put in array ghcnc. H_ND is the Hamiltonian due to magnetic dipoles on the nuclear sites.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
cwavef(2,npw*nspinor*ndat)=planewave coefficients of wavefunction. gs_ham <type(gs_hamiltonian_type)>=all data for the Hamiltonian to be applied my_nspinor=number of spinorial components of the wavefunctions (on current proc) ndat=number of FFT to do in //
OUTPUT
ghcnd(2,npw*my_nspinor*ndat)=matrix elements <G|H_ND|C>
SIDE EFFECTS
NOTES
Application of <k^prime|H|k> or <k|H|k^prime> not implemented!
PARENTS
getghc
CHILDREN
zhpmv
NOTES
This routine applies the Hamiltonian due to an array of magnetic dipoles located at the atomic nuclei to the input wavefunction. Strategy below is to take advantage of Hermiticity to store H_ND in triangular form and then use a BLAS call to zhpmv to apply to input vector in one shot.
SOURCE
47 #if defined HAVE_CONFIG_H 48 #include "config.h" 49 #endif 50 51 #include "abi_common.h" 52 53 subroutine getghcnd(cwavef,ghcnd,gs_ham,my_nspinor,ndat) 54 55 use defs_basis 56 use defs_abitypes 57 use m_errors 58 use m_profiling_abi 59 use m_xmpi 60 61 use m_hamiltonian, only : gs_hamiltonian_type 62 63 !This section has been created automatically by the script Abilint (TD). 64 !Do not modify the following lines by hand. 65 #undef ABI_FUNC 66 #define ABI_FUNC 'getghcnd' 67 !End of the abilint section 68 69 implicit none 70 71 !Arguments ------------------------------------ 72 !scalars 73 integer,intent(in) :: my_nspinor,ndat 74 type(gs_hamiltonian_type),intent(in),target :: gs_ham 75 !arrays 76 real(dp),intent(in) :: cwavef(2,gs_ham%npw_k*my_nspinor*ndat) 77 real(dp),intent(out) :: ghcnd(2,gs_ham%npw_k*my_nspinor*ndat) 78 79 !Local variables------------------------------- 80 !scalars 81 integer :: cwavedim 82 character(len=500) :: message 83 !arrays 84 complex(dpc),allocatable :: inwave(:),hggc(:) 85 86 ! ********************************************************************* 87 88 if (gs_ham%matblk /= gs_ham%natom) then 89 write(message,'(a,i4,a,i4)')' gs_ham%matblk = ',gs_ham%matblk,' but natom = ',gs_ham%natom 90 MSG_ERROR(message) 91 end if 92 if (ndat /= 1) then 93 write(message,'(a,i4,a)')' ndat = ',ndat,' but getghcnd requires ndat = 1' 94 MSG_ERROR(message) 95 end if 96 if (my_nspinor /= 1) then 97 write(message,'(a,i4,a)')' nspinor = ',my_nspinor,' but getghcnd requires nspinor = 1' 98 MSG_ERROR(message) 99 end if 100 if (any(abs(gs_ham%kpt_k(:)-gs_ham%kpt_kp(:))>tol8)) then 101 message=' not allowed for kpt(left)/=kpt(right)!' 102 MSG_BUG(message) 103 end if 104 105 cwavedim = gs_ham%npw_k*my_nspinor*ndat 106 ABI_ALLOCATE(hggc,(cwavedim)) 107 ABI_ALLOCATE(inwave,(cwavedim)) 108 109 inwave(1:gs_ham%npw_k) = cmplx(cwavef(1,1:gs_ham%npw_k),cwavef(2,1:gs_ham%npw_k),kind=dpc) 110 111 ! apply hamiltonian hgg to input wavefunction inwave, result in hggc 112 ! ZHPMV is a level-2 BLAS routine, does Matrix x Vector multiplication for double complex 113 ! objects, with the matrix as Hermitian in packed storage 114 call ZHPMV('L',cwavedim,cone,gs_ham%nucdipmom_k,inwave,1,czero,hggc,1) 115 116 ghcnd(1,1:gs_ham%npw_k) = real(hggc) 117 ghcnd(2,1:gs_ham%npw_k) = aimag(hggc) 118 119 ABI_DEALLOCATE(hggc) 120 ABI_DEALLOCATE(inwave) 121 122 end subroutine getghcnd