The Elk Code
 
Loading...
Searching...
No Matches
nesting.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 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 nesting
7use modmain
8use modomp
9implicit none
10! local variables
11integer iq,ik,jk,jkq,ivkq(3)
12integer ist,i1,i2,i3,nthd
13real(8) sm0,sm1,sm2,sm3
14real(8) vl(3),vc(3),x,t1
15! allocatable arrays
16real(8), allocatable :: nq(:)
17! external functions
18real(8), external :: sdelta
19! initialise universal variables
20call init0
21call init1
22call init2
23! read Fermi energy from file
24call readfermi
25! get the eigenvalues from file
26do ik=1,nkpt
27 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
28end do
29allocate(nq(nqpt))
30t1=1.d0/swidth
31sm0=0.d0
32call holdthd(nqpt,nthd)
33!$OMP PARALLEL DO DEFAULT(SHARED) &
34!$OMP PRIVATE(sm1,sm2,sm3,ik,jk) &
35!$OMP PRIVATE(ivkq,jkq,ist,x) &
36!$OMP REDUCTION(+:sm0) &
37!$OMP SCHEDULE(DYNAMIC) &
38!$OMP NUM_THREADS(nthd)
39do iq=1,nqpt
40!$OMP CRITICAL(nesting_)
41 write(*,'("Info(nesting): ",I6," of ",I6," q-points")') iq,nqpt
42!$OMP END CRITICAL(nesting_)
43 sm1=0.d0
44 do ik=1,nkptnr
45 jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
46 ivkq(:)=ivk(:,ik)+ivq(:,iq)
47 ivkq(:)=mod(ivkq(:),ngridk(:))
48 jkq=ivkik(ivkq(1),ivkq(2),ivkq(3))
49 sm2=0.d0
50 do ist=1,nstsv
51 x=(efermi-evalsv(ist,jk))*t1
52 sm2=sm2+sdelta(stype,x)*t1
53 end do
54 sm3=0.d0
55 do ist=1,nstsv
56 x=(efermi-evalsv(ist,jkq))*t1
57 sm3=sm3+sdelta(stype,x)*t1
58 end do
59 sm1=sm1+sm2*sm3
60 end do
61 nq(iq)=occmax*omegabz*wkptnr*sm1
62 sm0=sm0+omegabz*wqpt(iq)*nq(iq)
63end do
64!$OMP END PARALLEL DO
65call freethd(nthd)
66open(50,file='NEST3D.OUT',form='FORMATTED')
67write(50,'(3I6," : grid size")') ngridq(:)
68do i3=0,ngridq(3)-1
69 vl(3)=dble(i3)/dble(ngridq(3))
70 do i2=0,ngridq(2)-1
71 vl(2)=dble(i2)/dble(ngridq(2))
72 do i1=0,ngridq(1)-1
73 vl(1)=dble(i1)/dble(ngridq(1))
74 vc(:)=bvec(:,1)*vl(1)+bvec(:,2)*vl(2)+bvec(:,3)*vl(3)
75 iq=ivqiq(i1,i2,i3)
76 write(50,'(4G18.10)') vc(:),nq(iq)
77 end do
78 end do
79end do
80close(50)
81open(50,file='NESTING.OUT',form='FORMATTED')
82write(50,'(G18.10)') sm0
83close(50)
84write(*,*)
85write(*,'("Info(nesting):")')
86write(*,'(" Nesting function N(q) written to NEST3D.OUT for plotting")')
87write(*,*)
88write(*,'(" Total integrated nesting per unit volume written to NESTING.OUT")')
89deallocate(nq)
90end subroutine
91
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
real(8) efermi
Definition modmain.f90:904
integer nkptnr
Definition modmain.f90:463
character(256) filext
Definition modmain.f90:1300
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer, dimension(:,:), allocatable ivq
Definition modmain.f90:529
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
real(8) swidth
Definition modmain.f90:892
integer nqpt
Definition modmain.f90:525
integer nkpt
Definition modmain.f90:461
real(8) omegabz
Definition modmain.f90:22
real(8), dimension(:), allocatable wqpt
Definition modmain.f90:549
integer, dimension(3) ngridk
Definition modmain.f90:448
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) occmax
Definition modmain.f90:898
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
integer stype
Definition modmain.f90:888
integer, dimension(3) ngridq
Definition modmain.f90:515
real(8) wkptnr
Definition modmain.f90:477
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine nesting
Definition nesting.f90:7
subroutine readfermi
Definition readfermi.f90:10
real(8) function sdelta(stype, x)
Definition sdelta.f90:10