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