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, &
651 this%reduced%equiv, opt = .false.)
652 elseif (dos_method == option__dosmethod__tetrahedra_opt) then
653 this%simplex => simplex_t(dim, this%nik_axis, this%full%nshifts, this%full%shifts, this%full%red_point, &
654 this%reduced%equiv, opt = .true.)
655 end if
656
657 pop_sub(kpoints_init.read_mp)
658 end subroutine read_mp
659
660 ! ---------------------------------------------------------
662 subroutine read_path()
663 type(block_t) :: blk
664 integer :: nshifts, nkpoints, nhighsympoints, nsegments
665 integer :: icol, ik, idir, ncols
666 integer, allocatable :: resolution(:)
667 real(real64), allocatable :: highsympoints(:, :)
668 type(kpoints_grid_t) :: path_kpoints_grid
669
670 integer npoints
671 integer, allocatable :: symm_ops(:, :), num_symm_ops(:)
672
673 push_sub(kpoints_init.read_path)
674
675
676 !%Variable KPointsPath
677 !%Type block
678 !%Section Mesh::KPoints
679 !%Description
680 !% When this block is given, <i>k</i>-points are generated along a path
681 !% defined by the points of the list.
682 !% The points must be given in reduced coordinates.
683 !%
684 !% The first row of the block is a set of integers defining
685 !% the number of <i>k</i>-points for each segments of the path.
686 !% The number of columns should be equal to <tt>Dimensions</tt>,
687 !% and the k-points coordinate should be zero in finite directions.
688 !%
689 !% For example, the following input samples the BZ with 15 points:
690 !%
691 !% <tt>%KPointsPath
692 !% <br>&nbsp;&nbsp;10 | 5
693 !% <br>&nbsp;&nbsp; 0 | 0 | 0
694 !% <br>&nbsp;&nbsp; 0.5 | 0 | 0
695 !% <br>&nbsp;&nbsp; 0.5 | 0.5 | 0.5
696 !% <br>%</tt>
697 !%
698 !%End
699
700 if (parse_block(namespace, 'KPointsPath', blk) /= 0) then
701 write(message(1),'(a)') 'Internal error while reading KPointsPath.'
702 call messages_fatal(1, namespace=namespace)
703 end if
704
705 ! There is one high symmetry k-point per line
706 nsegments = parse_block_cols(blk, 0)
707 nhighsympoints = parse_block_n(blk) - 1
708 if (nhighsympoints /= nsegments + 1) then
709 write(message(1),'(a,i3,a,i3)') 'The first row of KPointsPath is not compatible with the number of specified k-points.'
710 call messages_fatal(1, namespace=namespace)
711 end if
712
713 safe_allocate(resolution(1:nsegments))
714 do icol = 1, nsegments
715 call parse_block_integer(blk, 0, icol-1, resolution(icol))
716 end do
717 !Total number of points in the segment
718 nkpoints = sum(resolution) + 1
719
720 safe_allocate(highsympoints(1:dim, 1:nhighsympoints))
721 do ik = 1, nhighsympoints
722 !Sanity check
723 ncols = parse_block_cols(blk, ik)
724 if (ncols /= dim) then
725 write(message(1),'(a,i8,a,i3,a,i3)') 'KPointsPath row ', ik, ' has ', ncols, ' columns but must have ', dim
726 call messages_fatal(1, namespace=namespace)
727 end if
728
729 do idir = 1, dim
730 call parse_block_float(blk, ik, idir-1, highsympoints(idir, ik))
731 end do
732 end do
733
734 call parse_block_end(blk)
735
736 !We do not have shifts
737 nshifts = 1
738 call kpoints_grid_init(dim, path_kpoints_grid, nkpoints, nshifts)
739
740 ! For the output of band-structures
741 safe_allocate(this%coord_along_path(1:nkpoints))
742
743 call kpoints_path_generate(dim, this%latt, nkpoints, nsegments, resolution, highsympoints, path_kpoints_grid%point, &
744 this%coord_along_path)
745
746 safe_deallocate_a(resolution)
747 safe_deallocate_a(highsympoints)
748
749 !Use zero weight for the path so that it can be used with any kind of calculation mode
750 !without interfering with the physical BZ integral.
751 if (this%method == kpoints_path) then
752 path_kpoints_grid%weight = m_one/path_kpoints_grid%npoints
753 else
754 path_kpoints_grid%weight = m_zero
755 this%nik_skip = this%nik_skip + path_kpoints_grid%npoints
756 end if
758
759 !The points have been generated in absolute coordinates
760 do ik = 1, path_kpoints_grid%npoints
761 call kpoints_to_reduced(this%latt, path_kpoints_grid%point(:, ik), path_kpoints_grid%red_point(:, ik))
762 end do
763
764 call kpoints_fold_to_1bz(path_kpoints_grid, this%latt)
765
766 !We need to copy the arrays containing the information on the symmetries
767 !Before calling kpoints_grid_addto
768 if (allocated(this%symmetry_ops)) then
769 npoints = this%reduced%npoints
770 safe_allocate(num_symm_ops(1:npoints))
771 safe_allocate(symm_ops(1:npoints, 1:maxval(this%num_symmetry_ops)))
772
773 num_symm_ops(1:npoints) = this%num_symmetry_ops(1:npoints)
774 symm_ops(1:npoints, 1:maxval(num_symm_ops)) = &
775 this%symmetry_ops(1:npoints, 1:maxval(num_symm_ops))
776 end if
777
778 call kpoints_grid_addto(this%full , path_kpoints_grid)
779 call kpoints_grid_addto(this%reduced, path_kpoints_grid)
780
781 if (allocated(this%symmetry_ops)) then
782 safe_deallocate_a(this%num_symmetry_ops)
783 safe_deallocate_a(this%symmetry_ops)
784 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
785 safe_allocate(this%symmetry_ops(1:this%reduced%npoints,1:maxval(num_symm_ops)))
786
787 this%num_symmetry_ops(1:npoints) = num_symm_ops(1:npoints)
788 this%symmetry_ops(1:npoints, 1:maxval(num_symm_ops)) = &
789 symm_ops(1:npoints, 1:maxval(num_symm_ops))
790 this%num_symmetry_ops(npoints+1:this%reduced%npoints) = 1
791 this%symmetry_ops(npoints+1:this%reduced%npoints, 1) = 1
792
793 safe_deallocate_a(num_symm_ops)
794 safe_deallocate_a(symm_ops)
795 else if(this%use_symmetries) then ! If symmetries are activated but no KPointsGrid specified
796 safe_allocate(this%num_symmetry_ops(1:this%reduced%npoints))
797 safe_allocate(this%symmetry_ops(1:this%reduced%npoints,1:1))
798
799 this%num_symmetry_ops(1:this%reduced%npoints) = 1
800 this%symmetry_ops(1:this%reduced%npoints, 1) = 1
801 end if
802
803
804 call kpoints_grid_end(path_kpoints_grid)
805
806 pop_sub(kpoints_init.read_path)
807 end subroutine read_path
808
809 ! ---------------------------------------------------------
810 subroutine read_user_kpoints()
811 type(block_t) :: blk
812 logical :: reduced
813 integer :: ik, idir
814 type(kpoints_grid_t) :: user_kpoints_grid
815
817
818 !%Variable KPoints
819 !%Type block
820 !%Section Mesh::KPoints
821 !%Description
822 !% This block defines an explicit set of <i>k</i>-points and their weights for
823 !% a periodic-system calculation. The first column is the weight
824 !% of each <i>k</i>-point and the following are the components of the <i>k</i>-point
825 !% vector. You only need to specify the components for the
826 !% periodic directions. Note that the <i>k</i>-points should be given in
827 !% Cartesian coordinates (not in reduced coordinates), in the units of inverse length.
828 !% The weights will be renormalized so they sum to 1 (and must be rational numbers).
829 !%
830 !% For example, if you want to include only the Gamma point, you can
831 !% use:
832 !%
833 !% <tt>%KPoints
834 !% <br>&nbsp;&nbsp;1.0 | 0 | 0 | 0
835 !% <br>%</tt>
836 !%
837 !%End
838
839 !%Variable KPointsReduced
840 !%Type block
841 !%Section Mesh::KPoints
842 !%Description
843 !% Same as the block <tt>KPoints</tt> but this time the input is given in reduced
844 !% coordinates, <i>i.e.</i>
845 !% what <tt>Octopus</tt> writes in a line in the ground-state standard output as
846 !%
847 !% <tt>#k = 1, k = ( 0.154000, 0.154000, 0.154000)</tt>.
848 !%End
849
850 reduced = .false.
851 if (parse_block(namespace, 'KPoints', blk) /= 0) then
852 if (parse_block(namespace, 'KPointsReduced', blk) == 0) then
853 reduced = .true.
854 else
855 ! This case should really never happen. But why not dying otherwise?!
856 write(message(1),'(a)') 'Internal error loading user-defined k-point list.'
857 call messages_fatal(1, namespace=namespace)
858 end if
859 end if
860
861! ! end the one initialized by KPointsGrid already
862! call kpoints_end(this)
863!
864 call kpoints_grid_init(dim, user_kpoints_grid, parse_block_n(blk), 1)
865
866 user_kpoints_grid%red_point = m_zero
867 user_kpoints_grid%point = m_zero
868 user_kpoints_grid%weight = m_zero
869 user_kpoints_grid%equiv = 0
870 user_kpoints_grid%shifts = m_zero
871
872
873 if (reduced) then
874 do ik = 1, user_kpoints_grid%npoints
875 call parse_block_float(blk, ik - 1, 0, user_kpoints_grid%weight(ik))
876 do idir = 1, dim
877 call parse_block_float(blk, ik - 1, idir, user_kpoints_grid%red_point(idir, ik))
878 end do
879 ! generate also the absolute coordinates
880 call kpoints_to_absolute(this%latt, user_kpoints_grid%red_point(:, ik), user_kpoints_grid%point(:, ik))
881 end do
882 else
883 do ik = 1, user_kpoints_grid%npoints
884 call parse_block_float(blk, ik - 1, 0, user_kpoints_grid%weight(ik))
885 do idir = 1, dim
886 call parse_block_float(blk, ik - 1, idir, user_kpoints_grid%point(idir, ik), unit_one/units_inp%length)
887 end do
888 ! generate also the reduced coordinates
889 call kpoints_to_reduced(this%latt, user_kpoints_grid%point(:, ik), user_kpoints_grid%red_point(:, ik))
890 end do
891 end if
892 call parse_block_end(blk)
893
894 this%nik_skip = 0
895 if (any(user_kpoints_grid%weight < m_epsilon)) then
896 call messages_experimental('K-points with zero weight')
897 message(1) = "Found k-points with zero weight. They are excluded from density calculation"
898 call messages_warning(1, namespace=namespace)
899 ! count k-points with zero weight and make sure the points are given in
900 ! a block after all regular k-points. This is for convenience, so they can be skipped
901 ! easily and not a big restraint for the user who has to provide the k-points
902 ! explicitly anyway.
903 do ik=1,user_kpoints_grid%npoints
904 if (user_kpoints_grid%weight(ik) < m_epsilon) then
905 ! check there are no points with positive weight following a zero weighted one
906 if (ik < user_kpoints_grid%npoints) then
907 if (user_kpoints_grid%weight(ik+1) > m_epsilon) then
908 message(1) = "K-points with zero weight must follow all regular k-points in a block"
909 call messages_fatal(1, namespace=namespace)
910 end if
911 end if
912 this%nik_skip = this%nik_skip + 1
913 ! set to machine zero
914 user_kpoints_grid%weight(ik) = m_zero
915 end if
916 end do
917 end if
918 ! renormalize weights
919 weight_sum = sum(user_kpoints_grid%weight)
920 if (weight_sum < m_epsilon) then
921 message(1) = "k-point weights must sum to a positive number."
922 call messages_fatal(1, namespace=namespace)
923 end if
924 user_kpoints_grid%weight = user_kpoints_grid%weight / weight_sum
925
926 ! for the moment we do not apply symmetries to user kpoints
927 ! call kpoints_grid_copy(this%full, this%reduced)
928
929 call kpoints_fold_to_1bz(user_kpoints_grid, this%latt)
930
931 call kpoints_grid_addto(this%full , user_kpoints_grid)
932 call kpoints_grid_addto(this%reduced, user_kpoints_grid)
933
934
935 write(message(1), '(a,i4,a)') 'Input: ', user_kpoints_grid%npoints, ' k-points were read from the input file'
936 call messages_info(1, namespace=namespace)
937
938 call kpoints_grid_end(user_kpoints_grid)
939
941 end subroutine read_user_kpoints
942
943 end subroutine kpoints_init
944
945 ! ---------------------------------------------------------
946 subroutine kpoints_end(this)
947 type(kpoints_t), intent(inout) :: this
948
949 push_sub(kpoints_end)
950
951 call kpoints_grid_end(this%full)
952 call kpoints_grid_end(this%reduced)
953
954 safe_deallocate_a(this%nik_axis)
955 safe_deallocate_a(this%niq_axis)
956 safe_deallocate_a(this%downsampling)
957 safe_deallocate_a(this%symmetry_ops)
958 safe_deallocate_a(this%num_symmetry_ops)
959 safe_deallocate_a(this%coord_along_path)
960 safe_deallocate_p(this%simplex)
961
962 pop_sub(kpoints_end)
963 end subroutine kpoints_end
964
965 ! ---------------------------------------------------------
966 subroutine kpoints_to_absolute(latt, kin, kout)
967 type(lattice_vectors_t), intent(in) :: latt
968 real(real64), intent(in) :: kin(:)
969 real(real64), intent(out) :: kout(:)
970
971 push_sub(kpoints_to_absolute)
972
973 kout = matmul(latt%klattice, kin)
974
975 pop_sub(kpoints_to_absolute)
976 end subroutine kpoints_to_absolute
977
978 ! ---------------------------------------------------------
979 subroutine kpoints_to_reduced(latt, kin, kout)
980 type(lattice_vectors_t), intent(in) :: latt
981 real(real64), intent(in) :: kin(:)
982 real(real64), intent(out) :: kout(:)
983
984 push_sub(kpoints_to_reduced)
985
986 kout = matmul(kin, latt%rlattice) / (m_two*m_pi)
987
988 pop_sub(kpoints_to_reduced)
989 end subroutine kpoints_to_reduced
990
991
992 ! ---------------------------------------------------------
993 subroutine kpoints_copy(kin, kout)
994 type(kpoints_t), intent(in) :: kin
995 type(kpoints_t), intent(inout) :: kout
996
997 push_sub(kpoints_copy)
998
999 call kpoints_end(kout)
1000
1001 kout%method = kin%method
1002
1003 call kpoints_grid_copy(kin%full, kout%full)
1004 call kpoints_grid_copy(kin%reduced, kout%reduced)
1005
1006 kout%use_symmetries = kin%use_symmetries
1007 kout%use_time_reversal = kin%use_time_reversal
1008
1009 safe_allocate(kout%nik_axis(1:kin%full%dim))
1010 safe_allocate(kout%niq_axis(1:kin%full%dim))
1011 safe_allocate(kout%downsampling(1:kin%full%dim))
1012 kout%nik_axis(:) = kin%nik_axis(:)
1013 kout%niq_axis(:) = kin%niq_axis(:)
1014 kout%downsampling(:) = kin%downsampling(:)
1015
1016 if (allocated(kin%coord_along_path)) then
1017 safe_allocate(kout%coord_along_path(1:kin%full%npoints))
1018 kout%coord_along_path(1:kin%full%npoints) = kin%coord_along_path(1:kin%full%npoints)
1019 end if
1020
1021 if (allocated(kin%symmetry_ops)) then
1022 safe_allocate(kout%num_symmetry_ops(1:kin%reduced%npoints))
1023 safe_allocate(kout%symmetry_ops(1:kin%reduced%npoints, 1:maxval(kin%num_symmetry_ops)))
1024 kout%num_symmetry_ops(1:kin%reduced%npoints) = kin%num_symmetry_ops(1:kin%reduced%npoints)
1025 kout%symmetry_ops(1:kin%reduced%npoints, 1:maxval(kin%num_symmetry_ops)) = &
1026 kin%symmetry_ops(1:kin%reduced%npoints, 1:maxval(kin%num_symmetry_ops))
1027 end if
1028
1029 pop_sub(kpoints_copy)
1030 end subroutine kpoints_copy
1031
1032
1033 ! ----------------------------------------------------------
1034 integer pure function kpoints_number(this) result(number)
1035 type(kpoints_t), intent(in) :: this
1036
1037 number = this%reduced%npoints
1038
1039 end function kpoints_number
1040
1042 ! ----------------------------------------------------------
1043 pure function kpoints_get_point(this, ik, absolute_coordinates) result(point)
1044 class(kpoints_t), intent(in) :: this
1045 integer, intent(in) :: ik
1046 logical, optional, intent(in) :: absolute_coordinates
1047 real(real64) :: point(1:this%full%dim)
1048
1049 if (optional_default(absolute_coordinates, .true.)) then
1050 point(1:this%full%dim) = this%reduced%point(1:this%full%dim, ik)
1051 else
1052 point(1:this%full%dim) = this%reduced%red_point(1:this%full%dim, ik)
1053 end if
1054
1055 end function kpoints_get_point
1056
1057
1058 ! ----------------------------------------------------------
1059 real(real64) pure function kpoints_get_weight(this, ik) result(weight)
1060 class(kpoints_t), intent(in) :: this
1061 integer, intent(in) :: ik
1062
1063 weight = this%reduced%weight(ik)
1064
1065 end function kpoints_get_weight
1066
1067
1068 ! ----------------------------------------------------------
1069 integer pure function kpoints_get_equiv(this, ik) result(equiv)
1070 class(kpoints_t), intent(in) :: this
1071 integer, intent(in) :: ik
1072
1073 equiv = this%reduced%equiv(ik)
1075 end function kpoints_get_equiv
1076
1077
1078 ! ----------------------------------------------------------
1086 subroutine kpoints_grid_generate(dim, naxis, nshifts, shift, kpoints, lk123)
1087 integer, intent(in) :: dim
1088 integer, intent(in) :: naxis(1:dim)
1089 integer, intent(in) :: nshifts
1090 real(real64), intent(in) :: shift(:, :)
1091 real(real64), intent(out) :: kpoints(:, :)
1092 integer, optional, intent(out) :: lk123(:, :)
1093 ! !< running from 0 to naxis(1:3).
1094
1095 real(real64) :: dx(dim), maxcoord
1096 integer :: ii, jj, divisor, ik, idir, npoints, is, ix(dim)
1097 integer, allocatable :: lk123_(:, :), idx(:)
1098 real(real64), allocatable :: nrm(:), shell(:), coords(:, :)
1099
1100 push_sub(kpoints_grid_generate)
1101
1102 dx = m_one/(m_two*naxis)
1103
1104 npoints = product(naxis)
1105
1106 safe_allocate(idx(1:npoints*nshifts))
1107 if (present(lk123)) then
1108 safe_allocate(lk123_(1:npoints*nshifts,1:dim))
1109 end if
1110
1111 do is = 1, nshifts
1112 do ii = 0, npoints - 1
1113 ik = npoints*is - ii
1114 jj = ii
1115 divisor = npoints
1116
1117 do idir = 1, dim
1118 divisor = divisor / naxis(idir)
1119 ix(idir) = jj / divisor + 1
1120 jj = mod(jj, divisor)
1121
1122 kpoints(idir, ik) = (m_two*ix(idir) - m_one*naxis(idir) + m_two*shift(idir,is))*dx(idir)
1123 !A default shift of +0.5 is including in case if (mod(naxis(idir), 2) /= 0)
1124
1125 !bring back point to first Brillouin zone, except for points at 1/2
1126 if (abs(kpoints(idir, ik) - m_half) > 1e-14_real64) then
1127 kpoints(idir, ik) = mod(kpoints(idir, ik) + m_half, m_one) - m_half
1128 end if
1130 end do
1131 if (present(lk123)) then
1132 lk123_(ik, :) = ix
1133 idx(ik) = ik
1134 end if
1135 end do
1136 end do
1137
1138 ! sort the k-points
1139
1140 safe_allocate(nrm(1:npoints*nshifts))
1141 safe_allocate(shell(1:npoints*nshifts))
1142 safe_allocate(coords(1:dim, 1:npoints*nshifts))
1143
1144 do ik = 1, npoints*nshifts
1145 shell(ik) = sum((kpoints(:, ik)/dx)**2)
1146 do idir = 1, dim
1147 coords(idir, ik) = kpoints(idir, ik)
1148 if (coords(idir, ik) < -1e-14_real64) coords(idir, ik) = coords(idir, ik) + m_one
1149 coords(idir, ik) = coords(idir, ik)/(dx(idir)*m_two)
1150 end do
1151 end do
1152
1153 nrm = m_zero
1154 maxcoord = m_one
1155 do idir = 1, dim
1156 maxcoord = maxcoord*max(m_one, maxval(coords(idir, 1:npoints*nshifts)))
1157 end do
1158
1159 do ik = 1, npoints*nshifts
1160 nrm(ik) = nrm(ik) + shell(ik)*maxcoord
1161 end do
1162
1163 coords = kpoints
1164 call robust_sort_by_abs(nrm, kpoints, idx)
1165 do ik = 1, npoints*nshifts
1166 kpoints(:,ik) = coords(:,idx(ik))
1167 end do
1168 if (present(lk123)) then
1169 do ik = 1, npoints*nshifts
1170 lk123(ik,:) = lk123_(idx(ik),:)
1171 end do
1172 safe_deallocate_a(lk123_)
1173 end if
1174 safe_deallocate_a(idx)
1175
1176 safe_deallocate_a(nrm)
1177 safe_deallocate_a(shell)
1178 safe_deallocate_a(coords)
1179
1180 pop_sub(kpoints_grid_generate)
1182
1183 ! --------------------------------------------------------------------------------------------
1185 subroutine kpoints_path_generate(dim, latt, nkpoints, nsegments, resolution, highsympoints, kpoints, coord)
1186 integer, intent(in) :: dim
1187 type(lattice_vectors_t), intent(in) :: latt
1188 integer, intent(in) :: nkpoints
1189 integer, intent(in) :: nsegments
1190 integer, intent(in) :: resolution(:)
1191 real(real64), intent(in) :: highsympoints(:,:)
1192 real(real64), intent(out) :: kpoints(1:dim, 1:nkpoints)
1193 real(real64), intent(out) :: coord(1:nkpoints)
1194
1195 integer :: is, ik, kpt_ind
1196 real(real64) :: length, total_length, accumulated_length
1197 real(real64) :: kpt1(dim), kpt2(dim), vec(dim)
1198
1199 push_sub(kpoints_path_generate)
1200
1201 total_length = m_zero
1202 assert(ubound(highsympoints, dim=2) >= nsegments+1)
1203 !We first compute the total length of the k-point path
1204 do is = 1, nsegments
1205 ! We need to work in abolute coordinates to get the correct path length
1206 call kpoints_to_absolute(latt, highsympoints(:,is), kpt1)
1207 call kpoints_to_absolute(latt, highsympoints(:,is+1), kpt2)
1208
1209 vec = kpt2 - kpt1
1210 length = norm2(vec)
1211 if (resolution(is) > 0) total_length = total_length + length
1212 end do
1213
1214 accumulated_length = m_zero
1215 kpt_ind = 0
1216 !Now we generate the points
1217 do is = 1, nsegments
1218 ! We need to work in abolute coordinates to get the correct path length
1219 call kpoints_to_absolute(latt, highsympoints(:, is), kpt1)
1220 call kpoints_to_absolute(latt, highsympoints(:, is+1), kpt2)
1221
1222 vec = kpt2 - kpt1
1223 length = norm2(vec)
1224 vec = vec/length
1225
1226 do ik = 1, resolution(is)
1227 kpt_ind = kpt_ind +1
1228 coord(kpt_ind) = accumulated_length + (ik-1)*length/resolution(is)
1229 kpoints(:, kpt_ind) = kpt1 + (ik-1)*length/resolution(is)*vec
1230 end do
1231 if (resolution(is) > 0) accumulated_length = accumulated_length + length
1232 end do
1233 !We add the last point
1234 kpt_ind = kpt_ind +1
1235 call kpoints_to_absolute(latt, highsympoints(:,nsegments+1), kpt1)
1236 coord(kpt_ind) = accumulated_length
1237 kpoints(:, kpt_ind) = kpt1
1238
1239 !The length of the total path is arbitrarily put to 1
1240 coord = coord/total_length
1241
1242 pop_sub(kpoints_path_generate)
1243 end subroutine kpoints_path_generate
1244
1245
1246 ! --------------------------------------------------------------------------------------------
1247 subroutine kpoints_grid_reduce(symm, time_reversal, nkpoints, dim, kpoints, weights, equiv, symm_ops, num_symm_ops)
1248 type(symmetries_t), intent(in) :: symm
1249 logical, intent(in) :: time_reversal
1250 integer, intent(inout) :: nkpoints
1251 integer, intent(in) :: dim
1252 real(real64), intent(inout) :: kpoints(1:dim,1:nkpoints)
1253 real(real64), intent(out) :: weights(1:nkpoints)
1254 integer, intent(out) :: equiv(1:nkpoints)
1255 integer, intent(out) :: symm_ops(:, :)
1256 integer, intent(out) :: num_symm_ops(:)
1257
1258 integer :: nreduced
1259 real(real64), allocatable :: reduced(:, :)
1260
1261 integer ik, iop, ik2
1262 real(real64) :: tran(dim), diff(dim)
1263 real(real64), allocatable :: kweight(:)
1264
1265 real(real64) :: prec
1266
1267 push_sub(kpoints_grid_reduce)
1268
1269 ! In case of really dense k-point grid, 1/nkpoints is might be smaller
1270 ! the symprec, causing problems
1271 ! Therefore we use PREC in the following
1272 prec = min(symprec, m_one/(nkpoints*100))
1273
1274 ! reduce to irreducible zone
1275
1276 ! kmap is used to mark reducible k-points and also to
1277 ! map reducible to irreducible k-points
1278
1279 safe_allocate(kweight(1:nkpoints))
1280 safe_allocate(reduced(1:dim, 1:nkpoints))
1281
1282 kweight = m_one / nkpoints
1283
1284 nreduced = 0
1285 equiv(:) = -1 ! for safety
1286
1287 num_symm_ops = 1
1288 symm_ops(:, 1) = symmetries_identity_index(symm)
1289
1290 do ik = 1, nkpoints
1291 if (kweight(ik) < prec) cycle
1292
1293 ! new irreducible point
1294 ! has reduced non-zero weight
1295 nreduced = nreduced + 1
1296 reduced(:, nreduced) = kpoints(:, ik)
1297 equiv(ik) = nreduced
1298
1299 !No need to check Gamma
1300 if (maxval(abs(kpoints(:, ik))) < m_epsilon) cycle
1301
1302 if (ik == nkpoints) cycle
1303
1304 ! operate with the symmetry operations
1305 do iop = 1, symmetries_number(symm)
1306 if (iop == symmetries_identity_index(symm) .and. &
1307 .not. time_reversal) cycle ! no need to check for the identity
1308
1309 call symmetries_apply_kpoint_red(symm, iop, reduced(1:dim, nreduced), tran)
1310 !We remove potential umklapp
1311 tran = tran - anint(tran + m_half*prec)
1312
1313 ! remove (mark) k-points related to irreducible reduced by symmetry
1314 do ik2 = ik + 1, nkpoints
1315 if (kweight(ik2) < prec) cycle
1316
1317 if (.not. iop == symmetries_identity_index(symm)) then ! no need to check for the identity
1318 diff = tran - kpoints(:, ik2)
1319 diff = diff - anint(diff)
1320
1321 ! both the transformed rk ...
1322 if (sum(abs(diff)) < prec) then
1323 kweight(ik) = kweight(ik) + kweight(ik2)
1324 kweight(ik2) = m_zero
1325 equiv(ik2) = nreduced
1326 weights(nreduced) = kweight(ik)
1327 num_symm_ops(nreduced) = num_symm_ops(nreduced) + 1
1328 symm_ops(nreduced, num_symm_ops(nreduced)) = iop
1329 cycle
1330 end if
1331 end if
1332
1333 if (time_reversal) then
1334 diff = tran + kpoints(:, ik2)
1335 diff = diff - anint(diff)
1336
1337 ! and its inverse
1338 if (sum(abs(diff)) < prec) then
1339 kweight(ik) = kweight(ik) + kweight(ik2)
1340 kweight(ik2) = m_zero
1341 equiv(ik2) = nreduced
1342 weights(nreduced) = kweight(ik)
1343 num_symm_ops(nreduced) = num_symm_ops(nreduced) + 1
1344 !We mark the symmetry+time-reversal operation as negative
1345 symm_ops(nreduced, num_symm_ops(nreduced)) = -iop
1346 end if
1347 end if
1348 end do
1349 end do
1350 end do
1351
1352 assert(sum(weights(1:nreduced)) - m_one < prec)
1353 assert(all(1 <= equiv(:)) .and. all(equiv(:) <= nreduced))
1354
1355 nkpoints = nreduced
1356 do ik = 1, nreduced
1357 kpoints(:, ik) = reduced(:, ik)
1358 end do
1359
1360 safe_deallocate_a(kweight)
1361 safe_deallocate_a(reduced)
1362
1363 pop_sub(kpoints_grid_reduce)
1364 end subroutine kpoints_grid_reduce
1365
1366
1367 ! --------------------------------------------------------------------------------------------
1368 subroutine kpoints_fold_to_1bz(grid, latt)
1369 type(kpoints_grid_t), intent(inout) :: grid
1370 type(lattice_vectors_t), intent(in) :: latt
1371
1372 integer :: ig1, ig2, ik
1373 real(real64) :: Gvec(grid%dim, 3**grid%dim), Gvec_cart(grid%dim, 3**grid%dim)
1374 real(real64) :: vec(grid%dim), kpt(grid%dim)
1375 real(real64) :: d, dmin
1376
1377 push_sub(kpoints_fold_to_1bz)
1378
1379 !We only need to compute the first G-vectors
1380 do ig1 = 1, grid%dim
1381 do ig2 = 1, 3**grid%dim
1382 gvec(ig1, ig2) = modulo((ig2 - 1)/grid%dim**(grid%dim - ig1), 3) - 1
1383 end do
1384 end do
1385
1386 do ig1 = 1, 3**grid%dim
1387 call kpoints_to_absolute(latt, gvec(:,ig1), gvec_cart(:,ig1))
1388 end do
1389
1390 do ik = 1, grid%npoints
1391
1392 dmin = 1e10_real64
1393 do ig1 = 1, 3**grid%dim
1394 vec = gvec_cart(:,ig1) - grid%point(:,ik)
1395 d = real(sum(vec**2), real32) !Conversion to simple precision
1396 ! To avoid numerical error problems
1397 if (d < dmin) then
1398 dmin = d
1399 ig2 = ig1
1400 end if
1401 end do
1402 kpt = grid%red_point(:,ik) - gvec(:,ig2)
1403 call kpoints_to_absolute(latt, kpt, grid%point1BZ(:, ik))
1404 end do
1405
1406 pop_sub(kpoints_fold_to_1bz)
1407 end subroutine kpoints_fold_to_1bz
1408
1409
1410 ! ---------------------------------------------------------
1411 subroutine kpoints_write_info(this, iunit, namespace, absolute_coordinates)
1412 class(kpoints_t), intent(in) :: this
1413 integer, optional, intent(in) :: iunit
1414 type(namespace_t), optional, intent(in) :: namespace
1415 logical, optional, intent(in) :: absolute_coordinates
1416
1417 integer :: ik, idir
1418 character(len=100) :: str_tmp
1419 character :: index
1420
1421 push_sub(kpoints_write_info)
1422
1423 call messages_print_with_emphasis('Brillouin zone sampling', iunit, namespace)
1424
1425 if (this%method == kpoints_monkh_pack) then
1426
1427 call messages_write('Dimensions of the k-point grid =')
1428 do idir = 1, this%full%dim
1429 call messages_write(this%nik_axis(idir), fmt = '(i3,1x)')
1430 end do
1431 call messages_new_line()
1432
1433 if (.not. all(this%downsampling(1:this%full%dim) == 1)) then
1434 call messages_write('Dimensions of the q-point grid =')
1435 do idir = 1, this%full%dim
1436 call messages_write(this%niq_axis(idir), fmt = '(i3,1x)')
1437 end do
1438 call messages_new_line()
1439 end if
1440
1441 call messages_write('Total number of k-points =')
1442 call messages_write(this%full%npoints)
1443 call messages_new_line()
1444
1445 call messages_write('Number of symmetry-reduced k-points =')
1446 call messages_write(this%reduced%npoints)
1447
1448 call messages_info(iunit=iunit, namespace=namespace)
1449
1450 else
1451
1452 call messages_write('Total number of k-points =')
1453 call messages_write(this%full%npoints)
1454 call messages_new_line()
1455 call messages_info(iunit=iunit, namespace=namespace)
1456
1457 end if
1458
1459 call messages_new_line()
1460 call messages_write('List of k-points:')
1461 call messages_info(iunit=iunit, namespace=namespace)
1462
1463 write(message(1), '(6x,a)') 'ik'
1464 do idir = 1, this%full%dim
1465 index = index2axis(idir)
1466 write(str_tmp, '(9x,2a)') 'k_', index
1467 message(1) = trim(message(1)) // trim(str_tmp)
1468 end do
1469 write(str_tmp, '(6x,a)') 'Weight'
1470 message(1) = trim(message(1)) // trim(str_tmp)
1471 message(2) = '---------------------------------------------------------'
1472 call messages_info(2, iunit)
1473
1474 do ik = 1, kpoints_number(this)
1475 write(message(1),'(i8,1x)') ik
1476 do idir = 1, this%full%dim
1477 if (optional_default(absolute_coordinates, .false.)) then
1478 write(str_tmp,'(f12.4)') this%reduced%point(idir, ik)
1479 else
1480 write(str_tmp,'(f12.4)') this%reduced%red_point(idir, ik)
1481 end if
1482 message(1) = trim(message(1)) // trim(str_tmp)
1483 end do
1484 write(str_tmp,'(f12.4)') kpoints_get_weight(this, ik)
1485 message(1) = trim(message(1)) // trim(str_tmp)
1486 call messages_info(1, iunit, namespace=namespace)
1487 end do
1488
1489 call messages_info(iunit=iunit, namespace=namespace)
1490
1491 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
1492
1493 pop_sub(kpoints_write_info)
1494 end subroutine kpoints_write_info
1495
1496
1497 ! ---------------------------------------------------------
1498 logical pure function kpoints_point_is_gamma(this, ik) result(is_gamma)
1499 class(kpoints_t), intent(in) :: this
1500 integer, intent(in) :: ik
1501
1502 is_gamma = (maxval(abs(kpoints_get_point(this, ik))) < m_epsilon)
1503
1504 end function kpoints_point_is_gamma
1505
1506 !--------------------------------------------------------
1507
1508 integer pure function kpoints_get_num_symmetry_ops(this, ik) result(num)
1509 type(kpoints_t), intent(in) :: this
1510 integer, intent(in) :: ik
1511
1512 if (this%use_symmetries) then
1513 num = this%num_symmetry_ops(ik)
1514 else
1515 num = 1
1516 end if
1517
1518 end function kpoints_get_num_symmetry_ops
1519 !--------------------------------------------------------
1520
1521 integer pure function kpoints_get_symmetry_ops(this, ik, index) result(iop)
1522 type(kpoints_t), intent(in) :: this
1523 integer, intent(in) :: ik
1524 integer, intent(in) :: index
1525
1526 if (this%use_symmetries) then
1527 iop = this%symmetry_ops(ik, index)
1528 else
1529 iop = 1
1530 end if
1531
1532 end function kpoints_get_symmetry_ops
1533
1534 !--------------------------------------------------------
1535 logical pure function kpoints_is_valid_symmetry(this, ik, index) result(valid)
1536 type(kpoints_t), intent(in) :: this
1537 integer, intent(in) :: ik
1538 integer, intent(in) :: index
1539
1540 integer :: iop
1541
1542 valid = .false.
1543 if (this%use_symmetries) then
1544 do iop = 1, this%num_symmetry_ops(ik)
1545 if (this%symmetry_ops(ik, iop) == index) then
1546 valid = .true.
1547 return
1548 end if
1549 end do
1550 else
1551 valid = .true.
1552 end if
1553
1554 end function kpoints_is_valid_symmetry
1555
1556
1557 !--------------------------------------------------------
1558 integer function kpoints_kweight_denominator(this)
1559 type(kpoints_t), intent(in) :: this
1560
1561 integer :: denom, nik, nik_skip
1562
1564
1565 nik = this%full%npoints
1566 nik_skip = this%nik_skip
1567
1568 if (this%method == kpoints_monkh_pack) then
1569 kpoints_kweight_denominator = this%full%npoints
1570 else
1572 ! NB largest reasonable value is: # k-points x 48. from space-group symmetries
1573 do denom = 1, 100000
1574 if (all(abs(int(this%full%weight(1:nik-nik_skip)*denom + 10_real64*m_epsilon) - &
1575 this%full%weight(1:nik-nik_skip)*denom) < 100_real64*m_epsilon)) then
1577 exit
1578 end if
1579 end do
1580 end if
1581
1583 end function kpoints_kweight_denominator
1584
1585 !--------------------------------------------------------
1586 logical pure function kpoints_have_zero_weight_path(this) result(have_zerow)
1587 class(kpoints_t), intent(in) :: this
1588
1589 if (this%nik_skip > 0) then
1590 have_zerow = .true.
1591 else
1592 have_zerow = .false.
1593 end if
1594
1596
1597 !--------------------------------------------------------
1598 integer pure function kpoints_get_kpoint_method(this)
1599 class(kpoints_t), intent(in) :: this
1600
1601 kpoints_get_kpoint_method = this%method
1602 end function kpoints_get_kpoint_method
1604 !--------------------------------------------------------
1605 real(real64) pure function kpoints_get_path_coord(this, ind) result(coord)
1606 type(kpoints_t), intent(in) :: this
1607 integer, intent(in) :: ind
1608
1609 coord = this%coord_along_path(ind)
1610 end function
1611
1612
1613
1614 !--------------------------------------------------------
1615 subroutine kpoints_check_symmetries(grid, symm, dim, time_reversal, namespace)
1616 type(kpoints_grid_t), intent(in) :: grid
1617 type(symmetries_t), intent(in) :: symm
1618 integer, intent(in) :: dim
1619 logical, intent(in) :: time_reversal
1620 type(namespace_t), intent(in) :: namespace
1621
1622 integer, allocatable :: kmap(:)
1623 real(real64) :: kpt(dim), diff(dim)
1624 integer :: nk, ik, ik2, iop
1625 type(distributed_t) :: kpt_dist
1626
1627 push_sub(kpoints_check_symmetries)
1628
1629 nk = grid%npoints
1631 !We distribute the k-points here for this routine, independently of the rest of the code
1632 call distributed_nullify(kpt_dist, nk)
1633#ifdef HAVE_MPI
1634 if (mpi_world%comm /= mpi_comm_undefined) then
1635 call distributed_init(kpt_dist, nk, mpi_comm_world, "kpt_check")
1636 end if
1637#endif
1638
1639 !A simple map to tell if the k-point as a matching symmetric point or not
1640 safe_allocate(kmap(kpt_dist%start:kpt_dist%end))
1641
1642 do iop = 1, symmetries_number(symm)
1643 if (iop == symmetries_identity_index(symm) .and. &
1644 .not. time_reversal) cycle
1645
1646 do ik = kpt_dist%start, kpt_dist%end
1647 kmap(ik) = ik
1648 end do
1649
1650 do ik = kpt_dist%start, kpt_dist%end
1651 !We apply the symmetry
1652 call symmetries_apply_kpoint_red(symm, iop, grid%red_point(:, ik), kpt)
1653 !We remove potential umklapp
1654 kpt = kpt - anint(kpt + m_half*symprec)
1655
1656 ! remove (mark) k-points which already have a symmetric point
1657 do ik2 = 1, nk
1658
1659 if (iop /= symmetries_identity_index(symm)) then
1660 diff = kpt - grid%red_point(:, ik2)
1661 diff = diff - anint(diff)
1662 !We found point corresponding to the symmetric kpoint
1663 if (sum(abs(diff)) < symprec) then
1664 kmap(ik) = -ik2
1665 exit
1666 end if
1667 end if
1668
1669 if (time_reversal) then
1670 diff = kpt + grid%red_point(:, ik2)
1671 diff = diff - anint(diff)
1672 !We found point corresponding to the symmetric kpoint
1673 if (sum(abs(diff)) < symprec) then
1674 kmap(ik) = -ik2
1675 exit
1676 end if
1677 end if
1678
1679 end do
1680 !In case we have not found a symnetric k-point...
1681 if (kmap(ik) == ik) then
1682 write(message(1),'(a,i5,a2,3(f7.3,a2),a)') "The reduced k-point ", ik, " (", &
1683 grid%red_point(1, ik), ", ", grid%red_point(2, ik), ", ", grid%red_point(3, ik), &
1684 ") ", "has no symmetric in the k-point grid for the following symmetry"
1685 write(message(2),'(i5,1x,a,2x,3(3i4,2x))') iop, ':', transpose(symm_op_rotation_matrix_red(symm%ops(iop)))
1686 message(3) = "Change your k-point grid or use KPointsUseSymmetries=no."
1687 call messages_fatal(3, namespace=namespace)
1688 end if
1689 end do
1690 end do
1691
1692 safe_deallocate_a(kmap)
1694 call distributed_end(kpt_dist)
1695
1697 end subroutine kpoints_check_symmetries
1698
1699 !--------------------------------------------------------
1700 logical function kpoints_is_compatible_downsampling(kpt, ik, iq) result(compatible)
1701 type(kpoints_t), intent(in) :: kpt
1702 integer, intent(in) :: ik
1703 integer, intent(in) :: iq
1704
1705 integer :: dim, idim
1706 real(real64) :: diff(kpt%reduced%dim), red(kpt%reduced%dim)
1707
1709
1710 compatible = .true.
1711 dim = kpt%reduced%dim
1712
1713 !No downsampling. We use all k-points
1714 if (all(kpt%downsampling == 1)) then
1716 return
1717 end if
1718
1719 assert(kpt%method == kpoints_monkh_pack)
1720
1721 diff = kpt%reduced%red_point(:, ik) - kpt%reduced%red_point(:, iq)
1722 do idim = 1, dim
1723 !We remove potential umklapp
1724 diff(idim) = diff(idim) - anint(diff(idim) + m_half*symprec)
1725 red(idim) = diff(idim)*kpt%nik_axis(idim)/real(kpt%downsampling(idim))
1726 if (abs(red(idim) - anint(red(idim))) > m_epsilon) then
1727 compatible = .false.
1729 return
1730 end if
1731 end do
1732
1735
1736 ! ---------------------------------------------------------
1737 integer function kpoints_nkpt_in_path(this) result(npath)
1738 class(kpoints_t), intent(in) :: this
1739
1740 npath = 0
1741 if (allocated(this%coord_along_path)) then
1742 npath = SIZE(this%coord_along_path)
1743 end if
1744
1745 end function kpoints_nkpt_in_path
1746
1747 ! ---------------------------------------------------------
1748 logical function kpoints_gamma_only(this) result(gamma_only)
1749 class(kpoints_t), intent(in) :: this
1750
1751 push_sub(kpoints_gamma_only)
1752
1753 gamma_only = kpoints_number(this) == 1 .and. kpoints_point_is_gamma(this, 1)
1754
1755 pop_sub(kpoints_gamma_only)
1756 end function kpoints_gamma_only
1757
1758 ! ---------------------------------------------------------
1759 subroutine kpoints_lattice_vectors_update(this, new_latt)
1760 class(kpoints_t), intent(inout) :: this
1761 type(lattice_vectors_t), intent(in) :: new_latt
1762
1763 integer :: ik
1764
1766
1767 this%latt = new_latt
1768
1769 do ik = 1, this%reduced%npoints
1770 call kpoints_to_absolute(this%latt, this%reduced%red_point(:, ik), this%reduced%point(:, ik))
1771 end do
1772
1773 do ik = 1, this%full%npoints
1774 call kpoints_to_absolute(this%latt, this%full%red_point(:, ik), this%full%point(:, ik))
1775 end do
1776
1778 end subroutine kpoints_lattice_vectors_update
1779
1780end module kpoints_oct_m
1781
1782!! Local Variables:
1783!! mode: f90
1784!! coding: utf-8
1785!! 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:758
subroutine read_user_kpoints()
Definition: kpoints.F90:906
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:1694
subroutine, public kpoints_grid_generate(dim, naxis, nshifts, shift, kpoints, lk123)
Generates the k-points grid.
Definition: kpoints.F90:1182
subroutine, public kpoints_end(this)
Definition: kpoints.F90:1042
logical pure function kpoints_have_zero_weight_path(this)
Definition: kpoints.F90:1682
subroutine kpoints_check_symmetries(grid, symm, dim, time_reversal, namespace)
Definition: kpoints.F90:1711
integer, parameter, public kpoints_monkh_pack
Definition: kpoints.F90:215
integer pure function, public kpoints_get_num_symmetry_ops(this, ik)
Definition: kpoints.F90:1604
subroutine kpoints_grid_addto(this, that)
Definition: kpoints.F90:282
real(real64) pure function, public kpoints_get_path_coord(this, ind)
Definition: kpoints.F90:1701
integer pure function kpoints_get_equiv(this, ik)
Definition: kpoints.F90:1165
logical pure function, public kpoints_is_valid_symmetry(this, ik, index)
Definition: kpoints.F90:1631
integer function, public kpoints_kweight_denominator(this)
Definition: kpoints.F90:1654
subroutine, public kpoints_copy(kin, kout)
Definition: kpoints.F90:1089
subroutine, public kpoints_path_generate(dim, latt, nkpoints, nsegments, resolution, highsympoints, kpoints, coord)
Generate the k-point along a path.
Definition: kpoints.F90:1281
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:1139
integer function, public kpoints_nkpt_in_path(this)
Definition: kpoints.F90:1833
logical function, public kpoints_is_compatible_downsampling(kpt, ik, iq)
Definition: kpoints.F90:1796
real(real64) pure function kpoints_get_weight(this, ik)
Definition: kpoints.F90:1155
integer pure function, public kpoints_get_symmetry_ops(this, ik, index)
Definition: kpoints.F90:1617
subroutine kpoints_grid_copy(bb, aa)
Definition: kpoints.F90:263
subroutine, public kpoints_fold_to_1bz(grid, latt)
Definition: kpoints.F90:1464
subroutine kpoints_grid_reduce(symm, time_reversal, nkpoints, dim, kpoints, weights, equiv, symm_ops, num_symm_ops)
Definition: kpoints.F90:1343
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:1594
integer, parameter, public kpoints_user
Definition: kpoints.F90:215
logical function kpoints_gamma_only(this)
Definition: kpoints.F90:1844
subroutine, public kpoints_lattice_vectors_update(this, new_latt)
Definition: kpoints.F90:1855
subroutine, public kpoints_to_reduced(latt, kin, kout)
Definition: kpoints.F90:1075
integer pure function, public kpoints_number(this)
Definition: kpoints.F90:1130
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1062
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:1507
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_experimental(name, namespace)
Definition: messages.F90:1039
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.
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)