The Elk Code
writew90amn.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 writew90amn
7 use modmain
8 use modw90
9 use modrandom
10 use modomp
11 implicit none
12 ! local variables
13 integer ik,ist,i,nthd
14 real(8) a,b
15 character(256) fname
16 ! allocatable arrays
17 complex(8), allocatable :: amn(:,:,:)
18 ! automatically generate the projections from the wavefunctions if required
19 if (projw90) then
20  allocate(amn(num_bands,num_wann,nkptnr))
21  call holdthd(nkptnr,nthd)
22 !$OMP PARALLEL DO DEFAULT(SHARED) &
23 !$OMP SCHEDULE(DYNAMIC) &
24 !$OMP NUM_THREADS(nthd)
25  do ik=1,nkptnr
26  call projkw90(ik,amn(:,:,ik))
27  end do
28 !$OMP END PARALLEL DO
29  call freethd(nthd)
30 end if
31 ! write the projections to file
32 fname=trim(seedname)//'.amn'
33 open(50,file=trim(fname),action='WRITE',form='FORMATTED')
34 write(50,'("Generated by Elk version ",I0,".",I0,".",I0)') version
35 write(50,'(3I6)') num_bands,nkptnr,num_wann
36 do ik=1,nkptnr
37  do i=1,num_wann
38  do ist=1,num_bands
39  if (projw90) then
40 ! calculated projections
41  a=amn(ist,i,ik)%re
42  b=amn(ist,i,ik)%im
43  else
44 ! random projections
45  a=0.5d0-randomu()
46  b=0.5d0-randomu()
47  end if
48  write(50,'(3I8,2G18.10)') ist,i,ik,a,b
49  end do
50  end do
51 end do
52 close(50)
53 write(*,*)
54 write(*,'("Info(writew90amn): created file ",A)') trim(fname)
55 if (projw90) deallocate(amn)
56 end subroutine
57 
integer num_bands
Definition: modw90.f90:22
logical projw90
Definition: modw90.f90:26
Definition: modomp.f90:6
integer num_wann
Definition: modw90.f90:20
integer nkptnr
Definition: modmain.f90:463
character(256) seedname
Definition: modw90.f90:14
Definition: modw90.f90:6
real(8) function randomu()
Definition: modrandom.f90:18
integer, dimension(3), parameter version
Definition: modmain.f90:1276
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine projkw90(ik, amn)
Definition: projkw90.f90:7
subroutine writew90amn
Definition: writew90amn.f90:7