The Elk Code
 
Loading...
Searching...
No Matches
writew90unk.f90
Go to the documentation of this file.
1
2! Copyright (C) 2024 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 writew90unk
7use modmain
8use modw90
9use modomp
10implicit none
11! local variables
12integer ispn,ik,ist,nthd
13integer is,ias,nrc,nrci,npc
14integer np,ngp(nspnfv),iu,i
15real(8) vc(3),kc
16character(256) fname
17! automatic arrays
18complex(8) ylmk(lmmaxo),sfack(natmtot)
19! allocatable arrays
20integer, allocatable :: igpig(:,:)
21real(8), allocatable :: vpl(:,:),jlkr(:,:)
22complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
23complex(8), allocatable :: expmt(:,:),wf(:,:,:)
24! total number of plot points
25np=np3d(1)*np3d(2)*np3d(3)
26! generate the 3D plotting points
27allocate(vpl(3,np))
28call plotpt3d(vpl)
29! parallel loop over non-reduced k-points
30call holdthd(nkptnr,nthd)
31!$OMP PARALLEL DEFAULT(SHARED) &
32!$OMP PRIVATE(igpig,jlkr,wfmt,wfir) &
33!$OMP PRIVATE(wf,expmt,ylmk,sfack) &
34!$OMP PRIVATE(ngp,vc,kc,ist,ispn,ias,is) &
35!$OMP PRIVATE(nrc,nrci,npc,fname,iu,i) &
36!$OMP NUM_THREADS(nthd)
37allocate(igpig(ngkmax,nspnfv))
38allocate(jlkr(njcmax,nspecies))
39allocate(wfmt(npcmtmax,natmtot,nspinor,num_bands))
40allocate(wfir(ngtot,nspinor,num_bands))
41allocate(wf(np,nspinor,num_bands))
42allocate(expmt(npcmtmax,natmtot))
43!$OMP DO SCHEDULE(DYNAMIC)
44do ik=1,nkptnr
45!$OMP CRITICAL(writew90unk_)
46 write(*,'("Info(writew90unk): ",I6," of ",I6," k-points")') ik,nkptnr
47!$OMP END CRITICAL(writew90unk_)
48! generate the second-variational wavefunctions
49 call genwfsvp(.false.,.false.,num_bands,idxw90,ngridg,igfft,vkl(:,ik),ngp, &
50 igpig,wfmt,ngtot,wfir)
51! generate the phase factor function exp(-ik.r) in the muffin-tins
52 vc(:)=-vkc(:,ik)
53 kc=sqrt(vc(1)**2+vc(2)**2+vc(3)**2)
54 call genjlgpr(1,kc,jlkr)
55 call genylmv(.true.,lmaxo,vc,ylmk)
56 call gensfacgp(1,vc,1,sfack)
57 call genexpmt(1,jlkr,ylmk,1,sfack,expmt)
58! split the wavefunctions into real and imaginary parts
59 do ist=1,num_bands
60 do ispn=1,nspinor
61 do ias=1,natmtot
62 is=idxis(ias)
63 nrc=nrcmt(is)
64 nrci=nrcmti(is)
65 npc=npcmt(is)
66! remove the explicit phase exp(ik.r) from the muffin-tin wavefunction
67 wfmt(1:npc,ias,ispn,ist)=wfmt(1:npc,ias,ispn,ist)*expmt(1:npc,ias)
68 call zfshtip(nrc,nrci,wfmt(:,ias,ispn,ist))
69 end do
70! generate the wavefunctions on a regular grid
71 call zfpts(np,vpl,wfmt(:,:,ispn,ist),wfir(:,ispn,ist),wf(:,ispn,ist))
72 end do
73 end do
74 if (spinpol) then
75 write(fname,'("UNK",I5.5,".NC")') ik
76 else
77 write(fname,'("UNK",I5.5,".1")') ik
78 end if
79 open(newunit=iu,file=trim(fname),form='UNFORMATTED',action='WRITE')
80 write(iu) np3d(1),np3d(2),np3d(3),ik,num_bands
81 do ist=1,num_bands
82 write(iu) (wf(i,1,ist),i=1,np)
83 if (spinpol) then
84 write(iu) (wf(i,2,ist),i=1,np)
85 end if
86 end do
87 close(iu)
88end do
89!$OMP END DO
90deallocate(igpig,jlkr,wfmt,wfir,wf,expmt)
91!$OMP END PARALLEL
92call freethd(nthd)
93write(*,*)
94write(*,'("Info(writew90unk): created the UNKkkkkk.s files")')
95deallocate(vpl)
96end subroutine
97
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition genexpmt.f90:7
subroutine genjlgpr(ngp, gpc, jlgpr)
Definition genjlgpr.f90:7
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition gensfacgp.f90:10
subroutine genwfsvp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition genwfsvp.f90:7
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition genylmv.f90:10
integer njcmax
Definition modmain.f90:1170
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
integer nkptnr
Definition modmain.f90:463
integer nspinor
Definition modmain.f90:267
logical spinpol
Definition modmain.f90:228
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
integer, dimension(3) np3d
Definition modmain.f90:1131
integer lmaxo
Definition modmain.f90:201
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
integer, dimension(:), allocatable idxw90
Definition modw90.f90:22
integer num_bands
Definition modw90.f90:20
pure subroutine plotpt3d(vpl)
Definition plotpt3d.f90:7
subroutine writew90unk
subroutine zfpts(np, vrl, zfmt, zfir, fp)
Definition zfpts.f90:7
subroutine zfshtip(nr, nri, zfmt)
Definition zfshtip.f90:7