The Elk Code
 
Loading...
Searching...
No Matches
torque.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017 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 torque
7use modmain
8implicit none
9! local variables
10integer idm
11real(8) torq(3)
12! allocatable arrays
13real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
14! external functions
15real(8), external :: rfint
16! initialise universal variables
17call init0
18if (.not.ncmag) then
19 torq(:)=0.d0
20 goto 10
21end if
22! read magnetisation and exchange-correlation magnetic field from file
23call readstate
24! compute m(r) x B_xc(r)
25allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
26call rvfcross(magmt,magir,bxcmt,bxcir,rvfmt,rvfir)
27! integrate to find the total torque
28do idm=1,ndmag
29 torq(idm)=rfint(rvfmt(:,:,idm),rvfir(:,idm))
30end do
3110 continue
32write(*,*)
33write(*,'("Info(torque):")')
34write(*,'(" Total torque exerted by B_xc on the magnetisation :")')
35write(*,'(3G18.10)') torq
36end subroutine
37
subroutine init0
Definition init0.f90:10
real(8), dimension(:,:,:), allocatable bxcmt
Definition modmain.f90:636
integer ngtot
Definition modmain.f90:390
logical ncmag
Definition modmain.f90:240
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
real(8), dimension(:,:), allocatable bxcir
Definition modmain.f90:636
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:), pointer, contiguous magir
Definition modmain.f90:616
integer npmtmax
Definition modmain.f90:216
integer ndmag
Definition modmain.f90:238
subroutine readstate
Definition readstate.f90:10
real(8) function rfint(rfmt, rfir)
Definition rfint.f90:7
subroutine rvfcross(rvfmt1, rvfir1, rvfmt2, rvfir2, rvfmt3, rvfir3)
Definition rvfcross.f90:10
subroutine torque
Definition torque.f90:7