The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine gentpmae
7use modmain
8implicit none
9! local variables
10integer na,n,i1,i2,i3,i
11integer isym,lspl
12real(8) v1(3),v2(3),t1
13! allocatable arrays
14real(8), allocatable :: vp(:,:)
15if (allocated(tpmae)) deallocate(tpmae)
16select case(npmae0)
17case(4:)
18! distribute points evenly on a sphere
20 allocate(tpmae(2,npmae))
21 call sphcover(npmae,tpmae)
22case(-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)
4710 continue
48 end do
49 end do
50 end do
51! convert vectors to spherical coordinates
52 allocate(tpmae(2,npmae))
53 do i=1,npmae
54 call sphcrd(vp(:,i),t1,tpmae(:,i))
55 end do
56 deallocate(vp)
57case(2)
58! use x- and z-directions
59 npmae=2
60 allocate(tpmae(2,npmae))
61 tpmae(1,1)=pi/2.d0
62 tpmae(2,1)=0.d0
63 tpmae(1,2)=0.d0
64 tpmae(2,2)=0.d0
65case(3)
66! use x-, y- and z-directions
67 npmae=3
68 allocate(tpmae(2,npmae))
69 tpmae(1,1)=pi/2.d0
70 tpmae(2,1)=0.d0
71 tpmae(1,2)=pi/2.d0
72 tpmae(2,2)=pi/2.d0
73 tpmae(1,3)=0.d0
74 tpmae(2,3)=0.d0
75case default
76 write(*,*)
77 write(*,'("Error(gentpmae): invalid npmae : ",I8)') npmae0
78 write(*,*)
79 stop
80end select
81end subroutine
82
subroutine gentpmae
Definition gentpmae.f90:7
real(8), parameter pi
Definition modmain.f90:1229
integer npmae
Definition modmain.f90:302
integer npmae0
Definition modmain.f90:302
real(8), dimension(:,:), allocatable tpmae
Definition modmain.f90:304
real(8) epslat
Definition modmain.f90:24
integer nsymcrys
Definition modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition modmain.f90:344
integer, dimension(maxsymcrys) lsplsymc
Definition modmain.f90:364
subroutine sphcover(n, tp)
Definition sphcover.f90:10
pure subroutine sphcrd(v, r, tp)
Definition sphcrd.f90:10