The Elk Code
 
Loading...
Searching...
No Matches
moment.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: moment
8! !INTERFACE:
9subroutine moment
10! !USES:
11use modmain
12use modtest
13! !DESCRIPTION:
14! Computes the muffin-tin, interstitial and total moments by integrating the
15! magnetisation.
16!
17! !REVISION HISTORY:
18! Created January 2005 (JKD)
19!EOP
20!BOC
21implicit none
22! local variables
23integer idm,is,ias
24! external functions
25real(8), external :: rfmtint
26! find the muffin-tin moments
27mommttot(:)=0.d0
28do idm=1,ndmag
29 do ias=1,natmtot
30 is=idxis(ias)
31 mommt(idm,ias)=rfmtint(nrmt(is),nrmti(is),wr2mt(:,is),magmt(:,ias,idm))
32 mommttot(idm)=mommttot(idm)+mommt(idm,ias)
33 end do
34end do
35! find the interstitial and total moments
36do idm=1,ndmag
37 momir(idm)=(omega/ngtot)*dot_product(magir(1:ngtot,idm),cfunir(1:ngtot))
38 momtot(idm)=mommttot(idm)+momir(idm)
39end do
40! total moment magnitude
41if (ncmag) then
42 momtotm=sqrt(momtot(1)**2+momtot(2)**2+momtot(3)**2)
43else
44 momtotm=abs(momtot(1))
45end if
46! write total moment magnitude to test file
47call writetest(450,'total moment magnitude',tol=1.d-3,rv=momtotm)
48end subroutine
49!EOC
50
subroutine moment
Definition moment.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(3) mommttot
Definition modmain.f90:746
integer ngtot
Definition modmain.f90:390
logical ncmag
Definition modmain.f90:240
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:,:), allocatable mommt
Definition modmain.f90:744
real(8) omega
Definition modmain.f90:20
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(3) momir
Definition modmain.f90:742
integer natmtot
Definition modmain.f90:40
real(8), dimension(:), allocatable cfunir
Definition modmain.f90:436
real(8), dimension(3) momtot
Definition modmain.f90:738
real(8), dimension(:,:), pointer, contiguous magir
Definition modmain.f90:616
real(8) momtotm
Definition modmain.f90:740
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
integer ndmag
Definition modmain.f90:238
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7