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