TABLE OF CONTENTS
- ABINIT/m_geometry
- defs_elphon/phdispl_cart2red
- m_geometry/normv_int_vector
- m_geometry/normv_int_vector_array
- m_geometry/normv_rdp_vector
- m_geometry/normv_rdp_vector_array
- m_geometry/vdotw_rc_vector
- m_geometry/vdotw_rr_vector
- m_geometry/wigner_seitz
ABINIT/m_geometry [ Modules ]
NAME
m_geometry
FUNCTION
This module contains basic tools to operate on vectors expressed in reduced coordinates.
COPYRIGHT
Copyright (C) 2008-2017 ABINIT group (MG) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
SOURCE
17 #if defined HAVE_CONFIG_H 18 #include "config.h" 19 #endif 20 21 #include "abi_common.h" 22 23 MODULE m_geometry 24 25 use defs_basis 26 use m_profiling_abi 27 use m_errors 28 29 implicit none 30 31 private 32 33 public :: normv ! Norm of vector(s) in reduced coordinates either in real or reciprocal space. 34 public :: vdotw ! Scalar product between two reduced vectors either in real or reciprocal space. 35 public :: wigner_seitz ! Find the grid of points falling inside the Wigner-Seitz cell. 36 public :: phdispl_cart2red ! Calculate the displacement vectors for all branches in reduced coordinates. 37 38 interface normv 39 module procedure normv_rdp_vector 40 module procedure normv_int_vector 41 !module procedure normv_int_vector_array ! WARNING for the time being, do not use these 2 procedures, 42 !module procedure normv_rdp_vector_array ! sunstudio12 is not able to resolve which sub should be called. 43 end interface normv 44 45 interface vdotw 46 module procedure vdotw_rr_vector 47 module procedure vdotw_rc_vector 48 end interface vdotw 49 50 CONTAINS !===========================================================
defs_elphon/phdispl_cart2red [ Functions ]
[ Top ] [ defs_elphon ] [ Functions ]
NAME
phdispl_cart2red
FUNCTION
Calculates the displacement vectors for all branches in reduced coordinates. $ displ_red = displ_cart \cdot gprimd $ for each phonon branch.
INPUTS
natom=Number of atoms. gprimd(3,3)Dimensional primitive translations for reciprocal space ($\textrm{bohr}^{-1}$) displ_cart(2,3*natom,3*natom)=Phonon displacement in Cartesian coordinates.
OUTPUT
displ_red(2,3*natom,3*natom)=Phonon displacement in reduded coordinates.
PARENTS
get_tau_k,m_ddb,m_ifc,mka2f,mkph_linwid,read_gkk
CHILDREN
SOURCE
622 subroutine phdispl_cart2red(natom,gprimd,displ_cart,displ_red) 623 624 625 !This section has been created automatically by the script Abilint (TD). 626 !Do not modify the following lines by hand. 627 #undef ABI_FUNC 628 #define ABI_FUNC 'phdispl_cart2red' 629 !End of the abilint section 630 631 implicit none 632 633 !Arguments ------------------------------------ 634 !scalars 635 integer,intent(in) :: natom 636 !arrays 637 real(dp),intent(in) :: gprimd(3,3) 638 real(dp),intent(in) :: displ_cart(2,3*natom,3*natom) 639 real(dp),intent(out) :: displ_red(2,3*natom,3*natom) 640 641 !Local variables------------------------- 642 !scalars 643 integer :: nbranch,jbranch,iatom,idir,ibranch,kdir,k1 644 645 ! ************************************************************************* 646 647 displ_red = zero 648 649 nbranch=3*natom 650 651 do jbranch=1,nbranch 652 ! 653 do iatom=1,natom 654 do idir=1,3 655 ibranch=idir+3*(iatom-1) 656 do kdir=1,3 657 k1 = kdir+3*(iatom-1) 658 ! WARNING: could be non-transpose of rprimd matrix : to be checked. 659 ! 23 june 2004: rprimd becomes gprimd 660 ! could be gprim and then multiply by acell... 661 ! Nope, checked and ok with gprimd 24 jun 2004 662 displ_red(1,ibranch,jbranch) = displ_red(1,ibranch,jbranch) + & 663 & gprimd(kdir,idir) * displ_cart(1,k1,jbranch) 664 665 displ_red(2,ibranch,jbranch) = displ_red(2,ibranch,jbranch) + & 666 & gprimd(kdir,idir) * displ_cart(2,k1,jbranch) 667 668 end do !kdir 669 end do !idir 670 end do !iatom 671 ! 672 end do !jbranch 673 674 end subroutine phdispl_cart2red
m_geometry/normv_int_vector [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
normv_int_vector
FUNCTION
Returns the norm of an integer 3D vector expressed in reduced coordinates. either in real or reciprocal space. In the later case the factor 2pi has to be included, due to the conventions used in abinit to define the reciprocal lattice.
INPUTS
OUTPUT
SOURCE
132 function normv_int_vector(xv,met,space) result(res) 133 134 135 !This section has been created automatically by the script Abilint (TD). 136 !Do not modify the following lines by hand. 137 #undef ABI_FUNC 138 #define ABI_FUNC 'normv_int_vector' 139 !End of the abilint section 140 141 implicit none 142 143 !Arguments ------------------------------------ 144 !scalars 145 real(dp) :: res 146 character(len=1),intent(in) :: space 147 !arrays 148 real(dp),intent(in) :: met(3,3) 149 integer,intent(in) :: xv(3) 150 151 ! ************************************************************************* 152 153 res = ( xv(1)*met(1,1)*xv(1) + xv(2)*met(2,2)*xv(2) + xv(3)*met(3,3)*xv(3) & 154 & +two*( xv(1)*met(1,2)*xv(2) + xv(1)*met(1,3)*xv(3) + xv(2)*met(2,3)*xv(3)) ) 155 156 select case (space) 157 case ('r','R') 158 res=SQRT(res) 159 case ('g','G') 160 res=two_pi*SQRT(res) 161 case default 162 MSG_BUG('Wrong value for space') 163 end select 164 165 end function normv_int_vector
m_geometry/normv_int_vector_array [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
normv_int_vector_array
FUNCTION
Returns the norm of an array of integer 3D vectors expressed in reduced coordinates. either in real or reciprocal space. In the later case the factor 2pi has to be included, due to the conventions used in abinit to define the reciprocal lattice.
INPUTS
OUTPUT
SOURCE
185 function normv_int_vector_array(xv,met,space) result(res) 186 187 188 !This section has been created automatically by the script Abilint (TD). 189 !Do not modify the following lines by hand. 190 #undef ABI_FUNC 191 #define ABI_FUNC 'normv_int_vector_array' 192 !End of the abilint section 193 194 implicit none 195 196 !Arguments ------------------------------------ 197 !scalars 198 character(len=1),intent(in) :: space 199 !arrays 200 real(dp),intent(in) :: met(3,3) 201 integer,intent(in) :: xv(:,:) 202 !this awful trick is needed to avoid problems with abilint 203 real(dp) :: res(SIZE(xv(1,:))) 204 !real(dp) :: res(SIZE(xv,DIM=2)) 205 206 ! ************************************************************************* 207 208 res(:) = ( xv(1,:)*met(1,1)*xv(1,:) + xv(2,:)*met(2,2)*xv(2,:) + xv(3,:)*met(3,3)*xv(3,:) & 209 & +two*(xv(1,:)*met(1,2)*xv(2,:) + xv(1,:)*met(1,3)*xv(3,:) + xv(2,:)*met(2,3)*xv(3,:)) ) 210 211 select case (space) 212 case ('r','R') 213 res(:)=SQRT(res(:)) 214 case ('g','G') 215 res(:)=two_pi*SQRT(res(:)) 216 case default 217 MSG_BUG('Wrong value for space') 218 end select 219 220 end function normv_int_vector_array
m_geometry/normv_rdp_vector [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
normv_rdp_vector
FUNCTION
Compute the norm of a vector expressed in reduced coordinates using the metric met. The result is multiplied by 2pi in case of a vector in reciprocal space to take into account the correct normalisation of the reciprocal lattice vectors
INPUTS
xv(3)=Vector in reduced coordinates met(3,3)=Metric tensor space=Character defining whether we are working in real (r|R) or reciprocal space (g|G)
OUTPUT
normv_rdp_vector=norm of xv
NOTES
The routine is able to deal both with a single vector as well as arrays of vectors. Versions for integer and real vectors are provided.
PARENTS
CHILDREN
SOURCE
80 function normv_rdp_vector(xv,met,space) result(res) 81 82 83 !This section has been created automatically by the script Abilint (TD). 84 !Do not modify the following lines by hand. 85 #undef ABI_FUNC 86 #define ABI_FUNC 'normv_rdp_vector' 87 !End of the abilint section 88 89 implicit none 90 91 !Arguments ------------------------------------ 92 !scalars 93 real(dp) :: res 94 character(len=1),intent(in) :: space 95 !arrays 96 real(dp),intent(in) :: met(3,3),xv(3) 97 98 ! ************************************************************************* 99 100 res = (xv(1)*met(1,1)*xv(1) + xv(2)*met(2,2)*xv(2) + xv(3)*met(3,3)*xv(3) & 101 & +two*(xv(1)*met(1,2)*xv(2) + xv(1)*met(1,3)*xv(3) + xv(2)*met(2,3)*xv(3)) ) 102 103 select case (space) 104 case ('r','R') 105 res=SQRT(res) 106 case ('g','G') 107 res=two_pi*SQRT(res) 108 case default 109 MSG_BUG('Wrong value for space') 110 end select 111 112 end function normv_rdp_vector
m_geometry/normv_rdp_vector_array [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
normv_rdp_vector_array
FUNCTION
Returns the norm of an array of real 3D vectors expressed in reduced coordinates. either in real or reciprocal space. In the later case the factor 2pi has to be included, due to the conventions used in abinit to define the reciprocal lattice.
INPUTS
OUTPUT
SOURCE
240 function normv_rdp_vector_array(xv,met,space) result(res) 241 242 243 !This section has been created automatically by the script Abilint (TD). 244 !Do not modify the following lines by hand. 245 #undef ABI_FUNC 246 #define ABI_FUNC 'normv_rdp_vector_array' 247 !End of the abilint section 248 249 implicit none 250 251 !Arguments ------------------------------------ 252 !scalars 253 character(len=1),intent(in) :: space 254 !arrays 255 real(dp),intent(in) :: met(3,3) 256 real(dp),intent(in) :: xv(:,:) 257 !this awful trick is needed to avoid problems with abilint 258 real(dp) :: res(SIZE(xv(1,:))) 259 !real(dp) :: res(SIZE(xv,DIM=2)) 260 261 ! ************************************************************************* 262 263 res(:) = ( xv(1,:)*met(1,1)*xv(1,:) + xv(2,:)*met(2,2)*xv(2,:) + xv(3,:)*met(3,3)*xv(3,:) & 264 & +two*(xv(1,:)*met(1,2)*xv(2,:) + xv(1,:)*met(1,3)*xv(3,:) + xv(2,:)*met(2,3)*xv(3,:)) ) 265 266 select case (space) 267 case ('r','R') 268 res(:)=SQRT(res(:)) 269 case ('g','G') 270 res(:)=two_pi*SQRT(res) 271 case default 272 MSG_BUG('Wrong value for space') 273 end select 274 275 end function normv_rdp_vector_array
m_geometry/vdotw_rc_vector [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
vdotw_rc_vector
FUNCTION
Compute the scalar product between two vectors expressed in reduced coordinates First vector is real, the second one is complex. The result is multiplied by (2pi)**2 in case of vectors in reciprocal space to take into account the correct normalisation of the reciprocal lattice vectors
INPUTS
xv(3),xw(3)=Vectors in reduced coordinates met(3,3)=Metric tensor space=Character defining whether we are working in real (r) or reciprocal space (g)
OUTPUT
res=complex scalar product of xv and xw
PARENTS
CHILDREN
SOURCE
368 function vdotw_rc_vector(xv,xw,met,space) result(res) 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 'vdotw_rc_vector' 375 !End of the abilint section 376 377 implicit none 378 379 !Arguments ------------------------------------ 380 !scalars 381 complex(dpc) :: res 382 character(len=1),intent(in) :: space 383 !arrays 384 real(dp),intent(in) :: met(3,3),xv(3) 385 complex(dpc),intent(in) :: xw(3) 386 387 ! ************************************************************************* 388 389 res = ( met(1,1)* xv(1)*xw(1) & 390 & +met(2,2)* xv(2)*xw(2) & 391 & +met(3,3)* xv(3)*xw(3) & 392 & +met(1,2)*(xv(1)*xw(2) + xv(2)*xw(1)) & 393 & +met(1,3)*(xv(1)*xw(3) + xv(3)*xw(1)) & 394 & +met(2,3)*(xv(2)*xw(3) + xv(3)*xw(2)) ) 395 396 select case (space) 397 case ('r','R') 398 return 399 case ('g','G') 400 res= res * (two_pi**2) 401 case default 402 MSG_BUG('Wrong value for space') 403 end select 404 405 end function vdotw_rc_vector
m_geometry/vdotw_rr_vector [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
vdotw_rr_vector
FUNCTION
Compute the scalar product between two vectors expressed in reduced coordinates The result is multiplied by (2pi)**2 in case of vectors in reciprocal space to take into account the correct normalisation of the reciprocal lattice vectors
INPUTS
xv(3),xw(3)=Vectors in reduced coordinates met(3,3)=Metric tensor space=Character defining whether we are working in real (r) or reciprocal space (g)
OUTPUT
res=scalar product of xv and xw
PARENTS
CHILDREN
SOURCE
303 function vdotw_rr_vector(xv,xw,met,space) result(res) 304 305 306 !This section has been created automatically by the script Abilint (TD). 307 !Do not modify the following lines by hand. 308 #undef ABI_FUNC 309 #define ABI_FUNC 'vdotw_rr_vector' 310 !End of the abilint section 311 312 implicit none 313 314 !Arguments ------------------------------------ 315 !scalars 316 real(dp) :: res 317 character(len=1),intent(in) :: space 318 !arrays 319 real(dp),intent(in) :: met(3,3),xv(3),xw(3) 320 321 ! ************************************************************************* 322 323 res = ( met(1,1)* xv(1)*xw(1) & 324 & +met(2,2)* xv(2)*xw(2) & 325 & +met(3,3)* xv(3)*xw(3) & 326 & +met(1,2)*(xv(1)*xw(2) + xv(2)*xw(1)) & 327 & +met(1,3)*(xv(1)*xw(3) + xv(3)*xw(1)) & 328 & +met(2,3)*(xv(2)*xw(3) + xv(3)*xw(2)) ) 329 330 select case (space) 331 case ('r','R') 332 return 333 case ('g','G') 334 res= res * (two_pi**2) 335 case default 336 MSG_BUG('Wrong value for space') 337 end select 338 339 end function vdotw_rr_vector
m_geometry/wigner_seitz [ Functions ]
[ Top ] [ m_geometry ] [ Functions ]
NAME
wigner_seitz
FUNCTION
Calculates a grid of points that falls inside of (and eventually on the surface of) the Wigner-Seitz supercell centered on the origin of the B lattice with primitive translations nmonkh(1)*a_1+nmonkh(2)*a_2+nmonkh(3)*a_3. Subroutine taken from the Wannier90 code. Modified by MG to fulfil abinit coding rules. API slightly changed the wrt wannier90 version.
COPYRIGHT
Copyright (C) 2007 Jonathan Yates, Arash Mostofi, Young-Su Lee, Nicola Marzari, Ivo Souza, David Vanderbilt. This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
center(3)=The Wigner-Seitz cell is centered on this point in reduced coordinates. rmet(3,3)=Real space metric ($\textrm{bohr}^{2}$). kptrlatt(3)=Values defining the supercell. prtvol=If different from 0 print out the points falling inside the W-S cell and the correponding weights. lmax(3)=see Notes below.
OUTPUT
npts=number of points falling inside the Wigner-Seitz cell irvec(3,npts)=Reduced coordinated of the points inside the W-S cell ndegen(npts)=Weigths associated to each point.
SIDE EFFECTS
In input irvec and ndegen are NULL pointers. They are allocated with the correct size inside the routine and returned to the caller.
NOTES
The Wannier functions live in a supercell of the real space unit cell. This supercell is mp_grid unit cells long in each direction The algorithm loops over grid points r on a unit cell that is 8 times larger than this primitive supercell. One of these points is in the W-S cell if it is closer to center(:) than any of the other points R where R are the translation vectors of the supercell. In the end npts contains the total number of grid points that have been found in the Wigner-Seitz cell The number of lattice vectors R along each direction of the supercell is defined by lmax.
PARENTS
CHILDREN
SOURCE
460 subroutine wigner_seitz(center,lmax,kptrlatt,rmet,npts,irvec,ndegen,prtvol) 461 462 463 !This section has been created automatically by the script Abilint (TD). 464 !Do not modify the following lines by hand. 465 #undef ABI_FUNC 466 #define ABI_FUNC 'wigner_seitz' 467 use interfaces_14_hidewrite 468 !End of the abilint section 469 470 implicit none 471 472 !Arguments ------------------------------------ 473 !scalars 474 integer,optional,intent(in) :: prtvol 475 integer,intent(out) :: npts 476 !arrays 477 integer,intent(in) :: kptrlatt(3,3),lmax(3) 478 integer,pointer :: irvec(:,:),ndegen(:) 479 real(dp),intent(in) :: center(3),rmet(3,3) 480 481 !Local variables------------------------------- 482 !scalars 483 integer :: in1,in2,in3,l1,l2,l3,ii,icount,n1,n2,n3 484 integer :: l0,l1_max,l2_max,l3_max,nl,verbose,mm1,mm2,mm3 485 real(dp) :: tot,dist_min 486 real(dp),parameter :: TOL_DIST=tol6 487 character(len=500) :: msg 488 !arrays 489 real(dp) :: diff(3) 490 real(dp),allocatable :: dist(:) 491 real(dp),allocatable :: swap2(:,:),swap1(:) 492 493 ! ************************************************************************* 494 495 verbose=0; if (PRESENT(prtvol)) verbose=prtvol 496 497 if (kptrlatt(1,2)/=0 .or. kptrlatt(2,1)/=0 .or. & 498 & kptrlatt(1,3)/=0 .or. kptrlatt(3,1)/=0 .or. & 499 & kptrlatt(2,3)/=0 .or. kptrlatt(3,2)/=0 ) then 500 MSG_ERROR('Off-diagonal elements of kptrlatt must be zero') 501 end if 502 503 n1=kptrlatt(1,1) 504 n2=kptrlatt(2,2) 505 n3=kptrlatt(3,3) 506 507 l1_max=lmax(1) 508 l2_max=lmax(2) 509 l3_max=lmax(3) 510 511 nl=(2*l1_max+1)*(2*l2_max+1)*(2*l3_max+1) 512 l0=1+l1_max*(1+(2*l2_max+1)**2+(2*l3_max+1)) ! Index of the origin. 513 ABI_MALLOC(dist,(nl)) 514 515 ! Allocate with maximum size 516 mm1=2*n1+1 517 mm2=2*n2+1 518 mm3=2*n3+1 519 ABI_MALLOC(irvec,(3,mm1*mm2*mm3)) 520 ABI_MALLOC(ndegen,(mm1*mm2*mm3)) 521 522 npts=0 523 do in1=-n1,n1 524 do in2=-n2,n2 525 do in3=-n3,n3 526 ! 527 ! Loop over the nl points R. R=0 corresponds to l1=l2=l3=1, or icount=l0 528 icount=0 529 do l1=-l1_max,l1_max 530 do l2=-l2_max,l2_max 531 do l3=-l3_max,l3_max 532 ! * Calculate |r-R-r_0|^2. 533 diff(1)= in1 -l1*n1 -center(1) 534 diff(2)= in2 -l2*n2 -center(2) 535 diff(3)= in3 -l3*n3 -center(3) 536 icount=icount+1 537 dist(icount)=DOT_PRODUCT(diff,MATMUL(rmet,diff)) 538 end do 539 end do 540 end do 541 542 dist_min=MINVAL(dist) 543 544 if (ABS(dist(l0)-dist_min)<TOL_DIST) then 545 npts=npts+1 546 ndegen(npts)=0 547 do ii=1,nl 548 if (ABS(dist(ii)-dist_min)<TOL_DIST) ndegen(npts)=ndegen(npts)+1 549 end do 550 irvec(1,npts)=in1 551 irvec(2,npts)=in2 552 irvec(3,npts)=in3 553 end if 554 end do !in3 555 end do !in2 556 end do !in1 557 558 if (verbose>=1) then 559 write(msg,'(a,i0)')' lattice points in Wigner-Seitz supercell: ',npts 560 call wrtout(std_out,msg,'COLL') 561 do ii=1,npts 562 write(msg,'(a,3(i3),a,i4)')' vector ', irvec(:,ii),' degeneracy: ', ndegen(ii) 563 call wrtout(std_out,msg,'COLL') 564 end do 565 end if 566 567 ! === Check the "sum rule" === 568 tot=zero 569 do ii=1,npts 570 tot=tot+one/ndegen(ii) 571 end do 572 if (ABS(tot-(n1*n2*n3))>tol8) then 573 write(msg,'(a,es16.8,a,i5)')'Something wrong in the generation of the mesh ',tot,' /= ',n1*n2*n3 574 MSG_ERROR(msg) 575 end if 576 577 ABI_FREE(dist) 578 579 ! === Reallocate the output with correct size === 580 ABI_MALLOC(swap1,(npts)) 581 swap1(:)=ndegen(1:npts) 582 ABI_FREE(ndegen) 583 ABI_MALLOC(ndegen,(npts)) 584 ndegen=swap1 585 ABI_FREE(swap1) 586 587 ABI_MALLOC(swap2,(3,npts)) 588 swap2(:,:)=irvec(1:3,1:npts) 589 ABI_FREE(irvec) 590 ABI_MALLOC(irvec,(3,npts)) 591 irvec=swap2 592 ABI_FREE(swap2) 593 594 end subroutine wigner_seitz