The Elk Code
reciplat.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: reciplat
8 ! !INTERFACE:
9 subroutine reciplat(avec,bvec,omega,omegabz)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! avec : lattice vectors (in,real(3,3))
12 ! bvec : reciprocal lattice vectors (out,real(3,3))
13 ! omega : unit cell volume (out,real)
14 ! omegabz : Brillouin zone volume (out,real)
15 ! !DESCRIPTION:
16 ! Generates the reciprocal lattice vectors from the real-space lattice vectors
17 ! \begin{align*}
18 ! {\bf b}_1&=\frac{2\pi}{s}({\bf a}_2\times{\bf a}_3)\\
19 ! {\bf b}_2&=\frac{2\pi}{s}({\bf a}_3\times{\bf a}_1)\\
20 ! {\bf b}_3&=\frac{2\pi}{s}({\bf a}_1\times{\bf a}_2)
21 ! \end{align*}
22 ! and finds the unit cell volume $\Omega=|s|$, where
23 ! $s={\bf a}_1\cdot({\bf a}_2\times{\bf a}_3)$, and the Brillouin zone volume
24 ! $\Omega_{\rm BZ}=(2\pi)^3/\Omega$.
25 !
26 ! !REVISION HISTORY:
27 ! Created September 2002 (JKD)
28 !EOP
29 !BOC
30 implicit none
31 ! arguments
32 real(8), intent(in) :: avec(3,3)
33 real(8), intent(out) :: bvec(3,3),omega,omegabz
34 ! local variables
35 real(8), parameter :: twopi=6.2831853071795864769d0
36 real(8) t1
37 call r3cross(avec(:,2),avec(:,3),bvec(:,1))
38 call r3cross(avec(:,3),avec(:,1),bvec(:,2))
39 call r3cross(avec(:,1),avec(:,2),bvec(:,3))
40 t1=avec(1,1)*bvec(1,1)+avec(2,1)*bvec(2,1)+avec(3,1)*bvec(3,1)
41 ! unit cell volume
42 omega=abs(t1)
43 if (omega < 1.d-6) then
44  write(*,*)
45  write(*,'("Error(reciplat) omega too small : ",G18.10)') omega
46  write(*,'(" Lattice vectors may be collinear")')
47  write(*,*)
48  stop
49 end if
50 bvec(:,:)=(twopi/t1)*bvec(:,:)
51 ! Brillouin zone volume
52 omegabz=(twopi**3)/omega
53 end subroutine
54 !EOC
55 
subroutine reciplat(avec, bvec, omega, omegabz)
Definition: reciplat.f90:10
pure subroutine r3cross(x, y, z)
Definition: r3cross.f90:10