The Elk Code
 
Loading...
Searching...
No Matches
gentaucr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 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 gentaucr
7use modmain
8use modomp
9implicit none
10! local variables
11integer ist,ispn,jspn
12integer is,ia,ias,nthd
13integer nr,nri,np,i,m
14! allocatable arrays
15complex(4), allocatable :: wfcr(:,:)
16complex(8), allocatable :: zfmt(:),gzfmt(:,:)
17taucr(1:npmtmax,1:natmtot,1:nspinor)=0.d0
18call holdthd(natmtot,nthd)
19!$OMP PARALLEL DEFAULT(SHARED) &
20!$OMP PRIVATE(wfcr,zfmt,gzfmt) &
21!$OMP PRIVATE(is,ia,nr,nri,np) &
22!$OMP PRIVATE(ist,m,ispn,jspn,i) &
23!$OMP NUM_THREADS(nthd)
24allocate(wfcr(npmtmax,2),zfmt(npmtmax),gzfmt(npmtmax,3))
25!$OMP DO SCHEDULE(DYNAMIC)
26do ias=1,natmtot
27 is=idxis(ias)
28 ia=idxia(ias)
29 nr=nrmt(is)
30 nri=nrmti(is)
31 np=npmt(is)
32 do ist=1,nstsp(is)
33 if (spcore(ist,is)) then
34 do m=-ksp(ist,is),ksp(ist,is)-1
35! generate the core wavefunction in spherical harmonics (pass in m-1/2)
36 call wavefcr(.true.,1,is,ia,ist,m,npmtmax,wfcr)
37 do ispn=1,2
38 if (spinpol) then
39 jspn=ispn
40 else
41 jspn=1
42 end if
43! compute the gradient of the wavefunction
44 zfmt(1:np)=wfcr(1:np,ispn)
45 call gradzfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),zfmt,npmtmax,gzfmt)
46 do i=1,3
47! convert gradient to spherical coordinates
48 call zbsht(nr,nri,gzfmt(:,i),zfmt)
49! add to total in muffin-tin
50 taucr(1:np,ias,jspn)=taucr(1:np,ias,jspn) &
51 +0.5d0*(dble(zfmt(1:np))**2+aimag(zfmt(1:np))**2)
52 end do
53 end do
54 end do
55 end if
56 end do
57! convert core tau to spherical harmonics
58 do ispn=1,nspinor
59 call rfshtip(nr,nri,taucr(:,ias,ispn))
60 end do
61end do
62!$OMP END DO
63deallocate(wfcr,zfmt,gzfmt)
64!$OMP END PARALLEL
65call freethd(nthd)
66end subroutine
67
subroutine gentaucr
Definition gentaucr.f90:7
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer nspinor
Definition modmain.f90:267
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
logical spinpol
Definition modmain.f90:228
real(8), dimension(:,:,:), allocatable taucr
Definition modmain.f90:674
integer, dimension(maxatoms *maxspecies) idxia
Definition modmain.f90:45
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
integer natmtot
Definition modmain.f90:40
integer, dimension(maxstsp, maxspecies) ksp
Definition modmain.f90:125
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rfshtip(nr, nri, rfmt)
Definition rfshtip.f90:7
pure subroutine wavefcr(tsh, lrstp, is, ia, ist, m, ld, wfcr)
Definition wavefcr.f90:7
subroutine zbsht(nr, nri, zfmt1, zfmt2)
Definition zbsht.f90:10