The Elk Code
 
Loading...
Searching...
No Matches
strainabg.f90
Go to the documentation of this file.
1
2! Copyright (C) 2016 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine strainabg
7use modmain
8use modgw
9use modvars
10implicit none
11integer is,ia,ig
12real(8) ta(3,3),tb(3,3),vc(3)
13! compute the strained lattice vectors
14if ((istrain >= 1).and.(istrain <= nstrain)) then
15 avec(1:3,1:3)=avec0(1:3,1:3)+deltast*strain(1:3,1:3,istrain)
16else if (tavref) then
17 avec(1:3,1:3)=avec0(1:3,1:3)
18else
19 return
20end if
21! generate the strained reciprocal lattice vectors and unit cell volume
23! determine the transformation matrix to the strained vectors
24call r3mm(avec,ainv,ta)
25call r3mm(bvec,binv,tb)
26! recalculate all required variables which depend on avec
27call r3minv(avec,ainv)
28call r3minv(bvec,binv)
29call r3mv(bvec,vqlss,vqcss)
30do is=1,nspecies
31 do ia=1,natoms(is)
32 call r3mv(avec,atposl(:,ia,is),atposc(:,ia,is))
33 end do
34end do
35call r3mv(bvec,vecql,vecqc)
38call symmetry
39! apply the transformation matrix to the G-vectors
40do ig=1,ngtot
41 vc(1:3)=vgc(1:3,ig)
42 call r3mv(tb,vc,vgc(:,ig))
43 gc(ig)=sqrt(vgc(1,ig)**2+vgc(2,ig)**2+vgc(3,ig)**2)
44end do
45end subroutine
46
Definition modgw.f90:6
real(8) deltast
Definition modmain.f90:1019
integer nstrain
Definition modmain.f90:1012
real(8), dimension(3, 3) avec0
Definition modmain.f90:12
integer ngtot
Definition modmain.f90:390
real(8), dimension(3) efieldc
Definition modmain.f90:312
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition modmain.f90:54
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(3) afieldc
Definition modmain.f90:325
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
real(8), dimension(3) vecql
Definition modmain.f90:1101
real(8) omega
Definition modmain.f90:20
real(8), dimension(3) vecqc
Definition modmain.f90:1101
real(8), dimension(3, 3) binv
Definition modmain.f90:18
real(8), dimension(3, 3) avec
Definition modmain.f90:12
real(8), dimension(3) efieldl
Definition modmain.f90:314
real(8), dimension(3) vqcss
Definition modmain.f90:295
logical tavref
Definition modmain.f90:1029
real(8) omegabz
Definition modmain.f90:22
real(8), dimension(3, 3, 9) strain
Definition modmain.f90:1016
integer istrain
Definition modmain.f90:1014
real(8), dimension(3, 3) ainv
Definition modmain.f90:14
real(8), dimension(3) afieldl
Definition modmain.f90:327
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
real(8), dimension(3) vqlss
Definition modmain.f90:293
integer nspecies
Definition modmain.f90:34
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition modmain.f90:51
subroutine r3minv(a, b)
Definition r3minv.f90:10
pure subroutine r3mm(a, b, c)
Definition r3mm.f90:10
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
subroutine reciplat(avec, bvec, omega, omegabz)
Definition reciplat.f90:10
subroutine strainabg
Definition strainabg.f90:7
subroutine symmetry
Definition symmetry.f90:7