The Elk Code
 
Loading...
Searching...
No Matches
genscss.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 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 genscss
7use modmain
8implicit none
9! local variables
10integer is,ia,na,i
11real(8) vc(3),cs,sn,t1
12! automatic arrays
13real(8) vsc(3,nqptnr)
14! find the smallest supercell which contains q-vector
15call findscq(iqss,nscss,vsc)
16! construct supercell atomic positions and magnetic fields
17do is=1,nspecies
18 na=0
19 do ia=1,natoms0(is)
20 do i=1,nscss
21 na=na+1
22 if (na > maxatoms) then
23 write(*,*)
24 write(*,'("Error(genscss): too many atoms in supercell : ",I8)') na
25 write(*,'(" for species ",I4)') is
26 write(*,'("Adjust maxatoms in modmain and recompile code")')
27 write(*,*)
28 stop
29 end if
30 vc(:)=vsc(:,i)+atposc0(:,ia,is)
31! new atomic position in lattice coordinates
32 call r3mv(ainv,vc,atposl(:,na,is))
33! rotate external B-field and fixed spin moment vector by angle q.r
34 t1=dot_product(vqc(:,iqss),vc(:))
35 cs=cos(t1); sn=sin(t1)
36 bfcmt0(1,na,is)=cs*bfcmt00(1,ia,is)-sn*bfcmt00(2,ia,is)
37 bfcmt0(2,na,is)=sn*bfcmt00(1,ia,is)+cs*bfcmt00(2,ia,is)
38 bfcmt0(3,na,is)=bfcmt00(3,ia,is)
39 mommtfix(1,na,is)=cs*mommtfix0(1,ia,is)-sn*mommtfix0(2,ia,is)
40 mommtfix(2,na,is)=sn*mommtfix0(1,ia,is)+cs*mommtfix0(2,ia,is)
41 mommtfix(3,na,is)=mommtfix0(3,ia,is)
42 end do
43 end do
44 natoms(is)=na
45end do
46end subroutine
47
subroutine findscq(iq, nsc, vsc)
Definition findscq.f90:7
subroutine genscss
Definition genscss.f90:7
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(3, maxatoms, maxspecies) mommtfix0
Definition modmain.f90:259
real(8), dimension(:,:), allocatable vqc
Definition modmain.f90:547
real(8), dimension(3, maxatoms, maxspecies) bfcmt00
Definition modmain.f90:275
integer iqss
Definition modmain.f90:297
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 nscss
Definition modmain.f90:299
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
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10