TABLE OF CONTENTS
- ABINIT/m_kptrank
- m_kptrank/copy_kptrank
- m_kptrank/destroy_kptrank
- m_kptrank/dump_kptrank
- m_kptrank/get_rank_1kpt
- m_kptrank/kptrank_index
- m_kptrank/kptrank_type
- m_kptrank/mkkptrank
ABINIT/m_kptrank [ Modules ]
NAME
m_kptrank
FUNCTION
This module deals with rank objects for hashing k-point vector lists
COPYRIGHT
Copyright (C) 2010-2018 ABINIT group (MVer) This file is distributed under the terms of the GNU General Public Licence, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
PARENTS
SOURCE
21 #if defined HAVE_CONFIG_H 22 #include "config.h" 23 #endif 24 25 #include "libtetra.h" 26 27 module m_kptrank 28 29 USE_MEMORY_PROFILING 30 USE_MSG_HANDLING 31 32 implicit none 33 34 private 35 36 double precision :: zero = 0.0d0, half = 0.5d0, one = 1.0d0, tol8 = 1.d-8, tol10 = 1.d-10, tol12 = 1.d-12
m_kptrank/copy_kptrank [ Functions ]
[ Top ] [ m_kptrank ] [ Functions ]
NAME
copy_kptrank
FUNCTION
Copy the object
INPUTS
OUTPUT
krank = object containing ranking and inverse ranking, to be deallocated
PARENTS
defs_elphon,elphon
CHILDREN
SOURCE
368 subroutine copy_kptrank (krank_in, krank_out) 369 370 371 !This section has been created automatically by the script Abilint (TD). 372 !Do not modify the following lines by hand. 373 #undef ABI_FUNC 374 #define ABI_FUNC 'copy_kptrank' 375 !End of the abilint section 376 377 implicit none 378 379 !Arguments ------------------------------------ 380 !scalars 381 type(kptrank_type), intent(in) :: krank_in 382 type(kptrank_type), intent(out) :: krank_out 383 384 ! ********************************************************************* 385 krank_out%max_linear_density = krank_in%max_linear_density 386 krank_out%max_rank = krank_in%max_rank 387 krank_out%npoints = krank_in%npoints 388 389 TETRA_ALLOCATE(krank_out%rank, (krank_out%npoints)) 390 krank_out%rank = krank_in%rank 391 392 TETRA_ALLOCATE(krank_out%invrank, (krank_out%max_rank)) 393 krank_out%invrank = krank_in%invrank 394 395 TETRA_ALLOCATE(krank_out%multipl, (krank_out%npoints)) 396 krank_out%multipl = krank_in%multipl 397 398 end subroutine copy_kptrank
m_kptrank/destroy_kptrank [ Functions ]
[ Top ] [ m_kptrank ] [ Functions ]
NAME
destroy_kptrank
FUNCTION
This routine deallocates the arrays in a kptrank_type structure
INPUTS
krank = object containing ranking and inverse ranking, to be deallocated
PARENTS
defs_elphon,get_full_kgrid,m_ddk,m_ebands,m_fstab,m_nesting,m_phgamma m_pptools,m_tetrahedron,mkfskgrid,mkqptequiv,order_fs_kpts,outelph read_el_veloc
CHILDREN
SOURCE
423 subroutine destroy_kptrank (krank) 424 425 426 !This section has been created automatically by the script Abilint (TD). 427 !Do not modify the following lines by hand. 428 #undef ABI_FUNC 429 #define ABI_FUNC 'destroy_kptrank' 430 !End of the abilint section 431 432 implicit none 433 434 !Arguments ------------------------------------ 435 !scalars 436 type(kptrank_type), intent(inout) :: krank 437 438 ! ********************************************************************* 439 440 if (allocated(krank%rank)) then 441 TETRA_DEALLOCATE(krank%rank) 442 end if 443 if (allocated(krank%invrank)) then 444 TETRA_DEALLOCATE(krank%invrank) 445 end if 446 if (allocated(krank%multipl)) then 447 TETRA_DEALLOCATE(krank%multipl) 448 end if 449 450 end subroutine destroy_kptrank
m_kptrank/dump_kptrank [ Functions ]
[ Top ] [ m_kptrank ] [ Functions ]
NAME
dump_kptrank
FUNCTION
This routine prints the arrays and dimensions of a kptrank_type structure
INPUTS
krank = object containing ranking and inverse ranking unout = unit for open file to print to
PARENTS
CHILDREN
SOURCE
473 subroutine dump_kptrank (krank, unout) 474 475 476 !This section has been created automatically by the script Abilint (TD). 477 !Do not modify the following lines by hand. 478 #undef ABI_FUNC 479 #define ABI_FUNC 'dump_kptrank' 480 !End of the abilint section 481 482 implicit none 483 484 !Arguments ------------------------------------ 485 !scalars 486 integer, intent(in) :: unout 487 !arrays 488 type(kptrank_type), intent(in) :: krank 489 490 ! ********************************************************************* 491 492 write(unout, *) 493 write(unout, '(a)') ' Dump of the contents of a kptrank_type structure with k-point rank information' 494 write(unout, '(a,I8)') ' max linear density of points in 3 directions: max_linear_density = ', krank%max_linear_density 495 write(unout, '(a,I8)') ' maximum rank for any point in grid: max_rank = ', krank%max_rank 496 write(unout, '(a,I8)') ' number of points in input grid: npoints = ', krank%npoints 497 write(unout, *) 498 write(unout, '(a)') ' invrank array = ' 499 write(unout, '(I4)') krank%invrank(:) 500 write(unout, '(a)') ' rank array = ' 501 write(unout, '(I4)') krank%rank(:) 502 write(unout, '(a)') ' multiplicity array = ' 503 write(unout, '(I4)') krank%multipl(:) 504 write(unout, *) 505 506 end subroutine dump_kptrank
m_kptrank/get_rank_1kpt [ Functions ]
[ Top ] [ m_kptrank ] [ Functions ]
NAME
get_rank_1kpt
FUNCTION
This routine calculates the rank for one kpt
INPUTS
kpt = coordinates of kpoints krank = rank object for the k-grid we are using
OUTPUT
rank = rank of the kpoint
PARENTS
elphon,get_full_kgrid,integrate_gamma,integrate_gamma_alt,k_neighbors m_ddk,m_kptrank,m_nesting,m_pptools,m_tetrahedron,mkfskgrid,mkqptequiv read_el_veloc,read_gkk
CHILDREN
SOURCE
229 subroutine get_rank_1kpt(kpt,rank,krank) 230 231 232 !This section has been created automatically by the script Abilint (TD). 233 !Do not modify the following lines by hand. 234 #undef ABI_FUNC 235 #define ABI_FUNC 'get_rank_1kpt' 236 !End of the abilint section 237 238 implicit none 239 240 !Arguments ------------------------------------ 241 !scalars 242 integer,intent(out) :: rank 243 type(kptrank_type), intent(in) :: krank 244 !arrays 245 double precision,intent(in) :: kpt(3) 246 247 !Local variables------------------------------- 248 !scalars 249 character(len=500) :: msg 250 !arrays 251 double precision :: redkpt(3) 252 253 ! ************************************************************************* 254 255 ! wrap to [0, 1[ -> replaced call to wrap2_zeroone inline, to encapsulate this module 256 if (kpt(1)>zero) then 257 redkpt(1)=mod((kpt(1)+tol12),one)-tol12 258 else 259 redkpt(1)=-mod(-(kpt(1)-one+tol12),one)+one-tol12 260 end if 261 if(abs(redkpt(1))<tol12)redkpt(1)=zero 262 263 if (kpt(2)>zero) then 264 redkpt(2)=mod((kpt(2)+tol12),one)-tol12 265 else 266 redkpt(2)=-mod(-(kpt(2)-one+tol12),one)+one-tol12 267 end if 268 if(abs(redkpt(2))<tol12)redkpt(2)=zero 269 270 if (kpt(3)>zero) then 271 redkpt(3)=mod((kpt(3)+tol12),one)-tol12 272 else 273 redkpt(3)=-mod(-(kpt(3)-one+tol12),one)+one-tol12 274 end if 275 if(abs(redkpt(3))<tol12)redkpt(3)=zero 276 277 278 279 ! rank = int(real(krank%max_linear_density)*(redkpt(3)+half+tol8 +& 280 !& real(krank%max_linear_density)*(redkpt(2)+half+tol8 +& 281 !& real(krank%max_linear_density)*(redkpt(1)+half+tol8)))) 282 rank = int(real(krank%max_linear_density)*(redkpt(1)+half+tol8 +& 283 & real(krank%max_linear_density)*(redkpt(2)+half+tol8 +& 284 & real(krank%max_linear_density)*(redkpt(3)+half+tol8)))) 285 286 if (rank > krank%max_rank) then 287 write(msg,'(a,i0)') ' rank should be inferior to ', krank%max_rank 288 TETRA_ERROR(msg) 289 end if 290 291 end subroutine get_rank_1kpt
m_kptrank/kptrank_index [ Functions ]
[ Top ] [ m_kptrank ] [ Functions ]
NAME
kptrank_index
FUNCTION
Return the index of the k-point `kpt` in the initial set. -1 if not found.
INPUTS
krank = rank object for the k-grid we are using kpt = coordinates of kpoints
OUTPUT
PARENTS
CHILDREN
SOURCE
316 integer function kptrank_index(krank, kpt) result(ikpt) 317 318 319 !This section has been created automatically by the script Abilint (TD). 320 !Do not modify the following lines by hand. 321 #undef ABI_FUNC 322 #define ABI_FUNC 'kptrank_index' 323 !End of the abilint section 324 325 implicit none 326 327 !Arguments ------------------------------------ 328 !scalars 329 type(kptrank_type), intent(in) :: krank 330 !arrays 331 double precision,intent(in) :: kpt(3) 332 333 !Local variables------------------------------- 334 !scalars 335 integer :: kpt_rank 336 337 ! ************************************************************************* 338 339 call get_rank_1kpt(kpt, kpt_rank, krank) 340 ikpt = -1 341 if (kpt_rank < krank%max_rank) ikpt = krank%invrank(kpt_rank) 342 343 end function kptrank_index
m_kptrank/kptrank_type [ Types ]
[ Top ] [ m_kptrank ] [ Types ]
NAME
kptrank_type
FUNCTION
structure to contain a rank/inverse rank pair of arrays, with dimensions
SOURCE
48 type,public :: kptrank_type 49 integer :: max_linear_density 50 integer :: max_rank 51 integer :: npoints 52 logical :: time_reversal 53 integer,allocatable :: invrank(:) 54 integer,allocatable :: rank(:) 55 integer,allocatable :: multipl(:) 56 end type kptrank_type 57 58 public :: mkkptrank ! Sets up the kpt ranks for comparing kpts 59 public :: get_rank_1kpt ! Calculates the rank for one kpt 60 public :: kptrank_index ! Return the index of the k-point `kpt` in the initial set. -1 if not found. 61 public :: copy_kptrank ! Copy the object 62 public :: destroy_kptrank ! Free memory 63 public :: dump_kptrank ! Prints the arrays and dimensions of a kptrank_type structure
m_kptrank/mkkptrank [ Functions ]
[ Top ] [ m_kptrank ] [ Functions ]
NAME
mkkptrank
FUNCTION
This routine sets up the kpt ranks for comparing kpts
INPUTS
npt = number of kpoints (eventually irreducible) kpt = coordinates of kpoints time_reversal = true or false to use time reversal symmetry. Default is true, but only important if nsym and symrec are present
OUTPUT
krank = object containing ranking and inverse ranking
PARENTS
get_full_kgrid,m_ddk,m_ebands,m_fstab,m_nesting,m_phgamma,m_pptools m_tetrahedron,mkfskgrid,mkqptequiv,order_fs_kpts,outelph,read_el_veloc
CHILDREN
SOURCE
94 subroutine mkkptrank (kpt,nkpt,krank,nsym,symrec, time_reversal) 95 96 97 !This section has been created automatically by the script Abilint (TD). 98 !Do not modify the following lines by hand. 99 #undef ABI_FUNC 100 #define ABI_FUNC 'mkkptrank' 101 !End of the abilint section 102 103 implicit none 104 105 !Arguments ------------------------------------ 106 !scalars 107 integer,intent(in) :: nkpt 108 integer,intent(in), optional :: nsym 109 logical,intent(in), optional :: time_reversal 110 !arrays 111 type(kptrank_type), intent(out) :: krank 112 double precision,intent(in) :: kpt(3,nkpt) 113 integer,intent(in), optional :: symrec(3,3, *) 114 115 !Local variables ------------------------- 116 !scalars 117 integer :: ikpt, isym, symkptrank, irank 118 integer :: timrev, itim 119 double precision :: smallestlen 120 character(len=500) :: msg 121 !arrays 122 double precision :: symkpt(3) 123 124 ! ********************************************************************* 125 126 ! find smallest linear length 127 smallestlen = one 128 do ikpt=1, nkpt 129 if (abs(kpt(1,ikpt)) > tol10) smallestlen = min(smallestlen, abs(kpt(1,ikpt))) 130 if (abs(kpt(2,ikpt)) > tol10) smallestlen = min(smallestlen, abs(kpt(2,ikpt))) 131 if (abs(kpt(3,ikpt)) > tol10) smallestlen = min(smallestlen, abs(kpt(3,ikpt))) 132 end do 133 134 krank%max_linear_density = int(one/smallestlen)+1 135 krank%max_rank = 2*krank%max_linear_density**3 136 krank%npoints = nkpt 137 138 TETRA_ALLOCATE(krank%rank, (nkpt)) 139 140 TETRA_ALLOCATE(krank%invrank, (krank%max_rank)) 141 krank%invrank(:) = -1 142 143 timrev = 2 144 krank%time_reversal = .true. 145 if (present(time_reversal)) then 146 if (.not. time_reversal) then 147 timrev = 1 148 end if 149 end if 150 151 !Ensure kpt(i)+one is positive, and the smallest 152 !difference between kpts should be larger than 1/100 153 !ie ngkpt < 100. 154 ! the following fills invrank for the k-points in the list provided (may be only the irred kpts) 155 do ikpt=1,nkpt 156 call get_rank_1kpt (kpt(:,ikpt), krank%rank(ikpt), krank) 157 158 if (krank%rank(ikpt) > krank%max_rank .or. krank%rank(ikpt) < 1) then 159 write(msg,'(a,2i0)')" max rank exceeded or < 1, ikpt, rank ", ikpt, krank%rank(ikpt) 160 TETRA_ERROR(msg) 161 end if 162 krank%invrank(krank%rank(ikpt)) = ikpt 163 end do 164 165 ! if symrec is provided, fill invrank with appropriate irred kpt indices 166 ! for symmetry completion: kptrank_t%invrank points to the irred k-point 167 ! equivalent to the k-point whose rank is provided 168 if (present(symrec)) then 169 if(.not. present(nsym)) then 170 write (msg,'(a)') "need both symrec and nsym arguments together" 171 TETRA_ERROR(msg) 172 end if 173 do ikpt=1,nkpt 174 ! itim == 1 for positive, and itim==2 gives Kramers opposite of k-point 175 ! favor the former by looping it last 176 do itim = timrev, 1, -1 177 do isym = 1, nsym 178 symkpt = (-1)**(timrev+1) * matmul(symrec(:,:,isym), kpt(:, ikpt)) 179 180 call get_rank_1kpt (symkpt(:), symkptrank, krank) 181 182 krank%invrank(symkptrank) = ikpt 183 end do 184 end do 185 end do 186 end if 187 188 TETRA_ALLOCATE(krank%multipl, (nkpt)) 189 190 krank%multipl = 0 191 ! find multiplicity of ikpt 192 do irank = 1, krank%max_rank 193 ikpt = krank%invrank(irank) 194 if (ikpt > 0) krank%multipl(ikpt) = krank%multipl(ikpt) + 1 195 end do 196 197 ! could add a check that all k in full grid are in the invrank... 198 199 200 end subroutine mkkptrank