Octopus
kpoints.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module kpoints_oct_m
22 use debug_oct_m
24 use global_oct_m
25 use, intrinsic :: iso_fortran_env
28 use mpi_oct_m
30 use parser_oct_m
32 use sort_oct_m
35 use unit_oct_m
37 use utils_oct_m
38
39 implicit none
40
41 private
42
43 public :: &
45 kpoints_t, &
66
68 ! Components are public by default
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
76 integer :: dim = 0
77 end type kpoints_grid_t
78
79 type kpoints_t
80 ! Components are public by default
81
82 type(lattice_vectors_t) :: latt
83
84 type(kpoints_grid_t) :: full
85 type(kpoints_grid_t) :: reduced
86
87 integer :: method = 0
88
89 logical :: use_symmetries = .false.
90 logical :: use_time_reversal = .false.
91 integer :: nik_skip = 0
92
94 integer, allocatable :: nik_axis(:)
95 integer, allocatable :: niq_axis(:)
96 integer, allocatable, private :: symmetry_ops(:, :)
97 integer, allocatable, private :: num_symmetry_ops(:)
98
100 real(real64), allocatable, private :: coord_along_path(:)
101
102 integer, allocatable :: downsampling(:)
103
104 type(symmetries_t), pointer :: symm => null()
105
106 contains
107 procedure :: have_zero_weight_path => kpoints_have_zero_weight_path
108 procedure :: get_kpoint_method => kpoints_get_kpoint_method
109 procedure :: gamma_only => kpoints_gamma_only
110 procedure :: get_weight => kpoints_get_weight
111 procedure :: get_point => kpoints_get_point
112 procedure :: write_info => kpoints_write_info
113 procedure :: nkpt_in_path => kpoints_nkpt_in_path
114 end type kpoints_t
115
116 integer, public, parameter :: &
117 KPOINTS_GAMMA = 2, &
118 kpoints_monkh_pack = 4, &
119 kpoints_user = 8, &
120 kpoints_path = 16
121
122contains
123
124 ! ---------------------------------------------------------
125 subroutine kpoints_grid_init(dim, this, npoints, nshifts)
126 integer, intent(in) :: dim
127 type(kpoints_grid_t), intent(out) :: this
128 integer, intent(in) :: npoints
129 integer, intent(in) :: nshifts
130
131 push_sub(kpoints_grid_init)
132
133 this%dim = dim
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))
141
142 pop_sub(kpoints_grid_init)
143 end subroutine kpoints_grid_init
144
145 ! ---------------------------------------------------------
146 subroutine kpoints_grid_end(this)
147 type(kpoints_grid_t), intent(inout) :: this
148
149 push_sub(kpoints_grid_end)
150
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)
156
157 pop_sub(kpoints_grid_end)
158 end subroutine kpoints_grid_end
159
160 ! ---------------------------------------------------------
161 subroutine kpoints_grid_copy(bb, aa)
162 type(kpoints_grid_t), intent(in) :: bb
163 type(kpoints_grid_t), intent(inout) :: aa
168 call kpoints_grid_init(bb%dim, aa, bb%npoints, bb%nshifts)
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)
174
176 end subroutine kpoints_grid_copy
178 ! ---------------------------------------------------------
179 subroutine kpoints_grid_addto(this, that)
180 type(kpoints_grid_t), intent(inout) :: this
181 type(kpoints_grid_t), intent(in) :: that
183 type(kpoints_grid_t) :: old_grid
185 push_sub(kpoints_grid_addto)
186
187 if (.not. allocated(that%point)) then
189 return
190 end if
191
192 if (.not. allocated(this%point)) then
193 call kpoints_grid_copy(that, this)
194 pop_sub(kpoints_grid_addto)
195 return
196 end if
198 call kpoints_grid_copy(this, old_grid)
199
201 call kpoints_grid_init(old_grid%dim, this, that%npoints + old_grid%npoints, old_grid%nshifts)
203 this%red_point = m_zero
204 this%point = m_zero
205 this%weight = m_zero
206 this%shifts = m_zero
207
208 ! Fill the the result with values form this
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)
214
215 ! Fill the result with that
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)
220
221 call kpoints_grid_end(old_grid)
222
223
224 pop_sub(kpoints_grid_addto)
225 end subroutine kpoints_grid_addto
226
227 ! ---------------------------------------------------------
228 subroutine kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
229 type(kpoints_t), intent(out) :: this
230 type(namespace_t), intent(in) :: namespace
231 type(symmetries_t), target, intent(in) :: symm
232 integer, intent(in) :: dim
233 integer, intent(in) :: periodic_dim
234 type(lattice_vectors_t), intent(in) :: latt
235
236 integer :: ik, idir, is
237 character(len=100) :: str_tmp
238 real(real64) :: weight_sum
239 logical :: default_timereversal, only_gamma
240
241 push_sub(kpoints_init)
242
243 this%latt = latt
244
245 safe_allocate(this%nik_axis(1:dim))
246 safe_allocate(this%niq_axis(1:dim))
247 safe_allocate(this%downsampling(1:dim))
248
249 this%nik_axis = 1
250
251 only_gamma = (periodic_dim == 0)
252
253 this%symm => symm
255 !%Variable KPointsUseSymmetries
256 !%Type logical
257 !%Default no
258 !%Section Mesh::KPoints
259 !%Description
260 !% This variable defines whether symmetries are taken into account
261 !% or not for the choice of <i>k</i>-points. If it is set to no, the <i>k</i>-point
262 !% sampling will range over the full Brillouin zone.
263 !%
264 !% When a perturbation is applied to the system, the full
265 !% symmetries of the system cannot be used. In this case you must
266 !% not use symmetries or use the <tt>SymmetryBreakDir</tt> to tell
267 !% Octopus the direction of the perturbation (for the moment this
268 !% has to be done by hand by the user, in the future it will be
269 !% automatic).
270 !%
271 !%End
272 call parse_variable(namespace, 'KPointsUseSymmetries', .false., this%use_symmetries)
273
274 !%Variable KPointsUseTimeReversal
275 !%Type logical
276 !%Section Mesh::KPoints
277 !%Description
278 !% If symmetries are used to reduce the number of <i>k</i>-points,
279 !% this variable defines whether time-reversal symmetry is taken
280 !% into account or not. If it is set to no, the <i>k</i>-point
281 !% sampling will not be reduced according to time-reversal
282 !% symmetry.
283 !%
284 !% The default is yes, unless symmetries are broken in one
285 !% direction by the <tt>SymmetryBreakDir</tt> block.
286 !%
287 !% Warning: For time propagation runs with an external field,
288 !% time-reversal symmetry should not be used.
289 !%
290 !%End
291 default_timereversal = this%use_symmetries .and. .not. symmetries_have_break_dir(symm)
292 call parse_variable(namespace, 'KPointsUseTimeReversal', default_timereversal, this%use_time_reversal)
293
294 !We determine the method used to define k-point
295 this%method = 0
296
297 if (only_gamma) then
298 this%method = kpoints_gamma
299 call read_mp(gamma_only = .true.)
300 pop_sub(kpoints_init)
301 return
302 end if
303
304 !Monkhorst Pack grid
305 if (parse_is_defined(namespace, 'KPointsGrid')) then
306 this%method = this%method + kpoints_monkh_pack
307
308 call read_mp(gamma_only = .false.)
309 end if
310
311 !User-defined kpoints
312 if (parse_is_defined(namespace, 'KPointsReduced').or. parse_is_defined(namespace, 'KPoints')) then
313 this%method = this%method + kpoints_user
314
315 if (this%use_symmetries) then
316 write(message(1), '(a)') "User-defined k-points are not compatible with KPointsUseSymmetries=yes."
317 call messages_warning(1, namespace=namespace)
318 end if
319
320 call read_user_kpoints()
321 end if
322
323 !User-defined k-points path
324 if (parse_is_defined(namespace, 'KPointsPath')) then
325 this%method = this%method + kpoints_path
326
327 if (this%use_symmetries) then
328 write(message(1), '(a)') "KPointsPath is not compatible with KPointsUseSymmetries=yes."
329 call messages_warning(1, namespace=namespace)
330 end if
331 call read_path()
332 end if
333
334
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."
338 call messages_warning(2, namespace=namespace)
339 this%method = kpoints_monkh_pack
340 call read_mp(gamma_only = .false.)
341 end if
342
343 !Printing the k-point list
344 if (bitand(this%method, kpoints_monkh_pack) /= 0) then
345
346 write(message(1),'(a)') ' '
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 ='
350 do idir = 1, dim
351 write(str_tmp,'(i5)') this%nik_axis(idir)
352 message(4) = trim(message(4)) // trim(str_tmp)
353 end do
354 call messages_info(4, namespace=namespace)
355
356 do is = 1, this%reduced%nshifts
357 write(message(1),'(a)') ' '
358 write(message(2),'(4x,a,i1,a)') 's', is, ' ='
359 do idir = 1, dim
360 write(str_tmp,'(f6.2)') this%reduced%shifts(idir,is)
361 message(2) = trim(message(2)) // trim(str_tmp)
362 end do
363 call messages_info(2, namespace=namespace)
364 end do
365 end if
366
367 write(message(1),'(a)') ' '
368 write(message(2),'(a)') ' index | weight | coordinates |'
369 call messages_info(2, namespace=namespace)
370
371 do ik = 1, this%reduced%npoints
372 write(str_tmp,'(i8,a,f12.6,a)') ik, " | ", this%reduced%weight(ik), " |"
373 message(1) = str_tmp
374 do idir = 1, dim
375 write(str_tmp,'(f12.6)') this%reduced%red_point(idir, ik)
376 message(1) = trim(message(1)) // trim(str_tmp)
377 end do
378 write(str_tmp,'(a)') " |"
379 message(1) = trim(message(1)) // trim(str_tmp)
380 call messages_info(1, namespace=namespace)
381 end do
382
383 write(message(1),'(a)') ' '
384 call messages_info(1, namespace=namespace)
385
386 pop_sub(kpoints_init)
387
388 contains
389
390 ! ---------------------------------------------------------
391 subroutine read_mp(gamma_only)
392 logical, intent(in) :: gamma_only
393
394 logical :: gamma_only_
395 integer :: ii, is, ncols, nshifts
396 type(block_t) :: blk
397 integer, allocatable :: symm_ops(:, :), num_symm_ops(:)
398 real(real64), allocatable :: shifts(:, :)
399
400
401 push_sub(kpoints_init.read_mp)
402
403 call messages_obsolete_variable(namespace, 'KPointsMonkhorstPack', 'KPointsGrid')
404
405 !%Variable KPointsGrid
406 !%Type block
407 !%Default <math>\Gamma</math>-point only
408 !%Section Mesh::KPoints
409 !%Description
410 !% When this block is given (and the <tt>KPoints</tt> block is not present),
411 !% <i>k</i>-points are distributed in a uniform grid, according to a modified
412 !% version of the Monkhorst-Pack scheme. For the original MP scheme, see
413 !% James D. Pack and Hendrik J. Monkhorst,
414 !% <i>Phys. Rev. B</i> <b>13</b>, 5188 (1976) and <i>Phys. Rev. B</i> <b>16</b>, 1748 (1977).
415 !%
416 !% The number of columns should be equal to <tt>Dimensions</tt>,
417 !% but the grid and shift numbers should be 1 and zero in finite directions.
418 !%
419 !% The first row of the block is a set of integers defining
420 !% the number of <i>k</i>-points to be used along each direction
421 !% in reciprocal space. The numbers refer to the whole Brillouin
422 !% zone, and the actual number of <i>k</i>-points is usually
423 !% reduced exploiting the symmetries of the system. By default
424 !% the grid will always include the <math>\Gamma</math>-point.
425 !%
426 !% Optional rows can be added to specify multiple shifts in the <i>k</i>-points (between 0.0 and 1.0),
427 !% in units of the Brillouin zone divided by the number in the first row.
428 !% Please note that without specifying any shift, an implicit shift of -0.5 is added for odd number
429 !% of k-points, such that the Gamma point is always included.
430 !%
431 !% For example, the following input samples the BZ with 100 points in the
432 !% <i>xy</i>-plane of reciprocal space:
433 !%
434 !% <tt>%KPointsGrid
435 !% <br>&nbsp;&nbsp;10 | 10 | 1
436 !% <br>%</tt>
437 !%
438 !% is equivalent to
439 !%
440 !% <tt>%KPointsGrid
441 !% <br> 10 | 10 | 1
442 !% <br> 0 | 0 | -0.5
443 !% <br>%</tt>
444 !%
445 !%End
446
447 gamma_only_ = gamma_only
448 if (.not. gamma_only_) then
449 gamma_only_ = (parse_block(namespace, 'KPointsGrid', blk) /= 0)
450 end if
451
452 this%nik_axis = 1
453
454 if (.not. gamma_only_) then
455 nshifts = max(parse_block_n(blk)-1,1)
456 else
457 nshifts = 1
458 end if
459 safe_allocate(shifts(1:dim,1:nshifts))
460 shifts = m_zero
461
462 if (.not. gamma_only_) then
463 ncols = parse_block_cols(blk, 0)
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
467 call messages_fatal(1, namespace=namespace)
468 else
469 write(message(2),'(a)') 'Continuing, but ignoring the additional values.'
470 call messages_warning(2, namespace=namespace)
471 end if
472 end if
473 do ii = 1, dim
474 call parse_block_integer(blk, 0, ii - 1, this%nik_axis(ii))
475 end do
476
477 if (any(this%nik_axis < 1)) then
478 message(1) = 'Input: KPointsGrid is not valid.'
479 call messages_fatal(1, namespace=namespace)
480 end if
481
482 if (parse_block_n(blk) > 1) then ! we have a shift, or even more
483 ncols = parse_block_cols(blk, 1)
484 if (ncols /= dim) then
485 write(message(1),'(a,i3,a,i3)') 'KPointsGrid shift has ', ncols, ' columns but must have ', dim
486 call messages_fatal(1, namespace=namespace)
487 end if
488 do is = 1, nshifts
489 do ii = 1, dim
490 call parse_block_float(blk, is, ii - 1, shifts(ii,is))
491 end do
492 end do
493 else
494 !We include a default 0.5 shift for odd number of k-points
495 do ii = 1, dim
496 if (mod(this%nik_axis(ii), 2) /= 0) then
497 shifts(ii,1) = m_half
498 end if
499 end do
500 end if
501
502 call parse_block_end(blk)
503 else
504 shifts(:, 1) = -m_half
505 end if
506
507 !%Variable QPointsGrid
508 !%Type block
509 !%Default KPointsGrid
510 !%Section Mesh::KPoints
511 !%Description
512 !% This block allows to define a q-point grid used for the calculation of the Fock operator
513 !% with k-points. The <i>q</i>-points are distributed in a uniform grid, as done for the
514 !% <tt>KPointsGrid</tt> variable.
515 !% See J. Chem Phys. 124, 154709 (2006) for details
516 !%
517 !% For each dimension, the number of q point must be a divider of the number of k point
518 !%
519 !% <tt>%QPointsGrid
520 !% <br>&nbsp;&nbsp;2 | 2 | 1
521 !% <br>%</tt>
522 !%
523 !% At the moment, this is not compatible with k-point symmetries.
524 !%
525 !%End
526 this%niq_axis(:) = this%nik_axis(:)
527 this%downsampling(:) = 1
528 if (parse_is_defined(namespace, 'QPointsGrid')) then
529 if (parse_block(namespace, 'QPointsGrid', blk) == 0) then
530 ncols = parse_block_cols(blk, 0)
531 if (ncols /= dim) then
532 write(message(1),'(a,i3,a,i3)') 'QPointsGrid first row has ', ncols, ' columns but must have ', dim
533 call messages_fatal(1, namespace=namespace)
534 end if
535 do ii = 1, dim
536 call parse_block_integer(blk, 0, ii - 1, this%niq_axis(ii))
537 end do
538
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.'
541 call messages_fatal(1, namespace=namespace)
542 end if
543
544 this%downsampling = this%nik_axis/this%niq_axis
545
546 if (any(this%downsampling /= 1)) then
547 call messages_not_implemented('QPointsGrid together with k-point symmetries', namespace=namespace)
548 end if
549
550 call parse_block_end(blk)
551 end if
552 end if
553
554 call kpoints_grid_init(dim, this%full, product(this%nik_axis)*nshifts, nshifts)
555
556 !We move the k-points into this%shifts
557 this%full%shifts(:, :) = shifts(:, :)
558 safe_deallocate_a(shifts)
559
560 call kpoints_grid_generate(dim, this%nik_axis, this%full%nshifts, this%full%shifts, this%full%red_point)
561
562 do ik = 1, this%full%npoints
563 call kpoints_to_absolute(this%latt, this%full%red_point(:, ik), this%full%point(:, ik))
564 end do
565
566 this%full%weight = m_one / this%full%npoints
567
568 if (this%use_symmetries) then
569 message(1) = "Checking if the generated full k-point grid is symmetric"
570 call messages_info(1, namespace=namespace)
571 call kpoints_check_symmetries(this%full, symm, dim, this%use_time_reversal, namespace)
572 end if
573
574 call kpoints_grid_copy(this%full, this%reduced)
575
576
577 if (this%use_symmetries) then
578
579 safe_allocate(num_symm_ops(1:this%full%npoints))
580
581 if (this%use_time_reversal) then
582 safe_allocate(symm_ops(1:this%full%npoints, 1:2*symmetries_number(symm)))
583 else
584 safe_allocate(symm_ops(1:this%full%npoints, 1:symmetries_number(symm)))
585 end if
586
587 call kpoints_grid_reduce(symm, this%use_time_reversal, &
588 this%reduced%npoints, dim, this%reduced%red_point, this%reduced%weight, symm_ops, num_symm_ops)
589
590 ! sanity checks
591 assert(maxval(num_symm_ops) >= 1)
592 if (this%use_time_reversal) then
593 assert(maxval(num_symm_ops) <= 2*symmetries_number(symm))
594 else
595 assert(maxval(num_symm_ops) <= symmetries_number(symm))
596 end if
597 ! the total number of symmetry operations in the list has to be equal to the number of k-points
598 assert(sum(num_symm_ops(1:this%reduced%npoints)) == this%full%npoints)
599
600 do ik = 1, this%reduced%npoints
601 assert(all(symm_ops(ik, 1:num_symm_ops(ik)) <= symm%nops))
602 end do
603
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))
606
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))
610
611 safe_deallocate_a(num_symm_ops)
612 safe_deallocate_a(symm_ops)
613
614 else
615
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
620
621 end if
622
623 do ik = 1, this%reduced%npoints
624 call kpoints_to_absolute(this%latt, this%reduced%red_point(:, ik), this%reduced%point(:, ik))
625 end do
626
627 call kpoints_fold_to_1bz(this%full, this%latt)
628 call kpoints_fold_to_1bz(this%reduced, this%latt)
629
630 pop_sub(kpoints_init.read_mp)
631 end subroutine read_mp
632
633 ! ---------------------------------------------------------
635 subroutine read_path()
636 type(block_t) :: blk
637 integer :: nshifts, nkpoints, nhighsympoints, nsegments
638 integer :: icol, ik, idir, ncols
639 integer, allocatable :: resolution(:)
640 real(real64), allocatable :: highsympoints(:, :)
641 type(kpoints_grid_t) :: path_kpoints_grid
642
643 integer npoints
644 integer, allocatable :: symm_ops(:, :), num_symm_ops(:)
645
646 push_sub(kpoints_init.read_path)
647
648
649 !%Variable KPointsPath
650 !%Type block
651 !%Section Mesh::KPoints
652 !%Description
653 !% When this block is given, <i>k</i>-points are generated along a path
654 !% defined by the points of the list.
655 !% The points must be given in reduced coordinates.
656 !%
657 !% The first row of the block is a set of integers defining
658 !% the number of <i>k</i>-points for each segments of the path.
659 !% The number of columns should be equal to <tt>Dimensions</tt>,
660 !% and the k-points coordinate should be zero in finite directions.
661 !%
662 !% For example, the following input samples the BZ with 15 points:
663 !%
664 !% <tt>%KPointsPath
665 !% <br>&nbsp;&nbsp;10 | 5
666 !% <br>&nbsp;&nbsp; 0 | 0 | 0
667 !% <br>&nbsp;&nbsp; 0.5 | 0 | 0
668 !% <br>&nbsp;&nbsp; 0.5 | 0.5 | 0.5
669 !% <br>%</tt>
670 !%
671 !%End
672
673 if (parse_block(namespace, 'KPointsPath', blk) /= 0) then
674 write(message(1),'(a)') 'Internal error while reading KPointsPath.'
675 call messages_fatal(1, namespace=namespace)
676 end if
677
678 ! There is one high symmetry k-point per line
679 nsegments = parse_block_cols(blk, 0)
680 nhighsympoints = parse_block_n(blk) - 1
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.'
683 call messages_fatal(1, namespace=namespace)
684 end if
685
686 safe_allocate(resolution(1:nsegments))
687 do icol = 1, nsegments
688 call parse_block_integer(blk, 0, icol-1, resolution(icol))
689 end do
690 !Total number of points in the segment
691 nkpoints = sum(resolution) + 1
692
693 safe_allocate(highsympoints(1:dim, 1:nhighsympoints))
694 do ik = 1, nhighsympoints
695 !Sanity check
696 ncols = parse_block_cols(blk, ik)
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
699 call messages_fatal(1, namespace=namespace)
700 end if
701
702 do idir = 1, dim
703 call parse_block_float(blk, ik, idir-1, highsympoints(idir, ik))
704 end do
705 end do
706
707 call parse_block_end(blk)
708
709 !We do not have shifts
710 nshifts = 1
711 call kpoints_grid_init(dim, path_kpoints_grid, nkpoints, nshifts)
712
713 ! For the output of band-structures
714 safe_allocate(this%coord_along_path(1:nkpoints))
715
716 call kpoints_path_generate(dim, this%latt, nkpoints, nsegments, resolution, highsympoints, path_kpoints_grid%point, &
717 this%coord_along_path)
718
719 safe_deallocate_a(resolution)
720 safe_deallocate_a(highsympoints)
721
722 !Use zero weight for the path so that it can be used with any kind of calculation mode
723 !without interfering with the physical BZ integral.
724 if (this%method == kpoints_path) then
725 path_kpoints_grid%weight = m_one/path_kpoints_grid%npoints
726 else
727 path_kpoints_grid%weight = m_zero
728 this%nik_skip = this%nik_skip + path_kpoints_grid%npoints
729 end if
730
731
732 !The points have been generated in absolute coordinates
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))
735 end do
736
737 call kpoints_fold_to_1bz(path_kpoints_grid, this%latt)
738
739 !We need to copy the arrays containing the information on the symmetries
740 !Before calling kpoints_grid_addto
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)))
745
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))
749 end if
750
751 call kpoints_grid_addto(this%full , path_kpoints_grid)
752 call kpoints_grid_addto(this%reduced, path_kpoints_grid)
753
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)))
759
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
765
766 safe_deallocate_a(num_symm_ops)
767 safe_deallocate_a(symm_ops)
768 else if(this%use_symmetries) then ! If symmetries are activated but no KPointsGrid specified
769 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
770 safe_allocate(this%symmetry_ops(1:this%reduced%npoints,1:1))
771
772 this%num_symmetry_ops(1:this%reduced%npoints) = 1
773 this%symmetry_ops(1:this%reduced%npoints, 1) = 1
774 end if
775
776
777 call kpoints_grid_end(path_kpoints_grid)
778
779 pop_sub(kpoints_init.read_path)
780 end subroutine read_path
781
782 ! ---------------------------------------------------------
783 subroutine read_user_kpoints()
784 type(block_t) :: blk
785 logical :: reduced
786 integer :: ik, idir
787 type(kpoints_grid_t) :: user_kpoints_grid
788
790
791 !%Variable KPoints
792 !%Type block
793 !%Section Mesh::KPoints
794 !%Description
795 !% This block defines an explicit set of <i>k</i>-points and their weights for
796 !% a periodic-system calculation. The first column is the weight
797 !% of each <i>k</i>-point and the following are the components of the <i>k</i>-point
798 !% vector. You only need to specify the components for the
799 !% periodic directions. Note that the <i>k</i>-points should be given in
800 !% Cartesian coordinates (not in reduced coordinates), in the units of inverse length.
801 !% The weights will be renormalized so they sum to 1 (and must be rational numbers).
802 !%
803 !% For example, if you want to include only the Gamma point, you can
804 !% use:
805 !%
806 !% <tt>%KPoints
807 !% <br>&nbsp;&nbsp;1.0 | 0 | 0 | 0
808 !% <br>%</tt>
809 !%
810 !%End
811
812 !%Variable KPointsReduced
813 !%Type block
814 !%Section Mesh::KPoints
815 !%Description
816 !% Same as the block <tt>KPoints</tt> but this time the input is given in reduced
817 !% coordinates, <i>i.e.</i>
818 !% what <tt>Octopus</tt> writes in a line in the ground-state standard output as
819 !%
820 !% <tt>#k = 1, k = ( 0.154000, 0.154000, 0.154000)</tt>.
821 !%End
822
823 reduced = .false.
824 if (parse_block(namespace, 'KPoints', blk) /= 0) then
825 if (parse_block(namespace, 'KPointsReduced', blk) == 0) then
826 reduced = .true.
827 else
828 ! This case should really never happen. But why not dying otherwise?!
829 write(message(1),'(a)') 'Internal error loading user-defined k-point list.'
830 call messages_fatal(1, namespace=namespace)
831 end if
832 end if
833
834! ! end the one initialized by KPointsGrid already
835! call kpoints_end(this)
836!
837 call kpoints_grid_init(dim, user_kpoints_grid, parse_block_n(blk), 1)
838
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
843
844
845 if (reduced) then
846 do ik = 1, user_kpoints_grid%npoints
847 call parse_block_float(blk, ik - 1, 0, user_kpoints_grid%weight(ik))
848 do idir = 1, dim
849 call parse_block_float(blk, ik - 1, idir, user_kpoints_grid%red_point(idir, ik))
850 end do
851 ! generate also the absolute coordinates
852 call kpoints_to_absolute(this%latt, user_kpoints_grid%red_point(:, ik), user_kpoints_grid%point(:, ik))
853 end do
854 else
855 do ik = 1, user_kpoints_grid%npoints
856 call parse_block_float(blk, ik - 1, 0, user_kpoints_grid%weight(ik))
857 do idir = 1, dim
858 call parse_block_float(blk, ik - 1, idir, user_kpoints_grid%point(idir, ik), unit_one/units_inp%length)
859 end do
860 ! generate also the reduced coordinates
861 call kpoints_to_reduced(this%latt, user_kpoints_grid%point(:, ik), user_kpoints_grid%red_point(:, ik))
862 end do
863 end if
864 call parse_block_end(blk)
865
866 this%nik_skip = 0
867 if (any(user_kpoints_grid%weight < m_epsilon)) then
868 call messages_experimental('K-points with zero weight')
869 message(1) = "Found k-points with zero weight. They are excluded from density calculation"
870 call messages_warning(1, namespace=namespace)
871 ! count k-points with zero weight and make sure the points are given in
872 ! a block after all regular k-points. This is for convenience, so they can be skipped
873 ! easily and not a big restraint for the user who has to provide the k-points
874 ! explicitly anyway.
875 do ik=1,user_kpoints_grid%npoints
876 if (user_kpoints_grid%weight(ik) < m_epsilon) then
877 ! check there are no points with positive weight following a zero weighted one
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"
881 call messages_fatal(1, namespace=namespace)
882 end if
883 end if
884 this%nik_skip = this%nik_skip + 1
885 ! set to machine zero
886 user_kpoints_grid%weight(ik) = m_zero
887 end if
888 end do
889 end if
890 ! renormalize weights
891 weight_sum = sum(user_kpoints_grid%weight)
892 if (weight_sum < m_epsilon) then
893 message(1) = "k-point weights must sum to a positive number."
894 call messages_fatal(1, namespace=namespace)
895 end if
896 user_kpoints_grid%weight = user_kpoints_grid%weight / weight_sum
897
898 ! for the moment we do not apply symmetries to user kpoints
899 ! call kpoints_grid_copy(this%full, this%reduced)
900
901 call kpoints_fold_to_1bz(user_kpoints_grid, this%latt)
902
903 call kpoints_grid_addto(this%full , user_kpoints_grid)
904 call kpoints_grid_addto(this%reduced, user_kpoints_grid)
905
906
907 write(message(1), '(a,i4,a)') 'Input: ', user_kpoints_grid%npoints, ' k-points were read from the input file'
908 call messages_info(1, namespace=namespace)
909
910 call kpoints_grid_end(user_kpoints_grid)
911
913 end subroutine read_user_kpoints
914
915 end subroutine kpoints_init
916
917 ! ---------------------------------------------------------
918 subroutine kpoints_end(this)
919 type(kpoints_t), intent(inout) :: this
920
921 push_sub(kpoints_end)
922
923 call kpoints_grid_end(this%full)
924 call kpoints_grid_end(this%reduced)
925
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)
932
933 pop_sub(kpoints_end)
934 end subroutine kpoints_end
935
936 ! ---------------------------------------------------------
937 subroutine kpoints_to_absolute(latt, kin, kout)
938 type(lattice_vectors_t), intent(in) :: latt
939 real(real64), intent(in) :: kin(:)
940 real(real64), intent(out) :: kout(:)
941
942 push_sub(kpoints_to_absolute)
943
944 kout = matmul(latt%klattice, kin)
945
946 pop_sub(kpoints_to_absolute)
947 end subroutine kpoints_to_absolute
948
949 ! ---------------------------------------------------------
950 subroutine kpoints_to_reduced(latt, kin, kout)
951 type(lattice_vectors_t), intent(in) :: latt
952 real(real64), intent(in) :: kin(:)
953 real(real64), intent(out) :: kout(:)
954
955 push_sub(kpoints_to_reduced)
956
957 kout = matmul(kin, latt%rlattice) / (m_two*m_pi)
958
959 pop_sub(kpoints_to_reduced)
960 end subroutine kpoints_to_reduced
961
962
963 ! ---------------------------------------------------------
964 subroutine kpoints_copy(kin, kout)
965 type(kpoints_t), intent(in) :: kin
966 type(kpoints_t), intent(inout) :: kout
967
968 push_sub(kpoints_copy)
969
970 call kpoints_end(kout)
971
972 kout%method = kin%method
973
974 call kpoints_grid_copy(kin%full, kout%full)
975 call kpoints_grid_copy(kin%reduced, kout%reduced)
976
977 kout%use_symmetries = kin%use_symmetries
978 kout%use_time_reversal = kin%use_time_reversal
979
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(:)
986
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)
990 end if
991
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))
998 end if
999
1000 pop_sub(kpoints_copy)
1001 end subroutine kpoints_copy
1002
1003
1004 ! ----------------------------------------------------------
1005 integer pure function kpoints_number(this) result(number)
1006 type(kpoints_t), intent(in) :: this
1007
1008 number = this%reduced%npoints
1009
1010 end function kpoints_number
1012
1013 ! ----------------------------------------------------------
1014 pure function kpoints_get_point(this, ik, absolute_coordinates) result(point)
1015 class(kpoints_t), intent(in) :: this
1016 integer, intent(in) :: ik
1017 logical, optional, intent(in) :: absolute_coordinates
1018 real(real64) :: point(1:this%full%dim)
1019
1020 if (optional_default(absolute_coordinates, .true.)) then
1021 point(1:this%full%dim) = this%reduced%point(1:this%full%dim, ik)
1022 else
1023 point(1:this%full%dim) = this%reduced%red_point(1:this%full%dim, ik)
1024 end if
1025
1026 end function kpoints_get_point
1027
1028
1029 ! ----------------------------------------------------------
1030 real(real64) pure function kpoints_get_weight(this, ik) result(weight)
1031 class(kpoints_t), intent(in) :: this
1032 integer, intent(in) :: ik
1033
1034 weight = this%reduced%weight(ik)
1035
1036 end function kpoints_get_weight
1037
1038
1039 ! ----------------------------------------------------------
1047 subroutine kpoints_grid_generate(dim, naxis, nshifts, shift, kpoints, lk123)
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(:, :)
1054 ! !< running from 0 to naxis(1:3).
1055
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(:, :)
1060
1061 push_sub(kpoints_grid_generate)
1062
1063 dx = m_one/(m_two*naxis)
1064
1065 npoints = product(naxis)
1066
1067 if (present(lk123)) then
1068 safe_allocate(lk123_(1:npoints*nshifts,1:dim))
1069 safe_allocate(idx(1:npoints*nshifts))
1070 end if
1071
1072 do is = 1, nshifts
1073 do ii = 0, npoints - 1
1074 ik = npoints*is - ii
1075 jj = ii
1076 divisor = npoints
1077
1078 do idir = 1, dim
1079 divisor = divisor / naxis(idir)
1080 ix(idir) = jj / divisor + 1
1081 jj = mod(jj, divisor)
1082
1083 kpoints(idir, ik) = (m_two*ix(idir) - m_one*naxis(idir) + m_two*shift(idir,is))*dx(idir)
1084 !A default shift of +0.5 is including in case if (mod(naxis(idir), 2) /= 0)
1085
1086 !bring back point to first Brillouin zone, except for points at 1/2
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
1089 end if
1090
1091 end do
1092 if (present(lk123)) then
1093 lk123_(ik, :) = ix
1094 idx(ik) = ik
1095 end if
1096 end do
1097 end do
1099 ! sort the k-points
1100
1101 safe_allocate(nrm(1:npoints*nshifts))
1102 safe_allocate(shell(1:npoints*nshifts))
1103 safe_allocate(coords(1:dim, 1:npoints*nshifts))
1104
1105 do ik = 1, npoints*nshifts
1106 shell(ik) = sum((kpoints(:, ik)/dx)**2)
1107 do idir = 1, dim
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)
1111 end do
1112 end do
1113
1114 nrm = m_zero
1115
1116 maxcoord = m_one
1117 do idir = 1, dim
1118 do ik = 1, npoints*nshifts
1119 nrm(ik) = nrm(ik) + coords(idir, ik)*maxcoord
1120 end do
1121 maxcoord = maxcoord*max(m_one, maxval(coords(idir, 1:npoints*nshifts)))
1122 end do
1124 do ik = 1, npoints*nshifts
1125 nrm(ik) = nrm(ik) + shell(ik)*maxcoord
1126 end do
1127
1128 if (present(lk123)) then
1129 call sort(nrm, idx)
1130 do ik = 1, npoints*nshifts
1131 lk123(ik,:) = lk123_(idx(ik),:)
1132 end do
1133 safe_deallocate_a(lk123_)
1134 safe_deallocate_a(idx)
1135 end if
1136
1137 call sort(nrm, kpoints)
1138
1139
1140 safe_deallocate_a(nrm)
1141 safe_deallocate_a(shell)
1142 safe_deallocate_a(coords)
1143
1144 pop_sub(kpoints_grid_generate)
1145 end subroutine kpoints_grid_generate
1146
1147 ! --------------------------------------------------------------------------------------------
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(:,:)
1156 real(real64), intent(out) :: kpoints(1:dim, 1:nkpoints)
1157 real(real64), intent(out) :: coord(1:nkpoints)
1158
1159 integer :: is, ik, kpt_ind
1160 real(real64) :: length, total_length, accumulated_length
1161 real(real64) :: kpt1(dim), kpt2(dim), vec(dim)
1162
1163 push_sub(kpoints_path_generate)
1164
1165 total_length = m_zero
1166 assert(ubound(highsympoints, dim=2) >= nsegments+1)
1167 !We first compute the total length of the k-point path
1168 do is = 1, nsegments
1169 ! We need to work in abolute coordinates to get the correct path length
1170 call kpoints_to_absolute(latt, highsympoints(:,is), kpt1)
1171 call kpoints_to_absolute(latt, highsympoints(:,is+1), kpt2)
1172
1173 vec = kpt2 - kpt1
1174 length = norm2(vec)
1175 if (resolution(is) > 0) total_length = total_length + length
1176 end do
1177
1178 accumulated_length = m_zero
1179 kpt_ind = 0
1180 !Now we generate the points
1181 do is = 1, nsegments
1182 ! We need to work in abolute coordinates to get the correct path length
1183 call kpoints_to_absolute(latt, highsympoints(:, is), kpt1)
1184 call kpoints_to_absolute(latt, highsympoints(:, is+1), kpt2)
1185
1186 vec = kpt2 - kpt1
1187 length = norm2(vec)
1188 vec = vec/length
1189
1190 do ik = 1, resolution(is)
1191 kpt_ind = kpt_ind +1
1192 coord(kpt_ind) = accumulated_length + (ik-1)*length/resolution(is)
1193 kpoints(:, kpt_ind) = kpt1 + (ik-1)*length/resolution(is)*vec
1194 end do
1195 if (resolution(is) > 0) accumulated_length = accumulated_length + length
1196 end do
1197 !We add the last point
1198 kpt_ind = kpt_ind +1
1199 call kpoints_to_absolute(latt, highsympoints(:,nsegments+1), kpt1)
1200 coord(kpt_ind) = accumulated_length
1201 kpoints(:, kpt_ind) = kpt1
1202
1203 !The length of the total path is arbitrarily put to 1
1204 coord = coord/total_length
1205
1206 pop_sub(kpoints_path_generate)
1207 end subroutine kpoints_path_generate
1208
1209
1210 ! --------------------------------------------------------------------------------------------
1211 subroutine kpoints_grid_reduce(symm, time_reversal, nkpoints, dim, kpoints, weights, symm_ops, num_symm_ops)
1212 type(symmetries_t), intent(in) :: symm
1213 logical, intent(in) :: time_reversal
1214 integer, intent(inout) :: nkpoints
1215 integer, intent(in) :: dim
1216 real(real64), intent(inout) :: kpoints(1:dim,1:nkpoints)
1217 real(real64), intent(out) :: weights(1:nkpoints)
1218 integer, intent(out) :: symm_ops(:, :)
1219 integer, intent(out) :: num_symm_ops(:)
1220
1221 integer :: nreduced
1222 real(real64), allocatable :: reduced(:, :)
1223
1224 integer ik, iop, ik2
1225 real(real64) :: tran(dim), diff(dim)
1226 real(real64), allocatable :: kweight(:)
1227
1228 real(real64) :: prec
1229
1230 push_sub(kpoints_grid_reduce)
1231
1232 ! In case of really dense k-point grid, 1/nkpoints is might be smaller
1233 ! the symprec, causing problems
1234 ! Therefore we use PREC in the following
1235 prec = min(symprec, m_one/(nkpoints*100))
1236
1237 ! reduce to irreducible zone
1238
1239 ! kmap is used to mark reducible k-points and also to
1240 ! map reducible to irreducible k-points
1241
1242 safe_allocate(kweight(1:nkpoints))
1243 safe_allocate(reduced(1:dim, 1:nkpoints))
1244
1245 kweight = m_one / nkpoints
1246
1247 nreduced = 0
1248
1249 num_symm_ops = 1
1250 symm_ops(:, 1) = symmetries_identity_index(symm)
1251
1252 do ik = 1, nkpoints
1253 if (kweight(ik) < prec) cycle
1254
1255 ! new irreducible point
1256 ! has reduced non-zero weight
1257 nreduced = nreduced + 1
1258 reduced(:, nreduced) = kpoints(:, ik)
1259
1260 !No need to check Gamma
1261 if (maxval(abs(kpoints(:, ik))) < m_epsilon) cycle
1262
1263 if (ik == nkpoints) cycle
1264
1265 ! operate with the symmetry operations
1266 do iop = 1, symmetries_number(symm)
1267 if (iop == symmetries_identity_index(symm) .and. &
1268 .not. time_reversal) cycle ! no need to check for the identity
1269
1270 call symmetries_apply_kpoint_red(symm, iop, reduced(1:dim, nreduced), tran)
1271 !We remove potential umklapp
1272 tran = tran - anint(tran + m_half*prec)
1273
1274 ! remove (mark) k-points related to irreducible reduced by symmetry
1275 do ik2 = ik + 1, nkpoints
1276 if (kweight(ik2) < prec) cycle
1277
1278 if (.not. iop == symmetries_identity_index(symm)) then ! no need to check for the identity
1279 diff = tran - kpoints(:, ik2)
1280 diff = diff - anint(diff)
1281
1282 ! both the transformed rk ...
1283 if (sum(abs(diff)) < prec) then
1284 kweight(ik) = kweight(ik) + kweight(ik2)
1285 kweight(ik2) = m_zero
1286 weights(nreduced) = kweight(ik)
1287 num_symm_ops(nreduced) = num_symm_ops(nreduced) + 1
1288 symm_ops(nreduced, num_symm_ops(nreduced)) = iop
1289 cycle
1290 end if
1291 end if
1292
1293 if (time_reversal) then
1294 diff = tran + kpoints(:, ik2)
1295 diff = diff - anint(diff)
1296
1297 ! and its inverse
1298 if (sum(abs(diff)) < prec) then
1299 kweight(ik) = kweight(ik) + kweight(ik2)
1300 kweight(ik2) = m_zero
1301 weights(nreduced) = kweight(ik)
1302 num_symm_ops(nreduced) = num_symm_ops(nreduced) + 1
1303 !We mark the symmetry+time-reversal operation as negative
1304 symm_ops(nreduced, num_symm_ops(nreduced)) = -iop
1305 end if
1306 end if
1307 end do
1308 end do
1309 end do
1310
1311 assert(sum(weights(1:nreduced)) - m_one < prec)
1312
1313 nkpoints = nreduced
1314 do ik = 1, nreduced
1315 kpoints(:, ik) = reduced(:, ik)
1316 end do
1317
1318 safe_deallocate_a(kweight)
1319 safe_deallocate_a(reduced)
1320
1321 pop_sub(kpoints_grid_reduce)
1322 end subroutine kpoints_grid_reduce
1323
1324
1325 ! --------------------------------------------------------------------------------------------
1326 subroutine kpoints_fold_to_1bz(grid, latt)
1327 type(kpoints_grid_t), intent(inout) :: grid
1328 type(lattice_vectors_t), intent(in) :: latt
1329
1330 integer :: ig1, ig2, ik
1331 real(real64) :: Gvec(grid%dim, 3**grid%dim), Gvec_cart(grid%dim, 3**grid%dim)
1332 real(real64) :: vec(grid%dim), kpt(grid%dim)
1333 real(real64) :: d, dmin
1334
1335 push_sub(kpoints_fold_to_1bz)
1336
1337 !We only need to compute the first G-vectors
1338 do ig1 = 1, grid%dim
1339 do ig2 = 1, 3**grid%dim
1340 gvec(ig1, ig2) = modulo((ig2 - 1)/grid%dim**(grid%dim - ig1), 3) - 1
1341 end do
1342 end do
1343
1344 do ig1 = 1, 3**grid%dim
1345 call kpoints_to_absolute(latt, gvec(:,ig1), gvec_cart(:,ig1))
1346 end do
1347
1348 do ik = 1, grid%npoints
1349
1350 dmin = 1e10_real64
1351 do ig1 = 1, 3**grid%dim
1352 vec = gvec_cart(:,ig1) - grid%point(:,ik)
1353 d = real(sum(vec**2), real32) !Conversion to simple precision
1354 ! To avoid numerical error problems
1355 if (d < dmin) then
1356 dmin = d
1357 ig2 = ig1
1358 end if
1359 end do
1360 kpt = grid%red_point(:,ik) - gvec(:,ig2)
1361 call kpoints_to_absolute(latt, kpt, grid%point1BZ(:, ik))
1362 end do
1363
1364 pop_sub(kpoints_fold_to_1bz)
1365 end subroutine kpoints_fold_to_1bz
1366
1367
1368 ! ---------------------------------------------------------
1369 subroutine kpoints_write_info(this, iunit, namespace, absolute_coordinates)
1370 class(kpoints_t), intent(in) :: this
1371 integer, optional, intent(in) :: iunit
1372 type(namespace_t), optional, intent(in) :: namespace
1373 logical, optional, intent(in) :: absolute_coordinates
1374
1375 integer :: ik, idir
1376 character(len=100) :: str_tmp
1377 character :: index
1378
1379 push_sub(kpoints_write_info)
1380
1381 call messages_print_with_emphasis('Brillouin zone sampling', iunit, namespace)
1382
1383 if (this%method == kpoints_monkh_pack) then
1384
1385 call messages_write('Dimensions of the k-point grid =')
1386 do idir = 1, this%full%dim
1387 call messages_write(this%nik_axis(idir), fmt = '(i3,1x)')
1388 end do
1389 call messages_new_line()
1390
1391 if (.not. all(this%downsampling(1:this%full%dim) == 1)) then
1392 call messages_write('Dimensions of the q-point grid =')
1393 do idir = 1, this%full%dim
1394 call messages_write(this%niq_axis(idir), fmt = '(i3,1x)')
1395 end do
1396 call messages_new_line()
1397 end if
1398
1399 call messages_write('Total number of k-points =')
1400 call messages_write(this%full%npoints)
1401 call messages_new_line()
1402
1403 call messages_write('Number of symmetry-reduced k-points =')
1404 call messages_write(this%reduced%npoints)
1405
1406 call messages_info(iunit=iunit, namespace=namespace)
1407
1408 else
1409
1410 call messages_write('Total number of k-points =')
1411 call messages_write(this%full%npoints)
1412 call messages_new_line()
1413 call messages_info(iunit=iunit, namespace=namespace)
1414
1415 end if
1416
1417 call messages_new_line()
1418 call messages_write('List of k-points:')
1419 call messages_info(iunit=iunit, namespace=namespace)
1420
1421 write(message(1), '(6x,a)') 'ik'
1422 do idir = 1, this%full%dim
1423 index = index2axis(idir)
1424 write(str_tmp, '(9x,2a)') 'k_', index
1425 message(1) = trim(message(1)) // trim(str_tmp)
1426 end do
1427 write(str_tmp, '(6x,a)') 'Weight'
1428 message(1) = trim(message(1)) // trim(str_tmp)
1429 message(2) = '---------------------------------------------------------'
1430 call messages_info(2, iunit)
1431
1432 do ik = 1, kpoints_number(this)
1433 write(message(1),'(i8,1x)') ik
1434 do idir = 1, this%full%dim
1435 if (optional_default(absolute_coordinates, .false.)) then
1436 write(str_tmp,'(f12.4)') this%reduced%point(idir, ik)
1437 else
1438 write(str_tmp,'(f12.4)') this%reduced%red_point(idir, ik)
1439 end if
1440 message(1) = trim(message(1)) // trim(str_tmp)
1441 end do
1442 write(str_tmp,'(f12.4)') kpoints_get_weight(this, ik)
1443 message(1) = trim(message(1)) // trim(str_tmp)
1444 call messages_info(1, iunit, namespace=namespace)
1445 end do
1446
1447 call messages_info(iunit=iunit, namespace=namespace)
1448
1449 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
1450
1451 pop_sub(kpoints_write_info)
1452 end subroutine kpoints_write_info
1453
1454
1455 ! ---------------------------------------------------------
1456 logical pure function kpoints_point_is_gamma(this, ik) result(is_gamma)
1457 class(kpoints_t), intent(in) :: this
1458 integer, intent(in) :: ik
1459
1460 is_gamma = (maxval(abs(kpoints_get_point(this, ik))) < m_epsilon)
1461
1463
1464 !--------------------------------------------------------
1465
1466 integer pure function kpoints_get_num_symmetry_ops(this, ik) result(num)
1467 type(kpoints_t), intent(in) :: this
1468 integer, intent(in) :: ik
1469
1470 if (this%use_symmetries) then
1471 num = this%num_symmetry_ops(ik)
1472 else
1473 num = 1
1474 end if
1475
1476 end function kpoints_get_num_symmetry_ops
1477 !--------------------------------------------------------
1478
1479 integer pure function kpoints_get_symmetry_ops(this, ik, index) result(iop)
1480 type(kpoints_t), intent(in) :: this
1481 integer, intent(in) :: ik
1482 integer, intent(in) :: index
1483
1484 if (this%use_symmetries) then
1485 iop = this%symmetry_ops(ik, index)
1486 else
1487 iop = 1
1488 end if
1489
1490 end function kpoints_get_symmetry_ops
1491
1492 !--------------------------------------------------------
1493 logical pure function kpoints_is_valid_symmetry(this, ik, index) result(valid)
1494 type(kpoints_t), intent(in) :: this
1495 integer, intent(in) :: ik
1496 integer, intent(in) :: index
1497
1498 integer :: iop
1499
1500 valid = .false.
1501 if (this%use_symmetries) then
1502 do iop = 1, this%num_symmetry_ops(ik)
1503 if (this%symmetry_ops(ik, iop) == index) then
1504 valid = .true.
1505 return
1506 end if
1507 end do
1508 else
1509 valid = .true.
1510 end if
1511
1512 end function kpoints_is_valid_symmetry
1513
1514
1515 !--------------------------------------------------------
1516 integer function kpoints_kweight_denominator(this)
1517 type(kpoints_t), intent(in) :: this
1518
1519 integer :: denom, nik, nik_skip
1520
1522
1523 nik = this%full%npoints
1524 nik_skip = this%nik_skip
1525
1526 if (this%method == kpoints_monkh_pack) then
1527 kpoints_kweight_denominator = this%full%npoints
1528 else
1530 ! NB largest reasonable value is: # k-points x 48. from space-group symmetries
1531 do denom = 1, 100000
1532 if (all(abs(int(this%full%weight(1:nik-nik_skip)*denom + 10_real64*m_epsilon) - &
1533 this%full%weight(1:nik-nik_skip)*denom) < 100_real64*m_epsilon)) then
1535 exit
1536 end if
1537 end do
1538 end if
1539
1541 end function kpoints_kweight_denominator
1542
1543 !--------------------------------------------------------
1544 logical pure function kpoints_have_zero_weight_path(this) result(have_zerow)
1545 class(kpoints_t), intent(in) :: this
1546
1547 if (this%nik_skip > 0) then
1548 have_zerow = .true.
1549 else
1550 have_zerow = .false.
1551 end if
1552
1554
1555 !--------------------------------------------------------
1556 integer pure function kpoints_get_kpoint_method(this)
1557 class(kpoints_t), intent(in) :: this
1558
1560 end function kpoints_get_kpoint_method
1561
1562 !--------------------------------------------------------
1563 real(real64) pure function kpoints_get_path_coord(this, ind) result(coord)
1564 type(kpoints_t), intent(in) :: this
1565 integer, intent(in) :: ind
1566
1567 coord = this%coord_along_path(ind)
1568 end function
1569
1570
1571
1572 !--------------------------------------------------------
1573 subroutine kpoints_check_symmetries(grid, symm, dim, time_reversal, namespace)
1574 type(kpoints_grid_t), intent(in) :: grid
1575 type(symmetries_t), intent(in) :: symm
1576 integer, intent(in) :: dim
1577 logical, intent(in) :: time_reversal
1578 type(namespace_t), intent(in) :: namespace
1579
1580 integer, allocatable :: kmap(:)
1581 real(real64) :: kpt(dim), diff(dim)
1582 integer :: nk, ik, ik2, iop
1583 type(distributed_t) :: kpt_dist
1584
1585 push_sub(kpoints_check_symmetries)
1587 nk = grid%npoints
1588
1589 !We distribute the k-points here for this routine, independently of the rest of the code
1590 call distributed_nullify(kpt_dist, nk)
1591#ifdef HAVE_MPI
1592 if (mpi_world%comm /= mpi_comm_undefined) then
1593 call distributed_init(kpt_dist, nk, mpi_comm_world, "kpt_check")
1594 end if
1595#endif
1596
1597 !A simple map to tell if the k-point as a matching symmetric point or not
1598 safe_allocate(kmap(kpt_dist%start:kpt_dist%end))
1599
1600 do iop = 1, symmetries_number(symm)
1601 if (iop == symmetries_identity_index(symm) .and. &
1602 .not. time_reversal) cycle
1603
1604 do ik = kpt_dist%start, kpt_dist%end
1605 kmap(ik) = ik
1606 end do
1607
1608 do ik = kpt_dist%start, kpt_dist%end
1609 !We apply the symmetry
1610 call symmetries_apply_kpoint_red(symm, iop, grid%red_point(:, ik), kpt)
1611 !We remove potential umklapp
1612 kpt = kpt - anint(kpt + m_half*symprec)
1613
1614 ! remove (mark) k-points which already have a symmetric point
1615 do ik2 = 1, nk
1616
1617 if (iop /= symmetries_identity_index(symm)) then
1618 diff = kpt - grid%red_point(:, ik2)
1619 diff = diff - anint(diff)
1620 !We found point corresponding to the symmetric kpoint
1621 if (sum(abs(diff)) < symprec) then
1622 kmap(ik) = -ik2
1623 exit
1624 end if
1625 end if
1626
1627 if (time_reversal) then
1628 diff = kpt + grid%red_point(:, ik2)
1629 diff = diff - anint(diff)
1630 !We found point corresponding to the symmetric kpoint
1631 if (sum(abs(diff)) < symprec) then
1632 kmap(ik) = -ik2
1633 exit
1634 end if
1635 end if
1636
1637 end do
1638 !In case we have not found a symnetric k-point...
1639 if (kmap(ik) == ik) then
1640 write(message(1),'(a,i5,a2,3(f7.3,a2),a)') "The reduced k-point ", ik, " (", &
1641 grid%red_point(1, ik), ", ", grid%red_point(2, ik), ", ", grid%red_point(3, ik), &
1642 ") ", "has no symmetric in the k-point grid for the following symmetry"
1643 write(message(2),'(i5,1x,a,2x,3(3i4,2x))') iop, ':', transpose(symm_op_rotation_matrix_red(symm%ops(iop)))
1644 message(3) = "Change your k-point grid or use KPointsUseSymmetries=no."
1645 call messages_fatal(3, namespace=namespace)
1646 end if
1647 end do
1648 end do
1650 safe_deallocate_a(kmap)
1651
1652 call distributed_end(kpt_dist)
1653
1655 end subroutine kpoints_check_symmetries
1657 !--------------------------------------------------------
1658 logical function kpoints_is_compatible_downsampling(kpt, ik, iq) result(compatible)
1659 type(kpoints_t), intent(in) :: kpt
1660 integer, intent(in) :: ik
1661 integer, intent(in) :: iq
1662
1663 integer :: dim, idim
1664 real(real64) :: diff(kpt%reduced%dim), red(kpt%reduced%dim)
1665
1667
1668 compatible = .true.
1669 dim = kpt%reduced%dim
1670
1671 !No downsampling. We use all k-points
1672 if (all(kpt%downsampling == 1)) then
1674 return
1675 end if
1676
1677 assert(kpt%method == kpoints_monkh_pack)
1678
1679 diff = kpt%reduced%red_point(:, ik) - kpt%reduced%red_point(:, iq)
1680 do idim = 1, dim
1681 !We remove potential umklapp
1682 diff(idim) = diff(idim) - anint(diff(idim) + m_half*symprec)
1683 red(idim) = diff(idim)*kpt%nik_axis(idim)/real(kpt%downsampling(idim))
1684 if (abs(red(idim) - anint(red(idim))) > m_epsilon) then
1685 compatible = .false.
1687 return
1688 end if
1689 end do
1690
1693
1694 ! ---------------------------------------------------------
1695 integer function kpoints_nkpt_in_path(this) result(npath)
1696 class(kpoints_t), intent(in) :: this
1697
1698 npath = 0
1699 if (allocated(this%coord_along_path)) then
1700 npath = SIZE(this%coord_along_path)
1701 end if
1702
1703 end function kpoints_nkpt_in_path
1704
1705 ! ---------------------------------------------------------
1706 logical function kpoints_gamma_only(this) result(gamma_only)
1707 class(kpoints_t), intent(in) :: this
1708
1709 push_sub(kpoints_gamma_only)
1710
1711 gamma_only = kpoints_number(this) == 1 .and. kpoints_point_is_gamma(this, 1)
1712
1713 pop_sub(kpoints_gamma_only)
1714 end function kpoints_gamma_only
1715
1716 ! ---------------------------------------------------------
1717 subroutine kpoints_lattice_vectors_update(this, new_latt)
1718 class(kpoints_t), intent(inout) :: this
1719 type(lattice_vectors_t), intent(in) :: new_latt
1720
1721 integer :: ik
1722
1724
1725 this%latt = new_latt
1726
1727 do ik = 1, this%reduced%npoints
1728 call kpoints_to_absolute(this%latt, this%reduced%red_point(:, ik), this%reduced%point(:, ik))
1729 end do
1730
1731 do ik = 1, this%full%npoints
1732 call kpoints_to_absolute(this%latt, this%full%red_point(:, ik), this%full%point(:, ik))
1733 end do
1734
1736 end subroutine kpoints_lattice_vectors_update
1737
1738end module kpoints_oct_m
1739
1740!! Local Variables:
1741!! mode: f90
1742!! coding: utf-8
1743!! End:
subroutine read_mp(gamma_only)
Definition: kpoints.F90:485
subroutine read_path()
Read the k-points path information and generate the k-points list.
Definition: kpoints.F90:729
subroutine read_user_kpoints()
Definition: kpoints.F90:877
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
integer pure function kpoints_get_kpoint_method(this)
Definition: kpoints.F90:1650
subroutine, public kpoints_grid_generate(dim, naxis, nshifts, shift, kpoints, lk123)
Generates the k-points grid.
Definition: kpoints.F90:1141
subroutine, public kpoints_end(this)
Definition: kpoints.F90:1012
logical pure function kpoints_have_zero_weight_path(this)
Definition: kpoints.F90:1638
subroutine kpoints_check_symmetries(grid, symm, dim, time_reversal, namespace)
Definition: kpoints.F90:1667
integer, parameter, public kpoints_monkh_pack
Definition: kpoints.F90:209
integer pure function, public kpoints_get_num_symmetry_ops(this, ik)
Definition: kpoints.F90:1560
subroutine kpoints_grid_addto(this, that)
Definition: kpoints.F90:273
real(real64) pure function, public kpoints_get_path_coord(this, ind)
Definition: kpoints.F90:1657
logical pure function, public kpoints_is_valid_symmetry(this, ik, index)
Definition: kpoints.F90:1587
integer function, public kpoints_kweight_denominator(this)
Definition: kpoints.F90:1610
subroutine, public kpoints_copy(kin, kout)
Definition: kpoints.F90:1058
subroutine, public kpoints_path_generate(dim, latt, nkpoints, nsegments, resolution, highsympoints, kpoints, coord)
Generate the k-point along a path.
Definition: kpoints.F90:1243
integer, parameter, public kpoints_path
Definition: kpoints.F90:209
pure real(real64) function, dimension(1:this%full%dim) kpoints_get_point(this, ik, absolute_coordinates)
Definition: kpoints.F90:1108
integer function, public kpoints_nkpt_in_path(this)
Definition: kpoints.F90:1789
logical function, public kpoints_is_compatible_downsampling(kpt, ik, iq)
Definition: kpoints.F90:1752
real(real64) pure function kpoints_get_weight(this, ik)
Definition: kpoints.F90:1124
integer pure function, public kpoints_get_symmetry_ops(this, ik, index)
Definition: kpoints.F90:1573
subroutine kpoints_grid_copy(bb, aa)
Definition: kpoints.F90:255
subroutine, public kpoints_fold_to_1bz(grid, latt)
Definition: kpoints.F90:1420
subroutine, public kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
Definition: kpoints.F90:322
subroutine, public kpoints_grid_end(this)
Definition: kpoints.F90:240
logical pure function, public kpoints_point_is_gamma(this, ik)
Definition: kpoints.F90:1550
integer, parameter, public kpoints_user
Definition: kpoints.F90:209
subroutine kpoints_grid_reduce(symm, time_reversal, nkpoints, dim, kpoints, weights, symm_ops, num_symm_ops)
Definition: kpoints.F90:1305
logical function kpoints_gamma_only(this)
Definition: kpoints.F90:1800
subroutine, public kpoints_lattice_vectors_update(this, new_latt)
Definition: kpoints.F90:1811
subroutine, public kpoints_to_reduced(latt, kin, kout)
Definition: kpoints.F90:1044
integer pure function, public kpoints_number(this)
Definition: kpoints.F90:1099
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
subroutine, public kpoints_grid_init(dim, this, npoints, nshifts)
Definition: kpoints.F90:219
subroutine kpoints_write_info(this, iunit, namespace, absolute_coordinates)
Definition: kpoints.F90:1463
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1113
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1085
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:543
logical pure function, public symmetries_have_break_dir(this)
Definition: symmetries.F90:581
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
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.
Definition: utils.F90:118
int true(void)