TABLE OF CONTENTS
- ABINIT/m_supercell
- m_effective_potential/getPBCIndexes_supercell
- m_supercell/distance_supercell
- m_supercell/freeze_displ_supercell
- m_supercell/init_supercell
- m_supercell/init_supercell_for_qpt
- m_supercell/mksupercell
- m_supercell/order_supercell_typat
- m_supercell/supercell_type
ABINIT/m_supercell [ Modules ]
NAME
m_supercell
FUNCTION
Module for using a supercell, in particular for phonon displacement freezing. Container type is defined, and destruction, print subroutines as well as the central init_supercell
COPYRIGHT
Copyright (C) 2010-2024 ABINIT group (MJV, DJA) 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 .
SOURCE
20 #if defined HAVE_CONFIG_H 21 #include "config.h" 22 #endif 23 24 #include "abi_common.h" 25 26 module m_supercell 27 28 use defs_basis 29 use m_errors 30 use m_abicore 31 32 use m_symtk, only : matr3inv 33 use m_copy, only : alloc_copy 34 use m_io_tools, only : open_file 35 use m_fstrings, only : int2char4, write_num 36 37 implicit none 38 39 private
m_effective_potential/getPBCIndexes_supercell [ Functions ]
[ Top ] [ m_effective_potential ] [ Functions ]
NAME
FUNCTION
Get the index of the cell by using PBC
INPUTS
index = index of the cell into the supercell ncell = number of total cell
OUTPUT
index = index of the cell into the supercell with PBC
SOURCE
612 subroutine getPBCIndexes_supercell(index,ncell) 613 614 implicit none 615 616 !Arguments --------------------------------------------- 617 integer, intent(inout) :: index(3) 618 integer, intent(in) :: ncell(3) 619 !Local variables --------------------------------------- 620 integer :: ii 621 ! ********************************************************************* 622 623 do ii=1,3 624 do while (index(ii) > ncell(ii)) 625 index(ii) = index(ii) - ncell(ii) 626 end do 627 do while (index(ii) <= 0) 628 index(ii) = index(ii) + ncell(ii) 629 end do 630 end do 631 632 end subroutine getPBCIndexes_supercell
m_supercell/distance_supercell [ Functions ]
[ Top ] [ m_supercell ] [ Functions ]
NAME
FUNCTION
compute the distance_supercell betwen 2 atoms in different cell
INPUTS
xcart1(3) = cartesian coordinates of the first atom xcart1(3) = cartesian coordinates of the second atom rprimd(3,3) = primitive lattice vectors cell1(3) = index of the cell of the first atom (for example -1 0 2) cell2(3) = index of the cell of the second atom (for example 0 0 2)
OUTPUT
distance_supercell = distance_supercell between the 2 atoms
SOURCE
690 function distance_supercell(xcart1,xcart2,rprimd,cell1,cell2) result(dist) 691 692 !Arguments ------------------------------------ 693 !scalar 694 !array 695 real(dp),intent(in):: rprimd(3,3) 696 real(dp),intent(in):: xcart1(3),xcart2(3) 697 integer,intent(in) :: cell1(3),cell2(3) 698 real(dp) :: dist 699 !Local variables ------------------------------- 700 real(dp) :: rpt1(3),rpt2(3) 701 integer :: mu 702 !! ************************************************************************* 703 do mu=1,3 704 rpt1(mu) = cell1(1)*rprimd(mu,1)+cell1(2)*rprimd(mu,2)+cell1(3)*rprimd(mu,3) 705 rpt2(mu) = cell2(1)*rprimd(mu,1)+cell2(2)*rprimd(mu,2)+cell2(3)*rprimd(mu,3) 706 end do 707 708 dist = ((xcart2(1)+rpt2(1)-xcart1(1)-rpt1(1))**2+& 709 & (xcart2(2)+rpt2(2)-xcart1(2)-rpt1(2))**2+& 710 & (xcart2(3)+rpt2(3)-xcart1(3)-rpt1(3))**2)**0.5 711 712 end function distance_supercell
m_supercell/freeze_displ_supercell [ Functions ]
[ Top ] [ m_supercell ] [ Functions ]
NAME
freeze_displ_supercell
FUNCTION
Freeze a specific displacement phonon field into the supercell scell
INPUTS
displ = phonon displacement vectors for this mode freeze_displ = desired amplitude for phonon displacement along displ. for thermal displacement use sqrt[ (1/2 + bose_einstein(freq,T)) / freq ] scell = supercell structure with reference atomic positions etc...
OUTPUT
scell = supercell structure: xcart will be updated with phonon displacement
SOURCE
347 subroutine freeze_displ_supercell (displ,freeze_displ,scell) 348 349 implicit none 350 351 !Arguments ------------------------------------ 352 !scalars 353 type(supercell_type), intent(inout) :: scell 354 real(dp), intent(in) :: freeze_displ 355 !arrays 356 real(dp), intent(in) :: displ(2,3*scell%natom_primcell) 357 358 !Local variables------------------------------- 359 !scalar 360 integer :: iatom, ipratom 361 complex(dpc) :: expqdotr, j=cmplx(zero,one) 362 complex(dpc) :: phase 363 complex(dpc) :: zdispl(3,scell%natom_primcell) 364 ! ************************************************************************* 365 366 zdispl = (cmplx(reshape(displ(1,:), (/3,scell%natom_primcell/)),& 367 & reshape(displ(2,:), (/3,scell%natom_primcell/)))) 368 369 ! fix gauge by imposing real displacement for first atom in first direction 370 ! multiply by normalized complex conjugate of first element 371 ! NB 6 March 2018: this may be imposing a positive (not just real) displacement for 1st atom along x!!! 372 ! That might be problematic below, though for the thermal displacement method freeze_displ swaps sign for each new mode 373 phase = cmplx(one,zero) 374 if (abs(zdispl(1,1)) > tol10) then 375 phase = conjg(zdispl(1,1)) / abs(zdispl(1,1)) 376 end if 377 378 do iatom = 1, scell%natom 379 expqdotr = exp(j*two_pi*(scell%qphon(1)*scell%uc_indexing(1,iatom) & 380 & +scell%qphon(2)*scell%uc_indexing(2,iatom) & 381 & +scell%qphon(3)*scell%uc_indexing(3,iatom))) 382 383 ! this is offset in zdispl vector due to primitive cell atom position 384 ipratom = scell%atom_indexing(iatom) 385 386 !add real part of displacement times Bloch phase 387 scell%xcart(:,iatom) = scell%xcart(:,iatom) & 388 & + freeze_displ * real(expqdotr * zdispl(:,ipratom) * phase) 389 390 ! scell%xcart(:,iatom) = scell%xcart(:,iatom) & 391 !& + freeze_displ * cos(qdotr) * displ(1,ipratom+1:ipratom+3) & 392 !& - freeze_displ * sin(qdotr) * displ(2,ipratom+1:ipratom+3) 393 end do 394 395 end subroutine freeze_displ_supercell
m_supercell/init_supercell [ Functions ]
[ Top ] [ m_supercell ] [ Functions ]
NAME
init_supercell
FUNCTION
Initialize scell structure, from unit cell vectors, and atoms, based on rlatt multiplicity matrix
INPUTS
natom_primcell = number of atoms in primitive cell rlatt(3,3) = multiplicity of primtive unit cells in supercell rprimd_primcell(3,3) = real space lattice vectors (bohr) typat_primcell(natom) = types of all atoms in primitive cell xcart_primcell(3,natom) = cartesian positions of atoms in primitive cell znucl = nuclear charges for all species ordering = if true, typat will be 1 1 1 1 1 1 2 2 2 2 2 2 3 3 3 3 .... if false, typat will be 1 2 3 4 1 2 3 4 1 2 3 4 1 2 3 4 ....
OUTPUT
scell = supercell structure to be initialized
SOURCE
197 subroutine init_supercell(natom_primcell, rlatt, rprimd_primcell, typat_primcell, xcart_primcell, znucl, scell, ordering) 198 199 implicit none 200 !Arguments ------------------------------------ 201 !scalars 202 integer, intent(in) :: natom_primcell 203 type(supercell_type), intent(out) :: scell 204 logical,optional,intent(in) :: ordering 205 !arrays 206 integer , intent(in) :: rlatt(3,3) 207 integer , intent(in) :: typat_primcell(natom_primcell) 208 real(dp), intent(in) :: znucl(:) 209 real(dp), intent(in) :: rprimd_primcell(3,3) 210 real(dp), intent(in) :: xcart_primcell(3,natom_primcell) 211 212 !local 213 !scalars 214 integer :: iatom_supercell, i1,i2,i3, iatom, icell 215 216 !arrays 217 218 scell%natom_primcell = natom_primcell 219 scell%rlatt = rlatt 220 scell%ncells = rlatt(1,1)*rlatt(2,2)*rlatt(3,3) 221 scell%rprimd(:,1) = rprimd_primcell(:,1) * rlatt(1,1) 222 scell%rprimd(:,2) = rprimd_primcell(:,2) * rlatt(2,2) 223 scell%rprimd(:,3) = rprimd_primcell(:,3) * rlatt(3,3) 224 225 scell%ntypat = size(znucl) 226 ABI_MALLOC(scell%znucl,(scell%ntypat)) 227 scell%znucl(:) = znucl(:) 228 229 !number of atoms in full supercell 230 scell%natom= natom_primcell*scell%ncells 231 ABI_MALLOC(scell%xcart,(3,scell%natom)) 232 ABI_MALLOC(scell%xcart_ref,(3,scell%natom)) 233 ABI_MALLOC(scell%typat,(scell%natom)) 234 ABI_MALLOC(scell%atom_indexing,(scell%natom)) 235 ABI_MALLOC(scell%uc_indexing,(3,scell%natom)) 236 ABI_MALLOC(scell%rvecs, (3, scell%ncells)) 237 238 iatom_supercell = 0 239 icell =0 240 do i1 = 1, rlatt(1,1) 241 do i2 = 1, rlatt(2,2) 242 do i3 = 1, rlatt(3,3) 243 244 icell=icell+1 245 scell%rvecs(:,icell)=(/i1-1, i2-1, i3-1/) 246 247 248 do iatom = 1, natom_primcell 249 iatom_supercell = iatom_supercell + 1 250 scell%uc_indexing(:,iatom_supercell) = (/i1-1,i2-1,i3-1/) 251 scell%xcart_ref(:,iatom_supercell) = xcart_primcell(:,iatom) & 252 & + matmul(rprimd_primcell,scell%uc_indexing(:,iatom_supercell)) 253 scell%atom_indexing(iatom_supercell) = iatom 254 scell%typat(iatom_supercell) = typat_primcell(iatom) 255 end do 256 end do 257 end do 258 end do 259 260 ABI_CHECK(iatom_supercell == scell%natom, "iatom_supercell /= scell%natom") 261 if(iatom_supercell /= scell%natom) then 262 write(std_out,*)"iatom_supercell /= scell%natom" 263 endif 264 265 266 267 scell%xcart = scell%xcart_ref 268 scell%qphon = zero 269 270 if (present(ordering)) then 271 if (ordering) then 272 call order_supercell_typat(scell) 273 end if 274 end if 275 276 end subroutine init_supercell
m_supercell/init_supercell_for_qpt [ Functions ]
[ Top ] [ m_supercell ] [ Functions ]
NAME
init_supercell_for_qpt
FUNCTION
Initialize scell structure, from unit cell vectors, and atoms, based on qpoint chosen
INPUTS
natom_primcell = number of atoms in primitive cell qphon(3) = phonon wavevector find smallest supercell which will accomodate phonon qphon = (1/2,1/2,1/2) rprimd_primcell(3,3) = real space lattice vectors (bohr) typat_primcell = types of atoms xcart_primcell(3,natom) = cartesian positions of atoms in primitive cell znucl = nuclear charges for all species ordering = if true, typat will be 1 1 1 1 1 1 2 2 2 2 2 2 3 3 3 3 .... if false, typat will be 1 2 3 4 1 2 3 4 1 2 3 4 1 2 3 4 ....
OUTPUT
scell = supercell structure to be initialized
SOURCE
111 subroutine init_supercell_for_qpt(natom_primcell, qphon, rprimd_primcell, & 112 & typat_primcell, xcart_primcell, znucl, scell, ordering) 113 114 implicit none 115 !Arguments ------------------------------------ 116 !scalars 117 integer, intent(in) :: natom_primcell 118 type(supercell_type), intent(out) :: scell 119 logical,optional,intent(in) :: ordering 120 !arrays 121 integer , intent(in) :: typat_primcell(natom_primcell) 122 real(dp), intent(in) :: qphon(3) 123 real(dp), intent(in) :: znucl(:) 124 real(dp), intent(in) :: rprimd_primcell(3,3) 125 real(dp), intent(in) :: xcart_primcell(3,natom_primcell) 126 127 !Local variables------------------------------- 128 !scalar 129 integer :: ii, maxsc, iscmult 130 real(dp) :: qbymult 131 !arrays 132 integer :: rlatt(3,3) ! number of primitive cells in each direction for the supercell 133 character(len=500) :: msg 134 ! ************************************************************************* 135 136 137 ! maximum number of unit cells in a given direction 138 maxsc = 10 139 140 ! find smallest supercell which will accomodate phonon. 141 ! FIXME: for the moment, just get smallest multiple along each direction, with an upper bound 142 rlatt = 0 143 rlatt(1,1) = -1 144 rlatt(2,2) = -1 145 rlatt(3,3) = -1 146 do ii=1,3 147 do iscmult=1,maxsc 148 qbymult = qphon(ii)*iscmult 149 if (abs(qbymult - int(qbymult)) < tol10) then 150 rlatt(ii,ii) = iscmult 151 exit 152 end if 153 end do 154 if (rlatt(ii,ii) == -1) then 155 write(msg,'(a,I4,a,I7,2a,3E20.10)')' No supercell found with less than ', & 156 & maxsc,' unit cells in direction ', & 157 & ii, ch10, ' qphon = ', qphon 158 ABI_ERROR(msg) 159 end if 160 end do 161 162 if (present(ordering)) then 163 call init_supercell(natom_primcell, rlatt, rprimd_primcell, typat_primcell, xcart_primcell, znucl, scell, ordering) 164 else 165 call init_supercell(natom_primcell, rlatt, rprimd_primcell, typat_primcell, xcart_primcell, znucl, scell) 166 end if 167 168 scell%qphon = qphon 169 170 end subroutine init_supercell_for_qpt
m_supercell/mksupercell [ Functions ]
[ Top ] [ m_supercell ] [ Functions ]
NAME
mksupercell
FUNCTION
computes atomic positons, magnetic ordering of supercell
INPUTS
magv_org (optional) magnetic ordering of atoms in primitive cell, ordering of atoms given als 1 and -1, if not given fm is assumed xred_org relative position of atoms in primitive cell rprimd_org unit cell dimensions of primitive cell natom=number of atoms in unit cell option= 1 output ion-ion distances / 2 output ordering of ion-ion distances / 3 output variables in varlist according to ion-ion distances * magnetic ordering
OUTPUT
magv_sc magnetic ordering of atoms in supercell xred_sc relative position of atoms in supercell rprimd_sc unit cell dimensions of supercell
SOURCE
788 subroutine mksupercell(xred_org,magv_org,rprimd_org,nat_org,nat_sc,xred_sc,magv_sc,rprimd_sc,ext,prtvol) 789 790 implicit none 791 792 !Arguments ------------------------------------ 793 !scalars 794 integer,intent(in) :: nat_org,nat_sc 795 integer,intent(in),optional :: prtvol 796 !arrays 797 real(dp),intent(in) :: rprimd_org(3,3) 798 integer,intent(in) :: ext(3) 799 real(dp),intent(in) :: xred_org(3,nat_org) 800 real(dp),intent(out) :: xred_sc(3,nat_sc) 801 real(dp),intent(out) :: magv_sc(nat_sc) 802 real(dp),intent(out) :: rprimd_sc(3,3) 803 integer,intent(in),optional :: magv_org(nat_org) 804 805 !Local variables------------------------------- 806 !scalars 807 integer :: prtvoll,ix,iy,iz,nprcl,iprcl,jdim,iatom 808 !arrays 809 real(dp) :: magvv_org(nat_org) 810 real(dp),allocatable :: transv(:,:,:) 811 812 ! ************************************************************************* 813 814 if (present(magv_org)) then 815 magvv_org=magv_org 816 else 817 magvv_org=(/ (1, iatom=1,nat_org) /) 818 end if 819 820 if (present(prtvol)) then 821 prtvoll=prtvol 822 else 823 prtvoll=1 824 end if 825 826 rprimd_sc=reshape((/ (rprimd_org(ix,:)*ext(ix) ,ix=1,3) /),(/3,3 /)) 827 nprcl=product(ext) 828 ABI_MALLOC(transv,(3,nat_org,nprcl)) 829 830 transv=reshape((/ (((((/ ix,iy,iz /),iatom=1,nat_org),ix=0,ext(1)-1),iy=0,ext(2)-1),iz=0,ext(3)-1) /), (/ 3, nat_org,nprcl/) ) 831 832 !DEBUG 833 !write(std_out,*)'mksupercell: xred_org ' ,xred_org 834 !END DEBUG 835 836 do iprcl=1,nprcl 837 xred_sc(:,1+(iprcl-1)*nat_org:iprcl*nat_org)=xred_org+transv(:,:,iprcl) 838 magv_sc(1+(iprcl-1)*nat_org:iprcl*nat_org)=magv_org 839 end do 840 841 842 do jdim=1,3 843 xred_sc(jdim,:)=xred_sc(jdim,:)/ext(jdim) 844 end do 845 846 !DEBUG 847 !write(std_out,*)'mksupercell: xred_sc ', xred_sc 848 !write(std_out,*)'mksupercell: magv_sc ', magv_sc 849 !END DEBUG 850 851 ABI_FREE(transv) 852 853 end subroutine mksupercell
m_supercell/order_supercell_typat [ Functions ]
[ Top ] [ m_supercell ] [ Functions ]
NAME
order_supercell_typat
FUNCTION
Re-order atoms in place for types
INPUTS
scell = supercell structure with reference atomic positions etc...
OUTPUT
scell = supercell structure: typat, xcart and so on will be updated
SOURCE
295 subroutine order_supercell_typat (scell) 296 297 implicit none 298 299 !Arguments ------------------------------------ 300 !scalars 301 type(supercell_type), intent(inout) :: scell 302 303 ! local tmp variables 304 integer :: itypat, iatom_supercell, iatom 305 type(supercell_type) :: scell_tmp 306 307 call copy_supercell (scell,scell_tmp) 308 309 iatom_supercell = 0 310 do itypat = 1, scell%ntypat 311 do iatom = 1, scell%natom 312 if (scell_tmp%typat(iatom) /= itypat) cycle 313 iatom_supercell = iatom_supercell + 1 314 scell%xcart(:,iatom_supercell) = scell_tmp%xcart(:,iatom) 315 scell%xcart_ref(:,iatom_supercell) = scell_tmp%xcart_ref(:,iatom) 316 scell%atom_indexing(iatom_supercell) = scell_tmp%atom_indexing(iatom) 317 scell%uc_indexing(:,iatom_supercell) = scell_tmp%uc_indexing(:,iatom) 318 scell%typat(iatom_supercell) = scell_tmp%typat(iatom) 319 end do 320 end do 321 322 call destroy_supercell(scell_tmp) 323 324 end subroutine order_supercell_typat
m_supercell/supercell_type [ Types ]
[ Top ] [ m_supercell ] [ Types ]
NAME
supercell_type
FUNCTION
structure for a supercell constructed from a basic rprimd and xcart, with indexing to original atoms The supercell may not be oriented the same way as the original cell, if you can reduce it by symmetry
SOURCE
53 type, public :: supercell_type 54 integer :: natom_primcell ! number of atoms in primitive cell 55 integer :: natom ! number of atoms in supercell 56 integer :: ntypat ! number of atom types 57 integer :: ncells ! number of unit cells in supercell 58 integer :: rlatt(3,3) ! matrix for multiplicity of supercell 59 real(dp) :: rprimd(3,3) ! new lattice vectors for supercell 60 real(dp) :: qphon(3) ! phonon q vector used to generate scell, if any 61 real(dp), allocatable :: xcart(:,:) ! (3, natom) positions of atoms 62 real(dp), allocatable :: xcart_ref(:,:) ! (3, natom) equilibrium positions of atoms 63 integer, allocatable :: atom_indexing(:) ! (natom) indexes original atom: 1..natom_primcell 64 integer, allocatable :: uc_indexing(:,:) ! (3, natom) indexes unit cell atom is in: 65 integer, allocatable :: typat(:) ! (3, natom) positions of atoms 66 real(dp), allocatable :: znucl(:) ! (ntypat) nuclear charges of species 67 integer, allocatable :: rvecs(:,:) ! supercell vectors 68 69 end type supercell_type 70 71 public :: init_supercell_for_qpt 72 public :: init_supercell 73 public :: freeze_displ_supercell 74 public :: prt_supercell_for_qpt 75 public :: prt_supercell 76 public :: copy_supercell 77 public :: distance_supercell 78 public :: findBound_supercell 79 public :: getPBCIndexes_supercell 80 public :: destroy_supercell 81 82 public :: mksupercell ! computes atomic positons, magnetic ordering of supercell