The Elk Code
 
Loading...
Searching...
No Matches
dbxcplot.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007 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 dbxcplot
7use modmain
8implicit none
9! local variables
10integer idm,is,ias,np
11! allocatable arrays
12real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
13real(8), allocatable :: rfmt(:,:),rfir(:)
14real(8), allocatable :: grfmt(:,:,:),grfir(:,:)
15! initialise universal variables
16call init0
17if (.not.spinpol) then
18 write(*,*)
19 write(*,'("Error(dbxcplot): spin-unpolarised magnetic field is zero")')
20 write(*,*)
21 stop
22end if
23! read magnetisation from file
24call readstate
25allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
26allocate(rfmt(npmtmax,natmtot),rfir(ngtot))
27allocate(grfmt(npmtmax,natmtot,3),grfir(ngtot,3))
28if (ncmag) then
29! non-collinear
30 rvfmt(:,:,:)=bxcmt(:,:,:)
31 rvfir(:,:)=bxcir(:,:)
32else
33! collinear
34 rvfmt(:,:,1:2)=0.d0
35 rvfir(:,1:2)=0.d0
36 rvfmt(:,:,3)=bxcmt(:,:,1)
37 rvfir(:,3)=bxcir(:,1)
38end if
39rfmt(:,:)=0.d0
40rfir(:)=0.d0
41do idm=1,3
42 call gradrf(rvfmt(:,:,idm),rvfir(:,idm),grfmt,grfir)
43 do ias=1,natmtot
44 is=idxis(ias)
45 np=npmt(is)
46 rfmt(1:np,ias)=rfmt(1:np,ias)+grfmt(1:np,ias,idm)
47 end do
48 rfir(:)=rfir(:)+grfir(:,idm)
49end do
50select case(task)
51case(91)
52 open(50,file='DBXC1D.OUT',form='FORMATTED')
53 open(51,file='DBXCLINES.OUT',form='FORMATTED')
54 call plot1d(50,51,1,rfmt,rfir)
55 close(50)
56 close(51)
57 write(*,*)
58 write(*,'("Info(dbxcplot):")')
59 write(*,'(" 1D divergence of exchange-correlation field written to &
60 &DBXC1D.OUT")')
61 write(*,'(" vertex location lines written to DBXCLINES.OUT")')
62case(92)
63 open(50,file='DBXC2D.OUT',form='FORMATTED')
64 call plot2d(.false.,50,1,rfmt,rfir)
65 close(50)
66 write(*,'("Info(dbxcplot):")')
67 write(*,'(" 2D divergence of exchange-correlation field written to &
68 &DBXC2D.OUT")')
69case(93)
70 open(50,file='DBXC3D.OUT',form='FORMATTED')
71 call plot3d(50,1,rfmt,rfir)
72 close(50)
73 write(*,'("Info(dbxcplot):")')
74 write(*,'(" 3D divergence of exchange-correlation field written to &
75 &DBXC3D.OUT")')
76end select
77deallocate(rvfmt,rvfir,rfmt,rfir,grfmt,grfir)
78end subroutine
79
subroutine dbxcplot
Definition dbxcplot.f90:7
subroutine gradrf(rfmt, rfir, grfmt, grfir)
Definition gradrf.f90:7
subroutine init0
Definition init0.f90:10
real(8), dimension(:,:,:), allocatable bxcmt
Definition modmain.f90:636
integer ngtot
Definition modmain.f90:390
logical ncmag
Definition modmain.f90:240
logical spinpol
Definition modmain.f90:228
real(8), dimension(:,:), allocatable bxcir
Definition modmain.f90:636
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer natmtot
Definition modmain.f90:40
integer task
Definition modmain.f90:1298
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
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