The Elk Code
readbfcr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2025 Wenhan Chen, J. K. Dewhurst and S. Sharma.
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 readbfcr
7 use modmain
8 use modulr
9 implicit none
10 ! local variables
11 integer idm,i1,i2,i3,ir
12 integer ngridq_(3),i1_,i2_,i3_
13 real(8) v(3)
14 ! automatic arrays
15 real(8) bfcr(nqpt,ndmag)
16 complex(8) zfft(nfqrz)
17 ! read the real-space external magnetic field in Cartesian coordinates from file
18 open(50,file='BFCR.OUT',form='FORMATTED')
19 read(50,*) ngridq_(1:3)
20 if (any(ngridq(1:3) /= ngridq_(1:3))) then
21  write(*,*)
22  write(*,'("Error(readbfcr): differing ngridq")')
23  write(*,'(" current : ",3I6)') ngridq
24  write(*,'(" BFCR.OUT : ",3I6)') ngridq_
25  write(*,*)
26  stop
27 end if
28 ir=0
29 do i3=1,ngridq(3)
30  do i2=1,ngridq(2)
31  do i1=1,ngridq(1)
32  ir=ir+1
33  read(50,*) i1_,i2_,i3_,v(:)
34  if ((i1 /= i1_).or.(i2 /= i2_).or.(i3 /= i3_)) then
35  write(*,*)
36  write(*,'("Error(readbfcr): differing i1, i2 or i3")')
37  write(*,'(" current : ",3I6)') i1,i2,i3
38  write(*,'(" BFCR.OUT : ",3I6)') i1_,i2_,i3_
39  write(*,*)
40  stop
41  end if
42  if (ncmag) then
43  bfcr(ir,1:3)=v(1:3)
44  else
45  bfcr(ir,1)=v(3)
46  end if
47  end do
48  end do
49 end do
50 close(50)
51 ! Fourier transform external magnetic field from real-space to Q-space
52 do idm=1,ndmag
53  call rzfftifc(3,ngridq,-1,bfcr(:,idm),zfft)
54  bfcq(idm,1:nfqrz)=zfft(1:nfqrz)
55 end do
56 end subroutine
57 
subroutine readbfcr
Definition: readbfcr.f90:7
integer, dimension(3) ngridq
Definition: modmain.f90:515
complex(8), dimension(:,:), allocatable bfcq
Definition: modulr.f90:71
subroutine rzfftifc(nd, n, sgn, r, z)
Definition: modulr.f90:6
logical ncmag
Definition: modmain.f90:240