The Elk Code
gentpmae.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2013 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
6
subroutine
gentpmae
7
use
modmain
8
implicit none
9
! local variables
10
integer
na,n,i1,i2,i3,i
11
integer
isym,lspl
12
real(8)
v1(3),v2(3),t1
13
! allocatable arrays
14
real(8)
,
allocatable
:: vp(:,:)
15
if
(
allocated
(
tpmae
))
deallocate
(
tpmae
)
16
select case
(
npmae0
)
17
case
(4:)
18
! distribute points evenly on a sphere
19
npmae
=
npmae0
20
allocate
(
tpmae
(2,
npmae
))
21
call
sphcover
(
npmae
,
tpmae
)
22
case
(-4:-1)
23
! use symmetry reduced cardinal directions
24
na=abs(
npmae0
)
25
n=(2*na+1)**3
26
allocate
(vp(3,n))
27
npmae
=0
28
do
i1=-na,na
29
v1(1)=dble(i1)
30
do
i2=-na,na
31
v1(2)=dble(i2)
32
do
i3=-na,na
33
v1(3)=dble(i3)
34
if
((i1 == 0).and.(i2 == 0).and.(i3 == 0)) cycle
35
do
isym=1,
nsymcrys
36
lspl=
lsplsymc
(isym)
37
v2(1:3)=
symlat
(1:3,1,lspl)*v1(1) &
38
+
symlat
(1:3,2,lspl)*v1(2) &
39
+
symlat
(1:3,3,lspl)*v1(3)
40
do
i=1,
npmae
41
t1=abs(vp(1,i)-v2(1))+abs(vp(2,i)-v2(2))+abs(vp(3,i)-v2(3))
42
if
(t1 <
epslat
)
goto
10
43
end do
44
end do
45
npmae
=
npmae
+1
46
vp(1:3,
npmae
)=v1(1:3)
47
10
continue
48
end do
49
end do
50
end do
51
allocate
(
tpmae
(2,
npmae
))
52
do
i=1,
npmae
53
! convert vector to Cartesian coordinates
54
call
r3mv
(
avec
,vp(:,i),v1)
55
! convert vector to spherical coordinates
56
call
sphcrd
(v1,t1,
tpmae
(:,i))
57
end do
58
deallocate
(vp)
59
case
(2)
60
! use x- and z-directions
61
npmae
=2
62
allocate
(
tpmae
(2,
npmae
))
63
tpmae
(1,1)=
pi
/2.d0
64
tpmae
(2,1)=0.d0
65
tpmae
(1,2)=0.d0
66
tpmae
(2,2)=0.d0
67
case
(3)
68
! use x-, y- and z-directions
69
npmae
=3
70
allocate
(
tpmae
(2,
npmae
))
71
tpmae
(1,1)=
pi
/2.d0
72
tpmae
(2,1)=0.d0
73
tpmae
(1,2)=
pi
/2.d0
74
tpmae
(2,2)=
pi
/2.d0
75
tpmae
(1,3)=0.d0
76
tpmae
(2,3)=0.d0
77
case default
78
write
(*,*)
79
write
(*,
'("Error(gentpmae): invalid npmae : ",I0)'
)
npmae0
80
write
(*,*)
81
stop
82
end select
83
end subroutine
84
modmain::tpmae
real(8), dimension(:,:), allocatable tpmae
Definition:
modmain.f90:304
modmain::nsymcrys
integer nsymcrys
Definition:
modmain.f90:358
modmain::symlat
integer, dimension(3, 3, 48) symlat
Definition:
modmain.f90:344
gentpmae
subroutine gentpmae
Definition:
gentpmae.f90:7
modmain::pi
real(8), parameter pi
Definition:
modmain.f90:1226
sphcrd
pure subroutine sphcrd(v, r, tp)
Definition:
sphcrd.f90:10
modmain
Definition:
modmain.f90:6
modmain::lsplsymc
integer, dimension(maxsymcrys) lsplsymc
Definition:
modmain.f90:364
modmain::avec
real(8), dimension(3, 3) avec
Definition:
modmain.f90:12
modmain::epslat
real(8) epslat
Definition:
modmain.f90:24
sphcover
subroutine sphcover(n, tp)
Definition:
sphcover.f90:10
r3mv
pure subroutine r3mv(a, x, y)
Definition:
r3mv.f90:10
modmain::npmae0
integer npmae0
Definition:
modmain.f90:302
modmain::npmae
integer npmae
Definition:
modmain.f90:302
gentpmae.f90
Generated by
1.8.14