The Elk Code
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 
6 subroutine writew90unk
7 use modmain
8 use modw90
9 use modomp
10 implicit none
11 ! local variables
12 integer ispn,ik,ist,nthd
13 integer is,ias,nrc,nrci,npc
14 integer np,ngp(nspnfv),iu,i
15 real(8) vc(3),kc
16 character(256) fname
17 ! automatic arrays
18 complex(8) ylmk(lmmaxo),sfack(natmtot)
19 ! allocatable arrays
20 integer, allocatable :: igpig(:,:)
21 real(8), allocatable :: vpl(:,:),jlkr(:,:)
22 complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
23 complex(8), allocatable :: expmt(:,:),wf(:,:,:)
24 ! total number of plot points
25 np=np3d(1)*np3d(2)*np3d(3)
26 ! generate the 3D plotting points
27 allocate(vpl(3,np))
28 call plotpt3d(vpl)
29 ! parallel loop over non-reduced k-points
30 call 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)
37 allocate(igpig(ngkmax,nspnfv))
38 allocate(jlkr(njcmax,nspecies))
39 allocate(wfmt(npcmtmax,natmtot,nspinor,num_bands))
40 allocate(wfir(ngtot,nspinor,num_bands))
41 allocate(wf(np,nspinor,num_bands))
42 allocate(expmt(npcmtmax,natmtot))
43 !$OMP DO SCHEDULE(DYNAMIC)
44 do 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)
88 end do
89 !$OMP END DO
90 deallocate(igpig,jlkr,wfmt,wfir,wf,expmt)
91 !$OMP END PARALLEL
92 call freethd(nthd)
93 write(*,*)
94 write(*,'("Info(writew90unk): created the UNKkkkkk.s files")')
95 deallocate(vpl)
96 end subroutine
97 
integer num_bands
Definition: modw90.f90:20
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer npcmtmax
Definition: modmain.f90:216
integer, dimension(3) ngridg
Definition: modmain.f90:386
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
integer ngtot
Definition: modmain.f90:390
logical spinpol
Definition: modmain.f90:228
integer ngkmax
Definition: modmain.f90:499
Definition: modomp.f90:6
integer lmaxo
Definition: modmain.f90:201
integer nkptnr
Definition: modmain.f90:463
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
subroutine writew90unk
Definition: writew90unk.f90:7
real(8), dimension(:,:), allocatable vkc
Definition: modmain.f90:473
pure subroutine plotpt3d(vpl)
Definition: plotpt3d.f90:7
subroutine zfpts(np, vrl, zfmt, zfir, fp)
Definition: zfpts.f90:7
Definition: modw90.f90:6
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine zfshtip(nr, nri, zfmt)
Definition: zfshtip.f90:7
integer, dimension(:), allocatable idxw90
Definition: modw90.f90:22
integer nspinor
Definition: modmain.f90:267
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition: genexpmt.f90:7
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine genjlgpr(ngp, gpc, jlgpr)
Definition: genjlgpr.f90:7
integer nspecies
Definition: modmain.f90:34
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer njcmax
Definition: modmain.f90:1173
integer, dimension(3) np3d
Definition: modmain.f90:1134
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
subroutine genwfsvp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition: genwfsvp.f90:7