The Elk Code
momentu.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 T. Mueller, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine momentu
7 use modmain
8 use modulr
9 use modtest
10 implicit none
11 ! local variables
12 integer idm,is,ias
13 integer nrc,nrci,ifq,ir
14 real(8) t1
15 ! automatic arrays
16 real(8) rfft(nqpt)
17 complex(8) zfft(nfqrz)
18 ! external functions
19 complex(8), external :: zfmtint
20 if (.not.spinpol) return
21 ! calculate muffin-tin moments
22 mommttot(:)=0.d0
23 do idm=1,ndmag
24  do ias=1,natmtot
25  is=idxis(ias)
26  nrc=nrcmt(is)
27  nrci=nrcmti(is)
28  do ifq=1,nfqrz
29  zfft(ifq)=zfmtint(nrc,nrci,wr2cmt(:,is),magqmt(:,ias,idm,ifq))
30  end do
31  mommt(idm,ias)=dble(zfft(1))
32  mommttot(idm)=mommttot(idm)+mommt(idm,ias)
33  call rzfftifc(3,ngridq,1,rfft,zfft)
34  mommtru(idm,ias,1:nqpt)=rfft(1:nqpt)
35  end do
36 end do
37 ! find the interstitial and total moments
38 t1=omega/dble(ngtc)
39 do idm=1,ndmag
40  do ifq=1,nfqrz
41  zfft(ifq)=sum(magqir(1:ngtc,idm,ifq)*cfrc(1:ngtc))
42  end do
43  momir(idm)=t1*dble(zfft(1))
44  momtot(idm)=mommttot(idm)+momir(idm)
45  call rzfftifc(3,ngridq,1,rfft,zfft)
46  do ir=1,nqpt
47  momirru(idm,ir)=t1*rfft(ir)
48  momtotru(idm,ir)=sum(mommtru(idm,1:natmtot,ir))+momirru(idm,ir)
49  end do
50 end do
51 ! total moment magnitude
52 if (ncmag) then
53  momtotm=sqrt(momtot(1)**2+momtot(2)**2+momtot(3)**2)
54 else
55  momtotm=abs(momtot(1))
56 end if
57 ! write the muffin-tin moments to test file
58 call writetest(770,'ULR muffin-tin moments',nv=ndmag*natmtot*nqpt,tol=2.d-2, &
59  rva=mommtru)
60 end subroutine
61 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
real(8), dimension(:,:), allocatable mommt
Definition: modmain.f90:744
real(8), dimension(3) momtot
Definition: modmain.f90:738
integer ngtc
Definition: modmain.f90:392
logical spinpol
Definition: modmain.f90:228
integer ndmag
Definition: modmain.f90:238
real(8) momtotm
Definition: modmain.f90:740
real(8), dimension(3) mommttot
Definition: modmain.f90:746
real(8) omega
Definition: modmain.f90:20
subroutine momentu
Definition: momentu.f90:7
complex(8), dimension(:,:,:,:), allocatable magqmt
Definition: modulr.f90:61
complex(8), dimension(:,:,:), allocatable magqir
Definition: modulr.f90:61
real(8), dimension(:,:), allocatable momirru
Definition: modulr.f90:58
real(8), dimension(3) momir
Definition: modmain.f90:742
integer, dimension(3) ngridq
Definition: modmain.f90:515
real(8), dimension(:,:), allocatable wr2cmt
Definition: modmain.f90:189
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:,:), allocatable mommtru
Definition: modulr.f90:58
subroutine rzfftifc(nd, n, sgn, r, z)
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:), allocatable cfrc
Definition: modmain.f90:438
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
Definition: modulr.f90:6
logical ncmag
Definition: modmain.f90:240
real(8), dimension(:,:), allocatable momtotru
Definition: modulr.f90:58