25 use,
intrinsic :: iso_fortran_env
69 real(real64),
allocatable :: point(:, :)
70 real(real64),
allocatable :: point1BZ(:, :)
71 real(real64),
allocatable :: red_point(:, :)
72 real(real64),
allocatable :: weight(:)
73 integer :: nshifts = 1
74 real(real64),
allocatable :: shifts(:,:)
75 integer :: npoints = 0
82 type(lattice_vectors_t) :: latt
84 type(kpoints_grid_t) :: full
85 type(kpoints_grid_t) :: reduced
89 logical :: use_symmetries = .false.
90 logical :: use_time_reversal = .false.
91 integer :: nik_skip = 0
94 integer,
allocatable :: nik_axis(:)
95 integer,
allocatable :: niq_axis(:)
96 integer,
allocatable,
private :: symmetry_ops(:, :)
97 integer,
allocatable,
private :: num_symmetry_ops(:)
100 real(real64),
allocatable,
private :: coord_along_path(:)
102 integer,
allocatable :: downsampling(:)
104 type(symmetries_t),
pointer :: symm => null()
116 integer,
public,
parameter :: &
126 integer,
intent(in) :: dim
127 type(kpoints_grid_t),
intent(out) :: this
128 integer,
intent(in) :: npoints
129 integer,
intent(in) :: nshifts
134 this%npoints = npoints
135 this%nshifts = nshifts
136 safe_allocate(this%red_point(1:dim, 1:npoints))
137 safe_allocate(this%point(1:dim, 1:npoints))
138 safe_allocate(this%point1bz(1:dim,1:npoints))
139 safe_allocate(this%weight(1:npoints))
140 safe_allocate(this%shifts(1:dim,1:nshifts))
147 type(kpoints_grid_t),
intent(inout) :: this
151 safe_deallocate_a(this%red_point)
152 safe_deallocate_a(this%point)
153 safe_deallocate_a(this%point1BZ)
154 safe_deallocate_a(this%weight)
155 safe_deallocate_a(this%shifts)
169 aa%weight(1:bb%npoints) = bb%weight(1:bb%npoints)
170 aa%point(1:bb%dim, 1:bb%npoints) = bb%point(1:bb%dim, 1:bb%npoints)
171 aa%point1BZ(1:bb%dim, 1:bb%npoints) = bb%point1BZ(1:bb%dim, 1:bb%npoints)
172 aa%red_point(1:bb%dim, 1:bb%npoints) = bb%red_point(1:bb%dim, 1:bb%npoints)
173 aa%shifts(1:bb%dim, 1:bb%nshifts) = bb%shifts(1:bb%dim, 1:bb%nshifts)
187 if (.not.
allocated(that%point))
then
192 if (.not.
allocated(this%point))
then
209 this%red_point(1:old_grid%dim, 1:old_grid%npoints)= old_grid%red_point(1:old_grid%dim, 1:old_grid%npoints)
210 this%point(1:old_grid%dim, 1:old_grid%npoints) = old_grid%point(1:old_grid%dim, 1:old_grid%npoints)
211 this%point1BZ(1:old_grid%dim, 1:old_grid%npoints) = old_grid%point1BZ(1:old_grid%dim, 1:old_grid%npoints)
212 this%weight(1:old_grid%npoints) = old_grid%weight(1:old_grid%npoints)
213 this%shifts(1:old_grid%dim, 1:old_grid%nshifts) = old_grid%shifts(1:old_grid%dim, 1:old_grid%nshifts)
216 this%red_point(1:old_grid%dim, old_grid%npoints+1:this%npoints)= that%red_point(1:that%dim, 1:that%npoints)
217 this%point(1:old_grid%dim, old_grid%npoints+1:this%npoints) = that%point(1:that%dim, 1:that%npoints)
218 this%point1BZ(1:old_grid%dim, old_grid%npoints+1:this%npoints) = that%point1BZ(1:that%dim, 1:that%npoints)
219 this%weight(old_grid%npoints+1:this%npoints) = that%weight(1:that%npoints)
228 subroutine kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
232 integer,
intent(in) :: dim
233 integer,
intent(in) :: periodic_dim
236 integer :: ik, idir, is
237 character(len=100) :: str_tmp
238 real(real64) :: weight_sum
239 logical :: default_timereversal, only_gamma
245 safe_allocate(this%nik_axis(1:dim))
246 safe_allocate(this%niq_axis(1:dim))
247 safe_allocate(this%downsampling(1:dim))
251 only_gamma = (periodic_dim == 0)
272 call parse_variable(namespace,
'KPointsUseSymmetries', .false., this%use_symmetries)
292 call parse_variable(namespace,
'KPointsUseTimeReversal', default_timereversal, this%use_time_reversal)
298 this%method = kpoints_gamma
308 call read_mp(gamma_only = .false.)
315 if (this%use_symmetries)
then
316 write(
message(1),
'(a)')
"User-defined k-points are not compatible with KPointsUseSymmetries=yes."
327 if (this%use_symmetries)
then
328 write(
message(1),
'(a)')
"KPointsPath is not compatible with KPointsUseSymmetries=yes."
335 if (this%method == 0)
then
336 write(
message(1),
'(a)')
"Unable to determine the method for defining k-points."
337 write(
message(2),
'(a)')
"Octopus will continue assuming a Monkhorst Pack grid."
340 call read_mp(gamma_only = .false.)
347 write(
message(2),
'(1x,i5,a)') this%reduced%npoints,
' k-points generated from parameters :'
348 write(
message(3),
'(1x,a)')
'---------------------------------------------------'
349 write(
message(4),
'(4x,a)')
'n ='
351 write(str_tmp,
'(i5)') this%nik_axis(idir)
356 do is = 1, this%reduced%nshifts
358 write(
message(2),
'(4x,a,i1,a)')
's', is,
' ='
360 write(str_tmp,
'(f6.2)') this%reduced%shifts(idir,is)
368 write(
message(2),
'(a)')
' index | weight | coordinates |'
371 do ik = 1, this%reduced%npoints
372 write(str_tmp,
'(i8,a,f12.6,a)') ik,
" | ", this%reduced%weight(ik),
" |"
375 write(str_tmp,
'(f12.6)') this%reduced%red_point(idir, ik)
378 write(str_tmp,
'(a)')
" |"
392 logical,
intent(in) :: gamma_only
394 logical :: gamma_only_
395 integer :: ii, is, ncols, nshifts
397 integer,
allocatable :: symm_ops(:, :), num_symm_ops(:)
398 real(real64),
allocatable :: shifts(:,:)
447 gamma_only_ = gamma_only
448 if (.not. gamma_only_)
then
449 gamma_only_ = (
parse_block(namespace,
'KPointsGrid', blk) /= 0)
454 if (.not. gamma_only_)
then
459 safe_allocate(shifts(1:dim,1:nshifts))
462 if (.not. gamma_only_)
then
464 if (ncols /= dim)
then
465 write(
message(1),
'(a,i3,a,i3)')
'KPointsGrid first row has ', ncols,
' columns but should have ', dim
466 if (ncols < dim)
then
469 write(
message(2),
'(a)')
'Continuing, but ignoring the additional values.'
477 if (any(this%nik_axis < 1))
then
478 message(1) =
'Input: KPointsGrid is not valid.'
484 if (ncols /= dim)
then
485 write(
message(1),
'(a,i3,a,i3)')
'KPointsGrid shift has ', ncols,
' columns but must have ', dim
496 if (mod(this%nik_axis(ii), 2) /= 0)
then
526 this%niq_axis(:) = this%nik_axis(:)
527 this%downsampling(:) = 1
529 if (
parse_block(namespace,
'QPointsGrid', blk) == 0)
then
531 if (ncols /= dim)
then
532 write(
message(1),
'(a,i3,a,i3)')
'QPointsGrid first row has ', ncols,
' columns but must have ', dim
539 if (any(this%nik_axis/this%niq_axis /= nint(this%nik_axis/real(this%niq_axis))))
then
540 message(1) =
'Input: QPointsGrid is not compatible with the KPointsGrid.'
544 this%downsampling = this%nik_axis/this%niq_axis
546 if (any(this%downsampling /= 1))
then
557 this%full%shifts = shifts
558 safe_deallocate_a(shifts)
560 call kpoints_grid_generate(dim, this%nik_axis, this%full%nshifts, this%full%shifts, this%full%red_point)
562 do ik = 1, this%full%npoints
566 this%full%weight =
m_one / this%full%npoints
568 if (this%use_symmetries)
then
569 message(1) =
"Checking if the generated full k-point grid is symmetric"
577 if (this%use_symmetries)
then
579 safe_allocate(num_symm_ops(1:this%full%npoints))
581 if (this%use_time_reversal)
then
588 this%reduced%npoints, dim, this%reduced%red_point, this%reduced%weight, symm_ops, num_symm_ops)
591 assert(maxval(num_symm_ops) >= 1)
592 if (this%use_time_reversal)
then
598 assert(sum(num_symm_ops(1:this%reduced%npoints)) == this%full%npoints)
600 do ik = 1, this%reduced%npoints
601 assert(all(symm_ops(ik, 1:num_symm_ops(ik)) <= symm%nops))
604 safe_allocate(this%symmetry_ops(1:this%reduced%npoints, 1:maxval(num_symm_ops)))
605 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
607 this%num_symmetry_ops(1:this%reduced%npoints) = num_symm_ops(1:this%reduced%npoints)
608 this%symmetry_ops(1:this%reduced%npoints, 1:maxval(num_symm_ops)) = &
609 symm_ops(1:this%reduced%npoints, 1:maxval(num_symm_ops))
611 safe_deallocate_a(num_symm_ops)
612 safe_deallocate_a(symm_ops)
616 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
617 safe_allocate(this%symmetry_ops(1:this%reduced%npoints, 1))
618 this%num_symmetry_ops(1:this%reduced%npoints) = 1
619 this%symmetry_ops(1:this%reduced%npoints, 1) = 1
623 do ik = 1, this%reduced%npoints
624 call kpoints_to_absolute(this%latt, this%reduced%red_point(:, ik), this%reduced%point(:, ik))
637 integer :: nshifts, nkpoints, nhighsympoints, nsegments
638 integer :: icol, ik, idir, ncols
639 integer,
allocatable :: resolution(:)
640 real(real64),
allocatable :: highsympoints(:,:)
644 integer,
allocatable :: symm_ops(:, :), num_symm_ops(:)
673 if (
parse_block(namespace,
'KPointsPath', blk) /= 0)
then
674 write(
message(1),
'(a)')
'Internal error while reading KPointsPath.'
681 if (nhighsympoints /= nsegments + 1)
then
682 write(
message(1),
'(a,i3,a,i3)')
'The first row of KPointsPath is not compatible with the number of specified k-points.'
686 safe_allocate(resolution(1:nsegments))
687 do icol = 1, nsegments
691 nkpoints = sum(resolution) + 1
693 safe_allocate(highsympoints(1:dim, 1:nhighsympoints))
694 do ik = 1, nhighsympoints
697 if (ncols /= dim)
then
698 write(
message(1),
'(a,i8,a,i3,a,i3)')
'KPointsPath row ', ik,
' has ', ncols,
' columns but must have ', dim
714 safe_allocate(this%coord_along_path(1:nkpoints))
716 call kpoints_path_generate(dim, this%latt, nkpoints, nsegments, resolution, highsympoints, path_kpoints_grid%point, &
717 this%coord_along_path)
719 safe_deallocate_a(resolution)
720 safe_deallocate_a(highsympoints)
725 path_kpoints_grid%weight =
m_one/path_kpoints_grid%npoints
727 path_kpoints_grid%weight =
m_zero
728 this%nik_skip = this%nik_skip + path_kpoints_grid%npoints
733 do ik = 1, path_kpoints_grid%npoints
734 call kpoints_to_reduced(this%latt, path_kpoints_grid%point(:, ik), path_kpoints_grid%red_point(:, ik))
741 if (
allocated(this%symmetry_ops))
then
742 npoints = this%reduced%npoints
743 safe_allocate(num_symm_ops(1:npoints))
744 safe_allocate(symm_ops(1:npoints, 1:maxval(this%num_symmetry_ops)))
746 num_symm_ops(1:npoints) = this%num_symmetry_ops(1:npoints)
747 symm_ops(1:npoints, 1:maxval(num_symm_ops)) = &
748 this%symmetry_ops(1:npoints, 1:maxval(num_symm_ops))
754 if (
allocated(this%symmetry_ops))
then
755 safe_deallocate_a(this%num_symmetry_ops)
756 safe_deallocate_a(this%symmetry_ops)
757 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
758 safe_allocate(this%symmetry_ops(1:this%reduced%npoints,1:maxval(num_symm_ops)))
760 this%num_symmetry_ops(1:npoints) = num_symm_ops(1:npoints)
761 this%symmetry_ops(1:npoints, 1:maxval(num_symm_ops)) = &
762 symm_ops(1:npoints, 1:maxval(num_symm_ops))
763 this%num_symmetry_ops(npoints+1:this%reduced%npoints) = 1
764 this%symmetry_ops(npoints+1:this%reduced%npoints, 1) = 1
766 safe_deallocate_a(num_symm_ops)
767 safe_deallocate_a(symm_ops)
768 else if(this%use_symmetries)
then
769 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
770 safe_allocate(this%symmetry_ops(1:this%reduced%npoints,1:1))
772 this%num_symmetry_ops(1:this%reduced%npoints) = 1
773 this%symmetry_ops(1:this%reduced%npoints, 1) = 1
824 if (
parse_block(namespace,
'KPoints', blk) /= 0)
then
825 if (
parse_block(namespace,
'KPointsReduced', blk) == 0)
then
829 write(
message(1),
'(a)')
'Internal error loading user-defined k-point list.'
839 user_kpoints_grid%red_point =
m_zero
840 user_kpoints_grid%point =
m_zero
841 user_kpoints_grid%weight =
m_zero
842 user_kpoints_grid%shifts =
m_zero
846 do ik = 1, user_kpoints_grid%npoints
852 call kpoints_to_absolute(this%latt, user_kpoints_grid%red_point(:, ik), user_kpoints_grid%point(:, ik))
855 do ik = 1, user_kpoints_grid%npoints
861 call kpoints_to_reduced(this%latt, user_kpoints_grid%point(:, ik), user_kpoints_grid%red_point(:, ik))
867 if (any(user_kpoints_grid%weight <
m_epsilon))
then
869 message(1) =
"Found k-points with zero weight. They are excluded from density calculation"
875 do ik=1,user_kpoints_grid%npoints
878 if (ik < user_kpoints_grid%npoints)
then
879 if (user_kpoints_grid%weight(ik+1) >
m_epsilon)
then
880 message(1) =
"K-points with zero weight must follow all regular k-points in a block"
884 this%nik_skip = this%nik_skip + 1
886 user_kpoints_grid%weight(ik) =
m_zero
891 weight_sum = sum(user_kpoints_grid%weight)
893 message(1) =
"k-point weights must sum to a positive number."
896 user_kpoints_grid%weight = user_kpoints_grid%weight / weight_sum
907 write(
message(1),
'(a,i4,a)')
'Input: ', user_kpoints_grid%npoints,
' k-points were read from the input file'
926 safe_deallocate_a(this%nik_axis)
927 safe_deallocate_a(this%niq_axis)
928 safe_deallocate_a(this%downsampling)
929 safe_deallocate_a(this%symmetry_ops)
930 safe_deallocate_a(this%num_symmetry_ops)
931 safe_deallocate_a(this%coord_along_path)
939 real(real64),
intent(in) :: kin(:)
940 real(real64),
intent(out) :: kout(:)
944 kout = matmul(latt%klattice, kin)
952 real(real64),
intent(in) :: kin(:)
953 real(real64),
intent(out) :: kout(:)
957 kout = matmul(kin, latt%rlattice) / (
m_two*
m_pi)
972 kout%method = kin%method
977 kout%use_symmetries = kin%use_symmetries
978 kout%use_time_reversal = kin%use_time_reversal
980 safe_allocate(kout%nik_axis(1:kin%full%dim))
981 safe_allocate(kout%niq_axis(1:kin%full%dim))
982 safe_allocate(kout%downsampling(1:kin%full%dim))
983 kout%nik_axis = kin%nik_axis
984 kout%niq_axis = kin%niq_axis
985 kout%downsampling = kin%downsampling
987 if (
allocated(kin%coord_along_path))
then
988 safe_allocate(kout%coord_along_path(1:kin%full%npoints))
989 kout%coord_along_path(1:kin%full%npoints) = kin%coord_along_path(1:kin%full%npoints)
992 if (
allocated(kin%symmetry_ops))
then
993 safe_allocate(kout%num_symmetry_ops(1:kin%reduced%npoints))
994 safe_allocate(kout%symmetry_ops(1:kin%reduced%npoints, 1:maxval(kin%num_symmetry_ops)))
995 kout%num_symmetry_ops(1:kin%reduced%npoints) = kin%num_symmetry_ops(1:kin%reduced%npoints)
996 kout%symmetry_ops(1:kin%reduced%npoints, 1:maxval(kin%num_symmetry_ops)) = &
997 kin%symmetry_ops(1:kin%reduced%npoints, 1:maxval(kin%num_symmetry_ops))
1005 integer pure function kpoints_number(this) result(number)
1008 number = this%reduced%npoints
1016 integer,
intent(in) :: ik
1017 logical,
optional,
intent(in) :: absolute_coordinates
1018 real(real64) :: point(1:this%full%dim)
1020 if (optional_default(absolute_coordinates, .
true.))
then
1021 point(1:this%full%dim) = this%reduced%point(1:this%full%dim, ik)
1023 point(1:this%full%dim) = this%reduced%red_point(1:this%full%dim, ik)
1032 integer,
intent(in) :: ik
1034 weight = this%reduced%weight(ik)
1048 integer,
intent(in) :: dim
1049 integer,
intent(in) :: naxis(1:dim)
1050 integer,
intent(in) :: nshifts
1051 real(real64),
intent(in) :: shift(:,:)
1052 real(real64),
intent(out) :: kpoints(:, :)
1053 integer,
optional,
intent(out) :: lk123(:,:)
1056 real(real64) :: dx(dim), maxcoord
1057 integer :: ii, jj, divisor, ik, idir, npoints, is, ix(dim)
1058 integer,
allocatable :: lk123_(:,:), idx(:)
1059 real(real64),
allocatable :: nrm(:), shell(:), coords(:, :)
1063 dx = m_one/(m_two*naxis)
1065 npoints = product(naxis)
1067 if (
present(lk123))
then
1068 safe_allocate(lk123_(1:npoints*nshifts,1:dim))
1069 safe_allocate(idx(1:npoints*nshifts))
1073 do ii = 0, npoints - 1
1074 ik = npoints*is - ii
1079 divisor = divisor / naxis(idir)
1080 ix(idir) = jj / divisor + 1
1081 jj = mod(jj, divisor)
1083 kpoints(idir, ik) = (m_two*ix(idir) - m_one*naxis(idir) + m_two*shift(idir,is))*dx(idir)
1087 if (abs(kpoints(idir, ik) - m_half) > 1e-14_real64)
then
1088 kpoints(idir, ik) = mod(kpoints(idir, ik) + m_half, m_one) - m_half
1092 if (
present(lk123))
then
1101 safe_allocate(nrm(1:npoints*nshifts))
1102 safe_allocate(shell(1:npoints*nshifts))
1103 safe_allocate(coords(1:dim, 1:npoints*nshifts))
1105 do ik = 1, npoints*nshifts
1106 shell(ik) = sum((kpoints(:, ik)/dx)**2)
1108 coords(idir, ik) = kpoints(idir, ik)
1109 if (coords(idir, ik) < m_zero) coords(idir, ik) = coords(idir, ik) + m_one
1110 coords(idir, ik) = coords(idir, ik)/(dx(idir)*m_two)
1118 do ik = 1, npoints*nshifts
1119 nrm(ik) = nrm(ik) + coords(idir, ik)*maxcoord
1121 maxcoord = maxcoord*max(m_one, maxval(coords(idir, 1:npoints*nshifts)))
1124 do ik = 1, npoints*nshifts
1125 nrm(ik) = nrm(ik) + shell(ik)*maxcoord
1128 if (
present(lk123))
then
1130 do ik = 1, npoints*nshifts
1131 lk123(ik,:) = lk123_(idx(ik),:)
1133 safe_deallocate_a(lk123_)
1134 safe_deallocate_a(idx)
1137 call sort(nrm, kpoints)
1140 safe_deallocate_a(nrm)
1141 safe_deallocate_a(shell)
1142 safe_deallocate_a(coords)
1149 subroutine kpoints_path_generate(dim, latt, nkpoints, nsegments, resolution, highsympoints, kpoints, coord)
1150 integer,
intent(in) :: dim
1151 type(lattice_vectors_t),
intent(in) :: latt
1152 integer,
intent(in) :: nkpoints
1153 integer,
intent(in) :: nsegments
1154 integer,
intent(in) :: resolution(:)
1155 real(real64),
intent(in) :: highsympoints(1:dim,1:nsegments)
1156 real(real64),
intent(out) :: kpoints(1:dim, 1:nkpoints)
1157 real(real64),
intent(out) :: coord(1:nkpoints)
1159 integer :: is, ik, kpt_ind
1160 real(real64) :: length, total_length, accumulated_length
1161 real(real64) :: kpt1(dim), kpt2(dim), vec(dim)
1165 total_length = m_zero
1167 do is = 1, nsegments
1174 if (resolution(is) > 0) total_length = total_length + length
1177 accumulated_length = m_zero
1180 do is = 1, nsegments
1189 do ik = 1, resolution(is)
1190 kpt_ind = kpt_ind +1
1191 coord(kpt_ind) = accumulated_length + (ik-1)*length/resolution(is)
1192 kpoints(:, kpt_ind) = kpt1 + (ik-1)*length/resolution(is)*vec
1194 if (resolution(is) > 0) accumulated_length = accumulated_length + length
1197 kpt_ind = kpt_ind +1
1199 coord(kpt_ind) = accumulated_length
1200 kpoints(:, kpt_ind) = kpt1
1203 coord = coord/total_length
1210 subroutine kpoints_grid_reduce(symm, time_reversal, nkpoints, dim, kpoints, weights, symm_ops, num_symm_ops)
1211 type(symmetries_t),
intent(in) :: symm
1212 logical,
intent(in) :: time_reversal
1213 integer,
intent(inout) :: nkpoints
1214 integer,
intent(in) :: dim
1215 real(real64),
intent(inout) :: kpoints(1:dim,1:nkpoints)
1216 real(real64),
intent(out) :: weights(1:nkpoints)
1217 integer,
intent(out) :: symm_ops(:, :)
1218 integer,
intent(out) :: num_symm_ops(:)
1221 real(real64),
allocatable :: reduced(:, :)
1223 integer ik, iop, ik2
1224 real(real64) :: tran(dim), diff(dim)
1225 real(real64),
allocatable :: kweight(:)
1227 real(real64) :: prec
1234 prec = min(symprec, m_one/(nkpoints*100))
1241 safe_allocate(kweight(1:nkpoints))
1242 safe_allocate(reduced(1:dim, 1:nkpoints))
1244 kweight = m_one / nkpoints
1249 symm_ops(:, 1) = symmetries_identity_index(symm)
1252 if (kweight(ik) < prec) cycle
1256 nreduced = nreduced + 1
1257 reduced(:, nreduced) = kpoints(:, ik)
1260 if (maxval(abs(kpoints(:, ik))) < m_epsilon) cycle
1262 if (ik == nkpoints) cycle
1265 do iop = 1, symmetries_number(symm)
1266 if (iop == symmetries_identity_index(symm) .and. &
1267 .not. time_reversal) cycle
1269 call symmetries_apply_kpoint_red(symm, iop, reduced(1:dim, nreduced), tran)
1271 tran = tran - anint(tran + m_half*prec)
1274 do ik2 = ik + 1, nkpoints
1275 if (kweight(ik2) < prec) cycle
1277 if (.not. iop == symmetries_identity_index(symm))
then
1278 diff = tran - kpoints(:, ik2)
1279 diff = diff - anint(diff)
1282 if (sum(abs(diff)) < prec)
then
1283 kweight(ik) = kweight(ik) + kweight(ik2)
1284 kweight(ik2) = m_zero
1285 weights(nreduced) = kweight(ik)
1286 num_symm_ops(nreduced) = num_symm_ops(nreduced) + 1
1287 symm_ops(nreduced, num_symm_ops(nreduced)) = iop
1292 if (time_reversal)
then
1293 diff = tran + kpoints(:, ik2)
1294 diff = diff - anint(diff)
1297 if (sum(abs(diff)) < prec)
then
1298 kweight(ik) = kweight(ik) + kweight(ik2)
1299 kweight(ik2) = m_zero
1300 weights(nreduced) = kweight(ik)
1301 num_symm_ops(nreduced) = num_symm_ops(nreduced) + 1
1303 symm_ops(nreduced, num_symm_ops(nreduced)) = -iop
1310 assert(sum(weights(1:nreduced)) - m_one < prec)
1314 kpoints(:, ik) = reduced(:, ik)
1317 safe_deallocate_a(kweight)
1318 safe_deallocate_a(reduced)
1327 type(lattice_vectors_t),
intent(in) :: latt
1329 integer :: ig1, ig2, ik
1330 real(real64) :: Gvec(grid%dim, 3**grid%dim), Gvec_cart(grid%dim, 3**grid%dim)
1331 real(real64) :: vec(grid%dim), kpt(grid%dim)
1332 real(real64) :: d, dmin
1337 do ig1 = 1, grid%dim
1338 do ig2 = 1, 3**grid%dim
1339 gvec(ig1, ig2) = modulo((ig2 - 1)/grid%dim**(grid%dim - ig1), 3) - 1
1343 do ig1 = 1, 3**grid%dim
1347 do ik = 1, grid%npoints
1350 do ig1 = 1, 3**grid%dim
1351 vec = gvec_cart(:,ig1) - grid%point(:,ik)
1352 d = real(sum(vec**2), real32)
1359 kpt = grid%red_point(:,ik) - gvec(:,ig2)
1370 integer,
optional,
intent(in) :: iunit
1371 type(namespace_t),
optional,
intent(in) :: namespace
1372 logical,
optional,
intent(in) :: absolute_coordinates
1375 character(len=100) :: str_tmp
1380 call messages_print_with_emphasis(
'Brillouin zone sampling', iunit, namespace)
1384 call messages_write(
'Dimensions of the k-point grid =')
1385 do idir = 1, this%full%dim
1386 call messages_write(this%nik_axis(idir), fmt =
'(i3,1x)')
1388 call messages_new_line()
1390 if (.not. all(this%downsampling(1:this%full%dim) == 1))
then
1391 call messages_write(
'Dimensions of the q-point grid =')
1392 do idir = 1, this%full%dim
1393 call messages_write(this%niq_axis(idir), fmt =
'(i3,1x)')
1395 call messages_new_line()
1398 call messages_write(
'Total number of k-points =')
1399 call messages_write(this%full%npoints)
1400 call messages_new_line()
1402 call messages_write(
'Number of symmetry-reduced k-points =')
1403 call messages_write(this%reduced%npoints)
1405 call messages_info(iunit=iunit, namespace=namespace)
1409 call messages_write(
'Total number of k-points =')
1410 call messages_write(this%full%npoints)
1411 call messages_new_line()
1412 call messages_info(iunit=iunit, namespace=namespace)
1416 call messages_new_line()
1417 call messages_write(
'List of k-points:')
1418 call messages_info(iunit=iunit, namespace=namespace)
1420 write(message(1),
'(6x,a)')
'ik'
1421 do idir = 1, this%full%dim
1422 index = index2axis(idir)
1423 write(str_tmp,
'(9x,2a)')
'k_', index
1424 message(1) = trim(message(1)) // trim(str_tmp)
1426 write(str_tmp,
'(6x,a)')
'Weight'
1427 message(1) = trim(message(1)) // trim(str_tmp)
1428 message(2) =
'---------------------------------------------------------'
1429 call messages_info(2, iunit)
1432 write(message(1),
'(i8,1x)') ik
1433 do idir = 1, this%full%dim
1434 if (optional_default(absolute_coordinates, .false.))
then
1435 write(str_tmp,
'(f12.4)') this%reduced%point(idir, ik)
1437 write(str_tmp,
'(f12.4)') this%reduced%red_point(idir, ik)
1439 message(1) = trim(message(1)) // trim(str_tmp)
1442 message(1) = trim(message(1)) // trim(str_tmp)
1443 call messages_info(1, iunit, namespace=namespace)
1446 call messages_info(iunit=iunit, namespace=namespace)
1448 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
1457 integer,
intent(in) :: ik
1465 integer pure function kpoints_get_num_symmetry_ops(this, ik) result(num)
1467 integer,
intent(in) :: ik
1469 if (this%use_symmetries)
then
1470 num = this%num_symmetry_ops(ik)
1478 integer pure function kpoints_get_symmetry_ops(this, ik, index) result(iop)
1480 integer,
intent(in) :: ik
1481 integer,
intent(in) :: index
1483 if (this%use_symmetries)
then
1484 iop = this%symmetry_ops(ik, index)
1492 logical pure function kpoints_is_valid_symmetry(this, ik, index) result(valid)
1494 integer,
intent(in) :: ik
1495 integer,
intent(in) :: index
1500 if (this%use_symmetries)
then
1501 do iop = 1, this%num_symmetry_ops(ik)
1502 if (this%symmetry_ops(ik, iop) == index)
then
1518 integer :: denom, nik, nik_skip
1522 nik = this%full%npoints
1523 nik_skip = this%nik_skip
1530 do denom = 1, 100000
1531 if (all(abs(int(this%full%weight(1:nik-nik_skip)*denom + 10_real64*m_epsilon) - &
1532 this%full%weight(1:nik-nik_skip)*denom) < 100_real64*m_epsilon))
then
1543 logical pure function kpoints_have_zero_weight_path(this) result(have_zerow)
1546 if (this%nik_skip > 0)
then
1549 have_zerow = .false.
1564 integer,
intent(in) :: ind
1566 coord = this%coord_along_path(ind)
1574 type(symmetries_t),
intent(in) :: symm
1575 integer,
intent(in) :: dim
1576 logical,
intent(in) :: time_reversal
1577 type(namespace_t),
intent(in) :: namespace
1579 integer,
allocatable :: kmap(:)
1580 real(real64) :: kpt(dim), diff(dim)
1581 integer :: nk, ik, ik2, iop
1582 type(distributed_t) :: kpt_dist
1589 call distributed_nullify(kpt_dist, nk)
1591 if (mpi_world%comm /= mpi_comm_undefined)
then
1592 call distributed_init(kpt_dist, nk, mpi_comm_world,
"kpt_check")
1597 safe_allocate(kmap(kpt_dist%start:kpt_dist%end))
1599 do iop = 1, symmetries_number(symm)
1600 if (iop == symmetries_identity_index(symm) .and. &
1601 .not. time_reversal) cycle
1603 do ik = kpt_dist%start, kpt_dist%end
1607 do ik = kpt_dist%start, kpt_dist%end
1609 call symmetries_apply_kpoint_red(symm, iop, grid%red_point(:, ik), kpt)
1611 kpt = kpt - anint(kpt + m_half*symprec)
1616 if (iop /= symmetries_identity_index(symm))
then
1617 diff = kpt - grid%red_point(:, ik2)
1618 diff = diff - anint(diff)
1620 if (sum(abs(diff)) < symprec)
then
1626 if (time_reversal)
then
1627 diff = kpt + grid%red_point(:, ik2)
1628 diff = diff - anint(diff)
1630 if (sum(abs(diff)) < symprec)
then
1638 if (kmap(ik) == ik)
then
1639 write(message(1),
'(a,i5,a2,3(f7.3,a2),a)')
"The reduced k-point ", ik,
" (", &
1640 grid%red_point(1, ik),
", ", grid%red_point(2, ik),
", ", grid%red_point(3, ik), &
1641 ") ",
"has no symmetric in the k-point grid for the following symmetry"
1642 write(message(2),
'(i5,1x,a,2x,3(3i4,2x))') iop,
':', transpose(symm_op_rotation_matrix_red(symm%ops(iop)))
1643 message(3) =
"Change your k-point grid or use KPointsUseSymmetries=no."
1644 call messages_fatal(3, namespace=namespace)
1649 safe_deallocate_a(kmap)
1651 call distributed_end(kpt_dist)
1659 integer,
intent(in) :: ik
1660 integer,
intent(in) :: iq
1662 integer :: dim, idim
1663 real(real64) :: diff(kpt%reduced%dim), red(kpt%reduced%dim)
1668 dim = kpt%reduced%dim
1671 if (all(kpt%downsampling == 1))
then
1678 diff = kpt%reduced%red_point(:, ik) - kpt%reduced%red_point(:, iq)
1681 diff(idim) = diff(idim) - anint(diff(idim) + m_half*symprec)
1682 red(idim) = diff(idim)*kpt%nik_axis(idim)/real(kpt%downsampling(idim))
1683 if (abs(red(idim) - anint(red(idim))) > m_epsilon)
then
1684 compatible = .false.
1698 if (
allocated(this%coord_along_path))
then
1699 npath =
SIZE(this%coord_along_path)
1718 type(lattice_vectors_t),
intent(in) :: new_latt
1724 this%latt = new_latt
1726 do ik = 1, this%reduced%npoints
1727 call kpoints_to_absolute(this%latt, this%reduced%red_point(:, ik), this%reduced%point(:, ik))
1730 do ik = 1, this%full%npoints
subroutine read_mp(gamma_only)
subroutine read_path()
Read the k-points path information and generate the k-points list.
subroutine read_user_kpoints()
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
integer pure function kpoints_get_kpoint_method(this)
subroutine, public kpoints_grid_generate(dim, naxis, nshifts, shift, kpoints, lk123)
Generates the k-points grid.
subroutine, public kpoints_end(this)
logical pure function kpoints_have_zero_weight_path(this)
subroutine kpoints_check_symmetries(grid, symm, dim, time_reversal, namespace)
integer, parameter, public kpoints_monkh_pack
integer pure function, public kpoints_get_num_symmetry_ops(this, ik)
subroutine kpoints_grid_addto(this, that)
real(real64) pure function, public kpoints_get_path_coord(this, ind)
logical pure function, public kpoints_is_valid_symmetry(this, ik, index)
integer function, public kpoints_kweight_denominator(this)
subroutine, public kpoints_copy(kin, kout)
subroutine, public kpoints_path_generate(dim, latt, nkpoints, nsegments, resolution, highsympoints, kpoints, coord)
Generate the k-point along a path.
integer, parameter, public kpoints_path
pure real(real64) function, dimension(1:this%full%dim) kpoints_get_point(this, ik, absolute_coordinates)
integer function, public kpoints_nkpt_in_path(this)
logical function, public kpoints_is_compatible_downsampling(kpt, ik, iq)
real(real64) pure function kpoints_get_weight(this, ik)
integer pure function, public kpoints_get_symmetry_ops(this, ik, index)
subroutine kpoints_grid_copy(bb, aa)
subroutine, public kpoints_fold_to_1bz(grid, latt)
subroutine, public kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
subroutine, public kpoints_grid_end(this)
logical pure function, public kpoints_point_is_gamma(this, ik)
integer, parameter, public kpoints_user
subroutine kpoints_grid_reduce(symm, time_reversal, nkpoints, dim, kpoints, weights, symm_ops, num_symm_ops)
logical function kpoints_gamma_only(this)
subroutine, public kpoints_lattice_vectors_update(this, new_latt)
subroutine, public kpoints_to_reduced(latt, kin, kout)
integer pure function, public kpoints_number(this)
subroutine, public kpoints_to_absolute(latt, kin, kout)
subroutine, public kpoints_grid_init(dim, this, npoints, nshifts)
subroutine kpoints_write_info(this, iunit, namespace, absolute_coordinates)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_experimental(name, namespace)
logical function, public parse_is_defined(namespace, name)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
This module is intended to contain "only mathematical" functions and procedures.
integer pure function, public symmetries_number(this)
logical pure function, public symmetries_have_break_dir(this)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
This module is intended to contain simple general-purpose utility functions and procedures.