File: momentu.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (34 lines) | stat: -rw-r--r-- 796 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34

! Copyright (C) 2018 T. Mueller, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine momentu
use modmain
use modulr
implicit none
! local variables
integer iq,idm,is,ias
! allocatable arrays
complex(8), allocatable :: zfft(:)
! external functions
complex(8) zfmtint
external zfmtint
if (.not.spinpol) return
allocate(zfft(nqpt))
do ias=1,natmtot
  is=idxis(ias)
  do idm=1,ndmag
    do iq=1,nqpt
      zfft(iq)=zfmtint(nrcmt(is),nrcmti(is),rcmt(:,is),r2cmt(:,is), &
       magqmt(:,ias,idm,iq))
    end do
    mommtu(idm,ias)=dble(zfft(1))
    call zfftifc(3,ngridq,1,zfft)
    mommtru(idm,ias,:)=dble(zfft(:))
  end do
end do
deallocate(zfft)
return
end subroutine