The Elk Code
aceplot.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2021 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 
6 subroutine aceplot
7 use modmain
8 use modphonon
9 use modbog
10 implicit none
11 ! local variables
12 integer ik,jk,iq,jq
13 integer i1,i2,i3,ist,i
14 real(8) ace,vn,xn
15 ! initialise universal variables
16 call init0
17 call init1
18 call init2
19 call initeph
20 !----------------------------------------------------------!
21 ! plot the fermionic anomalous correlation entropy !
22 !----------------------------------------------------------!
23 open(50,file='FACE3D.OUT',form='FORMATTED',action='WRITE')
24 write(50,'(3I6," : grid size")') ngridk(:)
25 do i3=0,ngridk(3)-1
26  do i2=0,ngridk(2)-1
27  do i1=0,ngridk(1)-1
28  ik=ivkiknr(i1,i2,i3)
29  jk=ivkik(i1,i2,i3)
30  ace=0.d0
31  do ist=1,nstsv
32  vn=vnorm(ist,jk)
33  if ((vn > 0.d0).and.(vn < 1.d0)) then
34  ace=ace+vn*log(vn)+(1.d0-vn)*log(1.d0-vn)
35  end if
36  end do
37  ace=-occmax*ace
38  write(50,'(4G18.10)') vkc(:,ik),ace
39  end do
40  end do
41 end do
42 close(50)
43 !--------------------------------------------------------!
44 ! plot the bosonic anomalous correlation entropy !
45 !--------------------------------------------------------!
46 open(50,file='BACE3D.OUT',form='FORMATTED',action='WRITE')
47 write(50,'(3I6," : grid size")') ngridq(:)
48 do i3=0,ngridq(3)-1
49  do i2=0,ngridq(2)-1
50  do i1=0,ngridq(1)-1
51  iq=ivqiqnr(i1,i2,i3)
52  jq=ivqiq(i1,i2,i3)
53  ace=0.d0
54  do i=1,nbph
55  xn=xnorm(i,jq)
56  if (xn > 0.d0) then
57  ace=ace+xn*log(xn)-(1.d0+xn)*log(1.d0+xn)
58  end if
59  end do
60  ace=-ace
61  write(50,'(4G18.10)') vqc(:,iq),ace
62  end do
63  end do
64 end do
65 close(50)
66 write(*,*)
67 write(*,'("Info(aceplot):")')
68 write(*,'(" 3D fermionic anomalous correlation entropy written to FACE3D.OUT")')
69 write(*,'(" 3D bosonic anomalous correlation entropy written to BACE3D.OUT")')
70 end subroutine
71 
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
real(8), dimension(:,:), allocatable vkc
Definition: modmain.f90:473
integer nstsv
Definition: modmain.f90:889
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
real(8), dimension(:,:), allocatable xnorm
Definition: modbog.f90:41
real(8) occmax
Definition: modmain.f90:901
Definition: modbog.f90:6
integer nbph
Definition: modphonon.f90:13
subroutine init2
Definition: init2.f90:7
integer, dimension(3) ngridk
Definition: modmain.f90:448
integer, dimension(:,:,:), allocatable ivqiq
Definition: modmain.f90:531
subroutine init1
Definition: init1.f90:10
integer, dimension(:,:,:), allocatable ivqiqnr
Definition: modmain.f90:533
integer, dimension(3) ngridq
Definition: modmain.f90:515
real(8), dimension(:,:), allocatable vnorm
Definition: modbog.f90:17
subroutine initeph
Definition: initeph.f90:7
subroutine init0
Definition: init0.f90:10
subroutine aceplot
Definition: aceplot.f90:7
integer, dimension(:,:,:), allocatable ivkiknr
Definition: modmain.f90:469