The Elk Code
 
Loading...
Searching...
No Matches
genscph.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 General Public License.
4! See the file COPYING for license details.
5
6subroutine genscph(p,dph)
7use modmain
8use modphonon
9implicit none
10! arguments
11integer, intent(in) :: p
12real(8), intent(in) :: dph
13! local variables
14integer is,ia,na,i
15real(8) vc(3),t1
16if ((p /= 0).and.(p /= 1)) then
17 write(*,*)
18 write(*,'("Error(genscph): phase (p) should be 0 or 1 : ",I8)') p
19 write(*,*)
20 stop
21end if
22! find the smallest supercell which contains the q-vector
24! construct supercell atomic positions and magnetic fields
25do is=1,nspecies
26 na=0
27 do ia=1,natoms0(is)
28 do i=1,nscph
29 na=na+1
30 if (na > maxatoms) then
31 write(*,*)
32 write(*,'("Error(genscph): too many atoms in supercell : ",I8)') na
33 write(*,'(" for species ",I4)') is
34 write(*,'("Adjust maxatoms in modmain and recompile code")')
35 write(*,*)
36 stop
37 end if
38 vc(:)=vscph(:,i)+atposc0(:,ia,is)
39! add small periodic displacement
40 if ((isph == is).and.(iaph == ia)) then
41 t1=dot_product(vqc(:,iqph),vscph(:,i))
42 if (p == 0) then
43 vc(ipph)=vc(ipph)+dph*cos(t1)
44 else
45 vc(ipph)=vc(ipph)+dph*sin(t1)
46 end if
47 end if
48! convert to new lattice coordinates
49 call r3mv(ainv,vc,atposl(:,na,is))
50 call r3frac(epslat,atposl(:,na,is))
51! set muffin-tin fields and fixed spin moments if required
52 if (spinpol) then
53 bfcmt0(:,na,is)=bfcmt00(:,ia,is)
54 mommtfix(:,na,is)=mommtfix0(:,ia,is)
55 end if
56 end do
57 end do
58 natoms(is)=na
59end do
60end subroutine
61
subroutine findscq(iq, nsc, vsc)
Definition findscq.f90:7
subroutine genscph(p, dph)
Definition genscph.f90:7
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical spinpol
Definition modmain.f90:228
real(8), dimension(3, maxatoms, maxspecies) mommtfix0
Definition modmain.f90:259
real(8), dimension(:,:), allocatable vqc
Definition modmain.f90:547
real(8) epslat
Definition modmain.f90:24
real(8), dimension(3, maxatoms, maxspecies) bfcmt00
Definition modmain.f90:275
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition modmain.f90:259
real(8), dimension(3, maxatoms, maxspecies) atposc0
Definition modmain.f90:54
real(8), dimension(3, 3) ainv
Definition modmain.f90:14
integer, dimension(maxspecies) natoms0
Definition modmain.f90:36
integer nspecies
Definition modmain.f90:34
integer, parameter maxatoms
Definition modmain.f90:32
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition modmain.f90:51
real(8), dimension(3, maxatoms, maxspecies) bfcmt0
Definition modmain.f90:275
integer iaph
Definition modphonon.f90:15
integer ipph
Definition modphonon.f90:15
integer iqph
Definition modphonon.f90:15
real(8), dimension(:,:), allocatable vscph
Definition modphonon.f90:46
integer isph
Definition modphonon.f90:15
integer nscph
Definition modphonon.f90:44
pure subroutine r3frac(eps, v)
Definition r3frac.f90:10
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10