The Elk Code
 
Loading...
Searching...
No Matches
wfcrplot.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 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 wfcrplot
7use modmain
8implicit none
9! local variables
10integer ist,is,ia,ias,ir
11character(256) fname
12! initialise universal variables
13call init0
14! read density and potentials from file
15call readstate
16! generate the core wavefunctions
17call gencore
18do is=1,nspecies
19 do ia=1,natoms(is)
20 ias=idxas(ia,is)
21 write(fname,'("WFCORE_S",I2.2,"_A",I4.4,".OUT")') is,ia
22 open(50,file=trim(fname),form='FORMATTED')
23 do ist=1,nstsp(is)
24 if (spcore(ist,is)) then
25 do ir=1,nrsp(is)
26 write(50,'(2G18.10)') rsp(ir,is),rwfcr(ir,1,ist,ias)
27 end do
28 write(50,*)
29 end if
30 end do
31 close(50)
32 end do
33end do
34write(*,*)
35write(*,'("Info(wfcrplot):")')
36write(*,'(" Core state wavefunctions written to WFCORE_Sss_Aaaaa.OUT")')
37write(*,'(" for all species and atoms")')
38end subroutine
39
subroutine gencore
Definition gencore.f90:10
subroutine init0
Definition init0.f90:10
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
real(8), dimension(:,:,:,:), allocatable rwfcr
Definition modmain.f90:936
integer, dimension(maxspecies) nrsp
Definition modmain.f90:107
integer nspecies
Definition modmain.f90:34
subroutine readstate
Definition readstate.f90:10
subroutine wfcrplot
Definition wfcrplot.f90:7