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