3
! Contains data structures and routines to deal with the kpoint-grid
4
! for the self-consistent calculation
5
! Other uses (bands, optical, polarization) have their own structures.
7
USE precision, only : dp
13
logical, public, save :: scf_kgrid_first_time = .true.
14
logical, public, save :: gamma_scf
15
integer, public, save :: maxk !
16
integer, public, save :: nkpnt ! Total number of k-points
17
real(dp) :: eff_kgrid_cutoff ! Effective kgrid_cutoff
19
real(dp), pointer, public, save :: kweight(:)
20
real(dp), pointer, public, save :: kpoint(:,:)
22
integer, dimension(3,3), save :: kscell = 0
23
real(dp), dimension(3), save :: kdispl = 0.0_dp
25
logical, save :: user_requested_mp = .false.
26
logical, save :: user_requested_cutoff = .false.
28
logical, save :: spiral = .false.
29
logical, save :: firm_displ = .false.
31
public :: setup_kpoint_grid
35
subroutine setup_Kpoint_grid( ucell )
36
USE parallel, only : Node
37
USE fdf, only : fdf_defined
38
USE m_find_kgrid, only : find_kgrid
45
real(dp) :: ucell(3,3)
51
if (scf_kgrid_first_time) then
52
nullify(kweight,kpoint)
54
spiral = fdf_defined('SpinSpiral')
57
call MPI_Bcast(spiral,1,MPI_logical,0,MPI_Comm_World,MPIerror)
59
call setup_scf_kscell(ucell, firm_displ)
61
scf_kgrid_first_time = .false.
64
if ( user_requested_mp ) then
65
! no need to set up the kscell again
67
! This was wrong in the old code
68
call setup_scf_kscell(ucell, firm_displ)
72
call find_kgrid(ucell,kscell,kdispl,firm_displ, &
74
nkpnt,kpoint,kweight, eff_kgrid_cutoff)
77
gamma_scf = (nkpnt == 1 .and. &
78
dot_product(kpoint(:,1),kpoint(:,1)) < 1.0e-20_dp)
80
if (Node .eq. 0) call siesta_write_k_points()
82
end subroutine setup_Kpoint_grid
84
!--------------------------------------------------------------------
85
subroutine setup_scf_kscell( cell, firm_displ )
87
! ***************** INPUT **********************************************
88
! real*8 cell(3,3) : Unit cell vectors in real space cell(ixyz,ivec)
89
! ***************** OUTPUT *********************************************
90
! logical firm_displ : User-specified displacements (firm)?
92
! The relevant fdf labels are kgrid_cutoff and kgrid_Monkhorst_Pack.
93
! If both are present, kgrid_Monkhorst_Pack has priority. If none is
94
! present, the cutoff default is zero, producing only the gamma point.
95
! Examples of fdf data specifications:
96
! kgrid_cutoff 50. Bohr
97
! %block kgrid_Monkhorst_Pack # Defines kscell and kdispl
98
! 4 0 0 0.50 # (kscell(i,1),i=1,3), kdispl(1)
99
! 0 4 0 0.50 # (kscell(i,2),i=1,3), kdispl(2)
100
! 0 0 4 0.50 # (kscell(i,3),i=1,3), kdispl(3)
101
! %endblock kgrid_Monkhorst_Pack
102
! **********************************************************************
106
use precision, only : dp
107
use parallel, only : Node
108
use m_minvec, only : minvec
117
real(dp), intent(in) :: cell(3,3)
118
logical, intent(out) :: firm_displ
121
integer i, iu, j, factor(3,3), expansion_factor
125
real(dp) scmin(3,3), vmod, cutoff
126
real(dp) ctransf(3,3)
129
real(dp), parameter :: defcut = 0.0_dp
130
integer, dimension(3,3), parameter :: unit_matrix = &
131
reshape ((/1,0,0,0,1,0,0,0,1/), (/3,3/))
135
mp_input = fdf_block('kgrid_Monkhorst_Pack',iu)
137
user_requested_mp = .true.
139
read(iu,*) (kscell(j,i),j=1,3), kdispl(i)
145
cutoff = fdf_physical('kgrid_cutoff',defcut,'Bohr')
146
if (cutoff /= defcut) then
147
!! write(6,"(a,f10.5)") "Kgrid cutoff input: ", cutoff
148
user_requested_cutoff = .true.
151
kdispl(1:3) = 0.0_dp ! Might be changed later
152
firm_displ = .false. ! In future we might add new options
153
! for user-specified displacements
155
! Find equivalent rounded unit-cell
156
call minvec( cell, scmin, ctransf )
161
vmod = sqrt(dot_product(scmin(1:3,j),scmin(1:3,j)))
162
factor(j,j) = int(2.0_dp*cutoff/vmod) + 1
163
expansion_factor = expansion_factor * factor(j,j)
165
! Generate actual supercell skeleton
166
kscell = matmul(ctransf, factor)
167
! Avoid confusing permutations
168
if (expansion_factor == 1) then
175
call MPI_Bcast(kscell(1,1),9,MPI_integer,0,MPI_Comm_World, MPIerror)
176
call MPI_Bcast(kdispl,3,MPI_double_precision,0,MPI_Comm_World, MPIerror)
177
call MPI_Bcast(firm_displ,1,MPI_logical,0,MPI_Comm_World, MPIerror)
178
call MPI_Bcast(user_requested_mp,1,MPI_logical,0, &
179
MPI_Comm_World, MPIerror)
180
call MPI_Bcast(user_requested_cutoff,1,MPI_logical,0, &
181
MPI_Comm_World, MPIerror)
184
end subroutine setup_scf_kscell
186
subroutine siesta_write_k_points()
187
USE siesta_options, only: writek
197
write(6,'(/,a)') 'siesta: k-point coordinates (Bohr**-1) and weights:'
198
write(6,'(a,i4,3f12.6,3x,f12.6)') &
199
('siesta: ', ik, (kpoint(ix,ik),ix=1,3), kweight(ik), &
202
call cmlAddProperty(xf=mainXML, property=kpoint, &
203
dictref='siesta:kpoint')
204
call cmlAddProperty(xf=mainXML, property=kweight, &
205
dictref='siesta:kweight')
208
call iokp( nkpnt, kpoint, kweight )
210
write(6,'(/a,i6)') 'siesta: k-grid: Number of k-points =', nkpnt
211
write(6,'(a,f10.3,a)') 'siesta: k-grid: Cutoff (effective) =', &
212
eff_kgrid_cutoff/Ang, ' Ang'
213
write(6,'(a)') 'siesta: k-grid: Supercell and displacements'
214
write(6,'(a,3i4,3x,f8.3)') 'siesta: k-grid: ', &
215
(kscell(i,1),i=1,3), kdispl(1)
216
write(6,'(a,3i4,3x,f8.3)') 'siesta: k-grid: ', &
217
(kscell(i,2),i=1,3), kdispl(2)
218
write(6,'(a,3i4,3x,f8.3)') 'siesta: k-grid: ', &
219
(kscell(i,3),i=1,3), kdispl(3)
221
call cmlStartPropertyList(mainXML, title='k-points')
222
call cmlAddProperty(xf=mainXML, property=nkpnt,dictref='siesta:nkpnt')
223
call cmlAddProperty(xf=mainXML, property=eff_kgrid_cutoff/Ang, &
224
dictref='siesta:kcutof', units='siestaUnits:angstrom')
225
call cmlEndPropertyList(mainXML)
226
call cmlAddProperty(xf=mainXML, property=kscell,dictref='siesta:kscell')
227
call cmlAddProperty(xf=mainXML, property=kdispl,dictref='siesta:kdispl')
229
end subroutine siesta_write_k_points
231
END MODULE Kpoint_grid