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