The Elk Code
 
Loading...
Searching...
No Matches
writemomqu.f90
Go to the documentation of this file.
1
2! Copyright (C) 2025 Wenhan Chen, J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine writemomqu
7use modmain
8use modulr
9implicit none
10! local variables
11integer idm,is,ia,ias,iq,ifq
12! automatic arrays
13complex(8) zfft(nqpt),mmtqu(ndmag,natmtot,nqpt)
14complex(8) mirqu(ndmag,nqpt),mtotqu(ndmag,nqpt)
15do idm=1,ndmag
16! muffin-tin moments
17 do ias=1,natmtot
18 is=idxis(ias)
19 zfft(1:nqpt)=mommtru(idm,ias,1:nqpt)
20 call zfftifc(3,ngridq,-1,zfft)
21 mmtqu(idm,ias,1:nqpt)=zfft(1:nqpt)
22 end do
23! interstitial moment
24 zfft(1:nqpt)=momirru(idm,1:nqpt)
25 call zfftifc(3,ngridq,-1,zfft)
26 mirqu(idm,1:nqpt)=zfft(1:nqpt)
27! total moment
28 zfft(1:nqpt)=momtotru(idm,1:nqpt)
29 call zfftifc(3,ngridq,-1,zfft)
30 mtotqu(idm,1:nqpt)=zfft(1:nqpt)
31end do
32open(50,file='MOMENTQU.OUT',form='FORMATTED')
33do iq=1,nqpt
34 ifq=iqfft(iq)
35 write(50,*)
36 write(50,'("Q-point number ",I6," of ",I6)') iq,nqpt
37 write(50,'("Q-point (lattice coordinates) :")')
38 write(50,'(3G18.10)') vql(:,iq)
39 write(50,'("Q-point (Cartesian coordinates) :")')
40 write(50,'(3G18.10)') vqc(:,iq)
41 write(50,'("Moments (complex):")')
42 write(50,'(" interstitial :")')
43 do idm=1,ndmag
44 write(50,'(2G18.10)') mirqu(idm,ifq)
45 end do
46 write(50,'(" muffin-tins")')
47 do is=1,nspecies
48 write(50,'(" species : ",I4," (",A,")")') is,trim(spsymb(is))
49 do ia=1,natoms(is)
50 ias=idxas(ia,is)
51 write(50,'(" atom ",I4," :")') ia
52 do idm=1,ndmag
53 write(50,'(2G18.10)') mmtqu(idm,ias,ifq)
54 end do
55 end do
56 end do
57 write(50,'(" total moment :")')
58 do idm=1,ndmag
59 write(50,'(2G18.10)') mtotqu(idm,ifq)
60 end do
61end do
62close(50)
63end subroutine
64
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(:,:), allocatable vqc
Definition modmain.f90:547
real(8), dimension(:,:), allocatable vql
Definition modmain.f90:545
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
integer, dimension(:), allocatable iqfft
Definition modmain.f90:537
integer nspecies
Definition modmain.f90:34
integer, dimension(3) ngridq
Definition modmain.f90:515
real(8), dimension(:,:), allocatable momirru
Definition modulr.f90:57
real(8), dimension(:,:), allocatable momtotru
Definition modulr.f90:57
real(8), dimension(:,:,:), allocatable mommtru
Definition modulr.f90:57
subroutine writemomqu
Definition writemomqu.f90:7
subroutine zfftifc(nd, n, sgn, z)