The Elk Code
checkfsm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2009 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 checkfsm
7 use modmain
8 implicit none
9 ! local variables
10 integer isym,lspn
11 integer is,ia,ja
12 real(8) sc(3,3),v(3),t1
13 if (fsmtype == 0) return
14 do isym=1,nsymcrys
15  lspn=lspnsymc(isym)
16 ! proper rotation matrix in Cartesian coordinates
17  sc(:,:)=dble(symlatd(lspn))*symlatc(:,:,lspn)
18 ! check invariance of global moment
19  if ((abs(fsmtype) == 1).or.(abs(fsmtype) == 3)) then
20  call r3mv(sc,momfix,v)
21  t1=sum(abs(momfix(:)-v(:)))
22  if (t1 > epslat) then
23  write(*,*)
24  write(*,'("Error(checkfsm): momfix not invariant under symmetry group")')
25  write(*,*)
26  stop
27  end if
28  end if
29 ! check invariance of muffin-tin moments
30  if ((abs(fsmtype) == 2).or.(abs(fsmtype) == 3)) then
31  do is=1,nspecies
32  do ia=1,natoms(is)
33 ! if any component is >= 1000 then do not fix the moment
34  t1=sum(abs(mommtfix(:,ia,is)))
35  if (t1 >= 1000.d0) cycle
36 ! equivalent atom
37  ja=ieqatom(ia,is,isym)
38  call r3mv(sc,mommtfix(:,ja,is),v)
39  t1=sum(abs(mommtfix(:,ia,is)-v(:)))
40  if (t1 > epslat) then
41  write(*,*)
42  write(*,'("Error(checkfsm): mommtfix not invariant under symmetry &
43  &group")')
44  write(*,'(" for species ",I4)') is
45  write(*,'(" and equivalent atoms ",2I4)') ia,ja
46  write(*,*)
47  stop
48  end if
49  end do
50  end do
51  end if
52 end do
53 end subroutine
54 
integer, dimension(maxsymcrys) lspnsymc
Definition: modmain.f90:366
integer nsymcrys
Definition: modmain.f90:358
integer, dimension(:,:,:), allocatable ieqatom
Definition: modmain.f90:368
real(8), dimension(3) momfix
Definition: modmain.f90:253
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition: modmain.f90:259
integer, dimension(48) symlatd
Definition: modmain.f90:346
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
real(8) epslat
Definition: modmain.f90:24
subroutine checkfsm
Definition: checkfsm.f90:7
integer nspecies
Definition: modmain.f90:34
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
integer fsmtype
Definition: modmain.f90:251