The Elk Code
rhomagq.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 rhomagq
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! local variables
12 integer ifq,idm,is,ias,nthd
13 ! partial Fourier transform of density to Q-space
15 ! convert density to spherical harmonics
16 call holdthd(nfqrz,nthd)
17 !$OMP PARALLEL DO DEFAULT(SHARED) &
18 !$OMP PRIVATE(ias,is) &
19 !$OMP SCHEDULE(DYNAMIC) &
20 !$OMP NUM_THREADS(nthd)
21 do ifq=1,nfqrz
22  do ias=1,natmtot
23  is=idxis(ias)
24  call zfshtip(nrcmt(is),nrcmti(is),rhoqmt(:,ias,ifq))
25  end do
26 end do
27 !$OMP END PARALLEL DO
28 call freethd(nthd)
29 if (.not.spinpol) return
30 ! partial Fourier transform of magnetisation to Q-space
32 ! convert magnetisation to spherical harmonics
33 call holdthd(nfqrz,nthd)
34 !$OMP PARALLEL DO DEFAULT(SHARED) &
35 !$OMP PRIVATE(idm,ias,is) &
36 !$OMP SCHEDULE(DYNAMIC) &
37 !$OMP NUM_THREADS(nthd)
38 do ifq=1,nfqrz
39  do idm=1,ndmag
40  do ias=1,natmtot
41  is=idxis(ias)
42  call zfshtip(nrcmt(is),nrcmti(is),magqmt(:,ias,idm,ifq))
43  end do
44  end do
45 end do
46 !$OMP END PARALLEL DO
47 call freethd(nthd)
48 end subroutine
49 
integer ngtc
Definition: modmain.f90:392
logical spinpol
Definition: modmain.f90:228
integer ndmag
Definition: modmain.f90:238
Definition: modomp.f90:6
complex(8), dimension(:,:,:,:), allocatable magqmt
Definition: modulr.f90:61
complex(8), dimension(:,:,:), allocatable magqir
Definition: modulr.f90:61
real(8), dimension(:,:,:,:), pointer, contiguous magrmt
Definition: modulr.f90:54
subroutine zfshtip(nr, nri, zfmt)
Definition: zfshtip.f90:7
subroutine rhomagq
Definition: rhomagq.f90:7
real(8), dimension(:,:,:), pointer, contiguous rhormt
Definition: modulr.f90:53
real(8), dimension(:,:), pointer, contiguous rhorir
Definition: modulr.f90:53
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer nfqrz
Definition: modmain.f90:539
complex(8), dimension(:,:,:), allocatable rhoqmt
Definition: modulr.f90:60
real(8), dimension(:,:,:), pointer, contiguous magrir
Definition: modulr.f90:54
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
subroutine rfzfftq(sgn, nf, ngt, rfmt, rfir, zfmt, zfir)
Definition: rfzfftq.f90:7
Definition: modulr.f90:6
complex(8), dimension(:,:), allocatable rhoqir
Definition: modulr.f90:60