The Elk Code
 
Loading...
Searching...
No Matches
wfplot.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine wfplot
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,ist
13real(8) x,t1
14! external functions
15real(8), external :: sdelta
16! initialise universal variables
17call init0
18call init1
19! read the density and potentials from file
20call readstate
21! find the new linearisation energies
22call linengy
23! generate the APW radial functions
24call genapwfr
25! generate the local-orbital radial functions
26call genlofr
27! set the occupation numbers
28if (any(task == [61,62,63])) then
29! plot of wavefunction modulus squared
30 ik=kstlist(1,1)
31 ist=kstlist(2,1)
32 if ((ik < 1).or.(ik > nkpt)) then
33 write(*,*)
34 write(*,'("Error(wfplot): k-point out of range : ",I8)') ik
35 write(*,*)
36 stop
37 end if
38 if ((ist < 1).or.(ist > nstsv)) then
39 write(*,*)
40 write(*,'("Error(wfplot): state out of range : ",I8)') ist
41 write(*,*)
42 stop
43 end if
44! select a particular wavefunction using its occupancy
45 occsv(:,:)=0.d0
46 occsv(ist,ik)=1.d0/wkpt(ik)
47! no symmetrisation required
48 nsymcrys=1
49 eqatoms(:,:,:)=.false.
50else
51! plotting an STM image by setting occupation numbers to be a delta function at
52! the Fermi energy
53 t1=1.d0/swidth
54 do ik=1,nkpt
55! get the eigenvalues from file
56 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
57 do ist=1,nstsv
58 x=(efermi-evalsv(ist,ik))*t1
59 occsv(ist,ik)=occmax*wkpt(ik)*sdelta(stype,x)*t1
60 end do
61 end do
62end if
63! compute the valence charge density with the new occupation numbers
64call rhomagv
65! write the wavefunction modulus squared plot to file
66if (mp_mpi) then
67 select case(task)
68 case(61)
69 open(50,file='WF1D.OUT',form='FORMATTED')
70 open(51,file='WFLINES.OUT',form='FORMATTED')
71 call plot1d(50,51,1,rhomt,rhoir)
72 close(50)
73 close(51)
74 write(*,*)
75 write(*,'("Info(wfplot):")')
76 write(*,'(" 1D wavefunction modulus squared written to WF1D.OUT")')
77 write(*,'(" vertex location lines written to WFLINES.OUT")')
78 case(62)
79 open(50,file='WF2D.OUT',form='FORMATTED')
80 call plot2d(.false.,50,1,rhomt,rhoir)
81 close(50)
82 write(*,*)
83 write(*,'("Info(wfplot):")')
84 write(*,'(" 2D wavefunction modulus squared written to WF2D.OUT")')
85 case(162)
86 open(50,file='STM2D.OUT',form='FORMATTED')
87 call plot2d(.false.,50,1,rhomt,rhoir)
88 close(50)
89 write(*,*)
90 write(*,'("Info(wfplot):")')
91 write(*,'(" 2D STM image written to STM2D.OUT")')
92 case(63)
93 open(50,file='WF3D.OUT',form='FORMATTED')
94 call plot3d(50,1,rhomt,rhoir)
95 close(50)
96 write(*,*)
97 write(*,'("Info(wfplot):")')
98 write(*,'(" 3D wavefunction modulus squared written to WF3D.OUT")')
99 end select
100 if (task /= 162) then
101 write(*,'(" for k-point ",I8," and state ",I6)') kstlist(1,1),kstlist(2,1)
102 end if
103end if
104end subroutine
105
subroutine genapwfr
Definition genapwfr.f90:10
subroutine genlofr
Definition genlofr.f90:10
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine linengy
Definition linengy.f90:10
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
real(8) efermi
Definition modmain.f90:904
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
character(256) filext
Definition modmain.f90:1300
real(8) swidth
Definition modmain.f90:892
integer nkpt
Definition modmain.f90:461
integer task
Definition modmain.f90:1298
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) occmax
Definition modmain.f90:898
logical, dimension(:,:,:), allocatable eqatoms
Definition modmain.f90:370
integer nsymcrys
Definition modmain.f90:358
integer stype
Definition modmain.f90:888
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
integer, dimension(2, maxkst) kstlist
Definition modmain.f90:926
logical mp_mpi
Definition modmpi.f90:17
subroutine plot1d(fnum1, fnum2, nf, rfmt, rfir)
Definition plot1d.f90:10
subroutine plot2d(tproj, fnum, nf, rfmt, rfir)
Definition plot2d.f90:10
subroutine plot3d(fnum, nf, rfmt, rfir)
Definition plot3d.f90:10
subroutine readstate
Definition readstate.f90:10
subroutine rhomagv
Definition rhomagv.f90:7
real(8) function sdelta(stype, x)
Definition sdelta.f90:10
subroutine wfplot
Definition wfplot.f90:7