The Elk Code
sfacmag.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2010 A. I. Baranov and F. Wagner.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: sfacmag
8 ! !INTERFACE:
9 subroutine sfacmag
10 ! !USES:
11 use modmain
12 use modpw
13 use modtest
14 ! !DESCRIPTION:
15 ! Outputs magnetic structure factors, i.e. the Fourier transform coefficients
16 ! of each component $j$ of magnetization ${\bf m}({\bf r})$,
17 ! $$ F_j({\bf H})=\int_{\Omega}d^3r\,m_j({\bf r})e^{i{\bf H}\cdot{\bf r}}, $$
18 ! to the files {\tt SFACMAG\_j.OUT}. The lattice coordinates $(h,k,l)$ of
19 ! $\bf H$-vectors in this file are transformed by the matrix {\tt vhmat}. See
20 ! also routines {\tt zftrf} and {\tt genhvec}.
21 !
22 ! !REVISION HISTORY:
23 ! Created July 2010 (Alexey I. Baranov)
24 ! Added multiplicity of the H-vectors, Oct. 2010 (Alexey I. Baranov)
25 !EOP
26 !BOC
27 implicit none
28 ! local variables
29 integer idm,ih,iv(3)
30 real(8) v(3),a,b,r
31 character(256) fname
32 ! allocatable arrays
33 complex(8), allocatable :: zmagh(:)
34 if (.not.spinpol) return
35 ! initialise the structure factor specific variables
36 call sfacinit
37 ! generate the magnetic structure factors
38 allocate(zmagh(nhvec))
39 do idm=1,ndmag
40  call zftrf(nhvec,ivh,vhc,magmt(:,:,idm),magir(:,idm),zmagh)
41  write(fname,'("SFACMAG_",I1.1,".OUT")') idm
42  open(50,file=trim(fname),form='FORMATTED')
43  write(50,*)
44  write(50,'("h k l indices transformed by vhmat matrix:")')
45  write(50,'(3G18.10)') vhmat(:,1)
46  write(50,'(3G18.10)') vhmat(:,2)
47  write(50,'(3G18.10)') vhmat(:,3)
48  write(50,*)
49  write(50,'(" h k l multipl. |H| Re(F)&
50  & Im(F) |F|")')
51  write(50,*)
52  do ih=1,nhvec
53 ! apply transformation matrix
54  v(:)=vhmat(:,1)*dble(ivh(1,ih)) &
55  +vhmat(:,2)*dble(ivh(2,ih)) &
56  +vhmat(:,3)*dble(ivh(3,ih))
57 ! in crystallography the forward Fourier transform of real-space density is
58 ! usually done with positive phase and without 1/omega prefactor
59  a=dble(zmagh(ih))*omega
60  b=-aimag(zmagh(ih))*omega
61  r=abs(zmagh(ih))*omega
62  iv(:)=nint(v(:))
63  if ((abs(v(1)-iv(1)) <= epslat).and. &
64  (abs(v(2)-iv(2)) <= epslat).and. &
65  (abs(v(3)-iv(3)) <= epslat)) then
66 ! integer hkl
67  write(50,'(4I7,4G16.8)') iv(:),mulh(ih),hc(ih),a,b,r
68  else
69 ! non-integer hkl
70  write(50,'(3F7.2,I7,4G16.8)') v(:),mulh(ih),hc(ih),a,b,r
71  end if
72  end do
73  close(50)
74 end do
75 write(*,*)
76 write(*,'("Info(sfacmag): magnetic structure factors written to &
77  &SFACMAG_j.OUT")')
78 write(*,'(" for magnetic components j = ",3I2)') (idm,idm=1,ndmag)
79 if (ndmag == 1) then
80  write(*,'(" (this corresponds to the z-component of the magnetisation)")')
81 end if
82 write(*,*)
83 write(*,'(" Energy window : ",2G18.10)') wsfac(:)
84 ! write the structure factors to test file
85 call writetest(196,'magnetic structure factors',nv=nhvec,tol=1.d-4,zva=zmagh(:))
86 deallocate(zmagh)
87 end subroutine
88 !EOC
89 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
logical spinpol
Definition: modmain.f90:228
integer ndmag
Definition: modmain.f90:238
real(8) omega
Definition: modmain.f90:20
real(8), dimension(3, 3) vhmat
Definition: modpw.f90:32
real(8), dimension(:,:), allocatable vhc
Definition: modpw.f90:28
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
real(8), dimension(2) wsfac
Definition: modmain.f90:1108
subroutine zftrf(npv, ivp, vpc, rfmt, rfir, zfp)
Definition: zftrf.f90:10
integer, dimension(:), allocatable mulh
Definition: modpw.f90:26
real(8) epslat
Definition: modmain.f90:24
subroutine sfacinit
Definition: sfacinit.f90:7
real(8), dimension(:), allocatable hc
Definition: modpw.f90:30
Definition: modpw.f90:6
subroutine sfacmag
Definition: sfacmag.f90:10
integer, dimension(:,:), allocatable ivh
Definition: modpw.f90:24
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
integer nhvec
Definition: modpw.f90:22