The Elk Code
 
Loading...
Searching...
No Matches
sfacrho.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: sfacrho
8! !INTERFACE:
9subroutine sfacrho
10! !USES:
11use modmain
12use modpw
13use modtest
14! !DESCRIPTION:
15! Outputs X-ray structure factors, i.e. the Fourier transform coefficients of
16! the total electron density
17! $$ F({\bf H})=\int_{\Omega}d^3r\,\rho({\bf r})e^{i{\bf H}\cdot{\bf r}}, $$
18! to the file {\tt SFACRHO.OUT}. The lattice coordinates $(h,k,l)$ of the
19! $\bf H$-vectors in this file are transformed by the matrix {\tt vhmat}. If
20! and energy window is set using the variable {\tt wsfac}, then only those
21! states within the window are used to compute the density. See also routines
22! {\tt zftrf} and {\tt genhvec}.
23!
24! !REVISION HISTORY:
25! Created July 2010 (Alexey I. Baranov)
26! Added multiplicity of the H-vectors, Oct. 2010 (Alexey I. Baranov)
27!EOP
28!BOC
29implicit none
30! local variables
31integer ih,iv(3)
32real(8) v(3),a,b,r
33! allocatable arrays
34complex(8), allocatable :: zrhoh(:)
35! initialise the structure factor specific variables
36call sfacinit
37! calculate the density structure factors
38allocate(zrhoh(nhvec))
39call zftrf(nhvec,ivh,vhc,rhomt,rhoir,zrhoh)
40open(50,file='SFACRHO.OUT',form='FORMATTED')
41write(50,*)
42write(50,'("h k l indices transformed by vhmat matrix:")')
43write(50,'(3G18.10)') vhmat(:,1)
44write(50,'(3G18.10)') vhmat(:,2)
45write(50,'(3G18.10)') vhmat(:,3)
46write(50,*)
47write(50,'(" h k l multipl. |H| Re(F)&
48 & Im(F) |F|")')
49write(50,*)
50do ih=1,nhvec
51! apply transformation matrix
52 v(:)=vhmat(:,1)*dble(ivh(1,ih)) &
53 +vhmat(:,2)*dble(ivh(2,ih)) &
54 +vhmat(:,3)*dble(ivh(3,ih))
55! in crystallography the forward Fourier transform of real-space density is
56! usually done with positive phase and without 1/omega prefactor
57 a=dble(zrhoh(ih))*omega
58 b=-aimag(zrhoh(ih))*omega
59 r=abs(zrhoh(ih))*omega
60 iv(:)=nint(v(:))
61 if ((abs(v(1)-iv(1)) <= epslat).and. &
62 (abs(v(2)-iv(2)) <= epslat).and. &
63 (abs(v(3)-iv(3)) <= epslat)) then
64! integer hkl
65 write(50,'(4I7,4G16.8)') iv(:),mulh(ih),hc(ih),a,b,r
66 else
67! non-integer hkl
68 write(50,'(3F7.2,I7,4G16.8)') v(:),mulh(ih),hc(ih),a,b,r
69 end if
70end do
71close(50)
72write(*,*)
73write(*,'("Info(sfacrho): density structure factors written to SFACRHO.OUT")')
74write(*,*)
75write(*,'(" Energy window : ",2G18.10)') wsfac(:)
76! write the structure factors to test file
77call writetest(195,'density structure factors',nv=nhvec,tol=1.d-5,zva=zrhoh(:))
78deallocate(zrhoh)
79end subroutine
80!EOC
81
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
real(8) omega
Definition modmain.f90:20
real(8), dimension(2) wsfac
Definition modmain.f90:1105
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
Definition modpw.f90:6
real(8), dimension(3, 3) vhmat
Definition modpw.f90:32
real(8), dimension(:,:), allocatable vhc
Definition modpw.f90:28
integer, dimension(:,:), allocatable ivh
Definition modpw.f90:24
integer nhvec
Definition modpw.f90:22
real(8), dimension(:), allocatable hc
Definition modpw.f90:30
integer, dimension(:), allocatable mulh
Definition modpw.f90:26
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine sfacinit
Definition sfacinit.f90:7
subroutine sfacrho
Definition sfacrho.f90:10
subroutine zftrf(npv, ivp, vpc, rfmt, rfir, zfp)
Definition zftrf.f90:10