TABLE OF CONTENTS
ABINIT/polcart [ Functions ]
NAME
polcart
FUNCTION
Transform polarization from reduced to cartesian coordinates, divide by ucvol and write the result to an output file
COPYRIGHT
Copyright (C) 2000-2018 ABINIT group (MVeithen) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
pel(3) = reduced coordinates of the electronic polarization pion(3) = reduced coordinates of the ionic polarization polunit = units used for the output of the polarization 1 : use atomic units 2 : use MKS units 3 : use both atomic and MKS units rprimd(3,3) = dimensional primitive translations (bohr) ucvol = volume of the primitive unit cell unit_out = unit for output of the results usepaw = 1 for PAW computation, zero else
OUTPUT
pel_cart(3) = cartesian coords of the electronic polarization in atomic units pelev(3)= expectation value polarization term (PAW only) in cartesian coordinates pion_cart(3)= cartesian coords of the ionic polarization in atomic units ptot_cart(3)= cartesian coords of the total polarization in atomic units
NOTES
The sum of the electronic and ionic Berry phase is folded into [-1,1] before it is transformed to cartesian coordinates. This means that in some cases, ptot_cart /= pel_cart + pion_cart
NOTES
- pel and pion do not take into account the factor 1/ucvol. At the opposite, this factor is taken into account in pel_cart and pion_cart - unit_out = 0 is allowed, in this case, there will be no output of the results
PARENTS
berryphase_new,relaxpol
CHILDREN
wrtout
SOURCE
58 #if defined HAVE_CONFIG_H 59 #include "config.h" 60 #endif 61 62 #include "abi_common.h" 63 64 subroutine polcart(red_ptot,pel,pel_cart,pelev,pion,pion_cart,polunit,& 65 & ptot_cart,rprimd,ucvol,unit_out,usepaw) 66 67 use defs_basis 68 use m_profiling_abi 69 70 !This section has been created automatically by the script Abilint (TD). 71 !Do not modify the following lines by hand. 72 #undef ABI_FUNC 73 #define ABI_FUNC 'polcart' 74 use interfaces_14_hidewrite 75 !End of the abilint section 76 77 implicit none 78 79 !Arguments ------------------------------------ 80 !scalars 81 integer,intent(in) :: polunit,unit_out,usepaw 82 real(dp),intent(in) :: ucvol 83 !arrays 84 real(dp),intent(in) :: red_ptot(3) !!REC 85 real(dp),intent(in) :: pel(3),pelev(3),pion(3),rprimd(3,3) 86 real(dp),intent(out) :: pel_cart(3),pion_cart(3),ptot_cart(3) 87 88 !Local variables ------------------------- 89 !scalars 90 integer :: idir 91 character(len=500) :: message 92 !arrays 93 real(dp) :: pel_mks(3),pelev_mks(3),pion_mks(3),ptot(3),ptot_mks(3) 94 95 ! *********************************************************************** 96 !!REC Note ptot has already been folded and kept onto same branch 97 !unless ptot=0d0, in which case ptot has not been computed yet 98 if( sum(abs(red_ptot(:))) < tol8 )then 99 ptot(:) = pel(:) + pion(:) 100 ! Fold ptot into [-1, 1] 101 do idir = 1, 3 102 ptot(idir) = ptot(idir) - 2_dp*nint(ptot(idir)/2_dp) 103 end do 104 else !!REC 105 ptot=red_ptot !!REC 106 end if !!REC 107 108 !Transform pel, pion and ptot to cartesian coordinates 109 pel_cart(:) = zero ; pion_cart(:) = zero ; ptot_cart(:) = zero 110 do idir = 1, 3 111 pel_cart(idir) = rprimd(idir,1)*pel(1) + rprimd(idir,2)*pel(2) + & 112 & rprimd(idir,3)*pel(3) 113 pion_cart(idir) = rprimd(idir,1)*pion(1) + rprimd(idir,2)*pion(2) + & 114 & rprimd(idir,3)*pion(3) 115 ptot_cart(idir) = rprimd(idir,1)*ptot(1) + rprimd(idir,2)*ptot(2) + & 116 & rprimd(idir,3)*ptot(3) 117 end do 118 119 !Divide by the unit cell volume 120 pel_cart(:) = pel_cart(:)/ucvol 121 pion_cart(:) = pion_cart(:)/ucvol 122 ptot_cart(:) = ptot_cart(:)/ucvol 123 !pelev is either zero (in NCPP case), or possibly non-zero (in PAW case) 124 !however, in the PAW case, it is already implicitly included in the computation 125 !of pel and should not be added in here. Below in the PAW case we report its 126 !value to the user, which is helpful for comparison to the USPP theory. 127 !13 June 2012 J Zwanziger 128 !note that pelev was AlREADY in cartesian frame. 129 !ptot_cart(:) = (ptot_cart(:)+pelev(:))/ucvol 130 131 !Write the results to unit_out (if /= 0) 132 !Use the coordinates specified by the value polunit 133 134 !Atomic units 135 136 if (((polunit == 1).or.(polunit == 3)).and.(unit_out /= 0)) then 137 138 ! write(message,'(7(a),3(e16.9,2x),a,a,3(e16.9,2x),a,a,3(e16.9,2x),a,a,3(e16.9,2x))')ch10,& 139 ! & ' Polarization in cartesian coordinates (a.u.):',ch10,& 140 ! & ' (the sum of the electronic and ionic Berry phase',& 141 ! & ' has been fold into [-1, 1])',ch10,& 142 ! & ' Electronic berry phase: ', (pel_cart(idir), idir = 1, 3), ch10,& 143 ! & ' Expectation value (PAW only): ', (pelev(idir)/ucvol, idir = 1, 3), ch10,& 144 ! & ' Ionic: ', (pion_cart(idir), idir = 1, 3), ch10, & 145 ! & ' Total: ', (ptot_cart(idir), idir = 1, 3) 146 ! call wrtout(unit_out,message,'COLL') 147 write(message,'(7(a),3(e16.9,2x),a,a,3(e16.9,2x))')ch10,& 148 & ' Polarization in cartesian coordinates (a.u.):',ch10,& 149 & ' (the sum of the electronic and ionic Berry phase',& 150 & ' has been folded into [-1, 1])',ch10,& 151 & ' Electronic berry phase: ', (pel_cart(idir), idir = 1, 3) 152 call wrtout(unit_out,message,'COLL') 153 if(usepaw==1) then 154 write(message,'(a,3(e16.9,2x))')& 155 & ' ...includes PAW on-site term: ', (pelev(idir)/ucvol, idir = 1, 3) 156 call wrtout(unit_out,message,'COLL') 157 end if 158 write(message,'(a,3(e16.9,2x),a,a,3(e16.9,2x))')& 159 & ' Ionic: ', (pion_cart(idir), idir = 1, 3), ch10, & 160 & ' Total: ', (ptot_cart(idir), idir = 1, 3) 161 call wrtout(unit_out,message,'COLL') 162 163 end if 164 165 !MKS units 166 167 if (((polunit == 2).or.(polunit == 3)).and.(unit_out /= 0)) then 168 169 pel_mks(:) = pel_cart(:)*(e_Cb)/(Bohr_Ang*1d-10)**2 170 pion_mks(:) = pion_cart(:)*(e_Cb)/(Bohr_Ang*1d-10)**2 171 pelev_mks(:) = pelev(:)/ucvol*(e_Cb)/(Bohr_Ang*1d-10)**2 172 ptot_mks(:) = (ptot_cart(:))*(e_Cb)/(Bohr_Ang*1d-10)**2 173 174 write(message,'(7(a),3(e16.9,2x),a,a,3(e16.9,2x),a,a,3(e16.9,2x),a,a,3(e16.9,2x))')ch10,& 175 & ' Polarization in cartesian coordinates (C/m^2):',ch10,& 176 & ' (the sum of the electronic and ionic Berry phase',& 177 & ' has been folded into [-1, 1])',ch10,& 178 & ' Electronic berry phase: ', (pel_mks(idir), idir = 1, 3) 179 call wrtout(unit_out,message,'COLL') 180 if(usepaw==1) then 181 write(message,'(a,3(e16.9,2x))')& 182 & ' ...includes PAW on-site term: ', (pelev_mks(idir), idir = 1, 3) 183 call wrtout(unit_out,message,'COLL') 184 end if 185 write(message,'(a,3(e16.9,2x),a,a,3(e16.9,2x))')& 186 & ' Ionic: ', (pion_mks(idir), idir = 1, 3), ch10, & 187 & ' Total: ', (ptot_mks(idir), idir = 1, 3) 188 call wrtout(unit_out,message,'COLL') 189 190 end if 191 192 end subroutine polcart