TABLE OF CONTENTS
- ABINIT/pbc_lotf
- pbc/dist_pbc_ext
- pbc/dist_pbc_int
- pbc_lotf/pbc_aa_contract
- pbc_lotf/pbc_bb_contract
- pbc_lotf/pbc_bb_proj
- pbc_lotf/pbc_init
- pbc_lotf/vecp
ABINIT/pbc_lotf [ Modules ]
NAME
pbc_lotf
FUNCTION
COPYRIGHT
Copyright (C) 2005-2024 ABINIT group (MMancini) 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
15 #if defined HAVE_CONFIG_H 16 #include "config.h" 17 #endif 18 19 #include "abi_common.h" 20 21 module pbc_lotf 22 23 use defs_basis 24 25 implicit none 26 private 27 28 real(dp),private :: aa(3,3),bb(3,3) 29 real(dp),public :: rd(3) ! these are pbc output 30 real(dp),public :: r2 ! these are pbc output 31 32 public :: & 33 dist_pbc, & 34 pbc_init, & 35 pbc_aa_contract, & 36 pbc_bb_contract, & 37 pbc_bb_proj 38 39 private :: & 40 vecp 41 42 43 interface dist_pbc 44 module procedure dist_pbc_ext 45 module procedure dist_pbc_int 46 end interface dist_pbc 47 48 contains
pbc/dist_pbc_ext [ Functions ]
NAME
dist_pbc_ext
FUNCTION
INPUTS
SOURCE
183 subroutine dist_pbc_ext(RI,RJ,r2,RD) 184 ! ONLY aa AND bb MATRICES ARE USED IN THIS VERSION 185 implicit none 186 187 !Arguments ------------------------ 188 real(dp),intent(out) :: r2 189 real(dp),intent(in),dimension(3) :: RI, RJ 190 real(dp),intent(out),dimension(3) :: RD 191 !Local --------------------------- 192 integer :: ii 193 real(dp),dimension(3) :: rad, radi 194 195 ! ************************************************************************* 196 197 !--at These are cartesian: 198 rad = RI - RJ 199 200 !--at These are monoclinic: 201 radi(:) = matmul(rad,bb) 202 203 do ii=1,3 204 if(radi(ii) < -half) then 205 rad = rad + aa(:,ii) 206 elseif(radi(ii) > half) then 207 rad = rad - aa(:,ii) 208 endif 209 enddo 210 211 r2 = dot_product(rad,rad) 212 RD = rad 213 end subroutine dist_pbc_ext
pbc/dist_pbc_int [ Functions ]
NAME
dist_pbc_int
FUNCTION
INPUTS
SOURCE
226 subroutine dist_pbc_int(RI,RJ) 227 ! ONLY aa AND bb MATRICES ARE USED IN THIS VERSION 228 implicit none 229 230 !Arguments ------------------------ 231 real(dp),intent(in),dimension(3) :: RI, RJ 232 !Local --------------------------- 233 integer :: ii 234 real(dp),dimension(3) :: rad, radi 235 236 ! ************************************************************************* 237 238 !--at These are cartesian: 239 rad = RI - RJ 240 241 !--at These are monoclinic: 242 radi(:) = matmul(rad,bb) 243 244 do ii=1,3 245 if(radi(ii) < -half) then 246 rad = rad + aa(:,ii) 247 elseif(radi(ii) > half) then 248 rad = rad - aa(:,ii) 249 endif 250 enddo 251 252 r2 = dot_product(rad,rad) 253 RD = rad 254 end subroutine dist_pbc_int 255 256 end module pbc_lotf
pbc_lotf/pbc_aa_contract [ Functions ]
[ Top ] [ pbc_lotf ] [ Functions ]
NAME
pbc_aa_contract
FUNCTION
Compute the contraction of aa
INPUTS
SOURCE
113 function pbc_aa_contract() 114 115 implicit none 116 117 !Arguments ------------------------ 118 real(dp) :: pbc_aa_contract(3) 119 pbc_aa_contract = sqrt(sum(aa(:,:)**two,dim=1)) 120 return 121 end function pbc_aa_contract
pbc_lotf/pbc_bb_contract [ Functions ]
[ Top ] [ pbc_lotf ] [ Functions ]
NAME
pbc_bb_contract
FUNCTION
Compute the contraction of bb
INPUTS
SOURCE
134 function pbc_bb_contract() 135 136 implicit none 137 138 !Arguments ------------------------ 139 real(dp) :: pbc_bb_contract(3) 140 141 ! ************************************************************************* 142 143 pbc_bb_contract = sqrt(sum(bb(:,:)**two,dim=1)) 144 return 145 end function pbc_bb_contract
pbc_lotf/pbc_bb_proj [ Functions ]
[ Top ] [ pbc_lotf ] [ Functions ]
NAME
pbc_bb_proj
FUNCTION
Compute the application of a vector on bb
INPUTS
vi(3)=real vector
SOURCE
159 function pbc_bb_proj(vi) 160 161 implicit none 162 163 !Arguments ------------------------ 164 real(dp),intent(in) :: vi(3) 165 real(dp) :: pbc_bb_proj(3) 166 167 ! ************************************************************************* 168 169 pbc_bb_proj = matmul(vi,bb) 170 return 171 end function pbc_bb_proj
pbc_lotf/pbc_init [ Functions ]
[ Top ] [ pbc_lotf ] [ Functions ]
NAME
pbc_init
FUNCTION
INPUTS
SOURCE
83 subroutine pbc_init(rprimd) 84 85 implicit none 86 87 !Arguments ------------------------ 88 real(dp),intent(in) :: rprimd(3,3) 89 !Local --------------------------- 90 real(dp) :: avol 91 92 aa(:,:) = rprimd(:,:) 93 call vecp(aa(1,2),aa(1,3),bb(1,1)) !--b^c 94 call vecp(aa(1,3),aa(1,1),bb(1,2)) !--c^a 95 call vecp(aa(1,1),aa(1,2),bb(1,3)) !--a^b 96 avol = dot_product(aa(:,3),bb(:,3)) !--c.(a^b) 97 98 bb = (one/avol)*bb 99 end subroutine pbc_init
pbc_lotf/vecp [ Functions ]
[ Top ] [ pbc_lotf ] [ Functions ]
NAME
vecp
FUNCTION
INPUTS
SOURCE
60 subroutine vecp(a,b,c) 61 62 real(dp),intent(in) :: a(3),b(3) 63 real(dp),intent(out) :: c(3) 64 65 ! ************************************************************************* 66 67 c(1) = a(2) * b(3) - b(2) * a(3) 68 c(2) = a(3) * b(1) - b(3) * a(1) 69 c(3) = a(1) * b(2) - b(1) * a(2) 70 end subroutine vecp