The Elk Code
 
Loading...
Searching...
No Matches
genbs.f90
Go to the documentation of this file.
1
2! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6subroutine genbs
7use modmain
8use modomp
9implicit none
10! local variables
11integer idm,is,ia,ias
12integer nrc,nrci,npc,nthd
13real(8) cb,t1
14if (.not.spinpol) return
15! coupling constant of the external field (g_e/4c)
16cb=gfacte/(4.d0*solsc)
17!------------------------------------!
18! muffin-tin Kohn-Sham field !
19!------------------------------------!
20call holdthd(natmtot+ndmag,nthd)
21!$OMP PARALLEL DEFAULT(SHARED) &
22!$OMP PRIVATE(ias,is,ia,nrc) &
23!$OMP PRIVATE(nrci,npc,idm,t1) &
24!$OMP NUM_THREADS(nthd)
25!$OMP DO SCHEDULE(DYNAMIC)
26do ias=1,natmtot
27 is=idxis(ias)
28 ia=idxia(ias)
29 nrc=nrcmt(is)
30 nrci=nrcmti(is)
31 npc=npcmt(is)
32! exchange-correlation magnetic field in spherical coordinates
33 do idm=1,ndmag
34 call rfmtftoc(nrc,nrci,bxcmt(:,ias,idm),bsmt(:,ias,idm))
35 call rbshtip(nrc,nrci,bsmt(:,ias,idm))
36 end do
37! add the external magnetic field
38 t1=cb*(bfcmt(3,ia,is)+bfieldc(3))
39 bsmt(1:npc,ias,ndmag)=bsmt(1:npc,ias,ndmag)+t1
40 if (ncmag) then
41 do idm=1,2
42 t1=cb*(bfcmt(idm,ia,is)+bfieldc(idm))
43 bsmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)+t1
44 end do
45 end if
46end do
47!$OMP END DO NOWAIT
48!-----------------------------------------------!
49! interstitial Kohn-Sham magnetic field !
50!-----------------------------------------------!
51!$OMP DO SCHEDULE(DYNAMIC)
52do idm=1,ndmag
53 if (ncmag) then
54 t1=cb*bfieldc(idm)
55 else
56 t1=cb*bfieldc(3)
57 end if
58 bsir(1:ngtot,idm)=bxcir(1:ngtot,idm)+t1
59end do
60!$OMP END DO
61!$OMP END PARALLEL
62call freethd(nthd)
63! add the magnetic dipole field if required
64if (tbdip) call bdipole
65! multiply interstitial part by characteristic function and store on coarse grid
66do idm=1,ndmag
67 call rfirftoc(bsir(:,idm),bsirc(:,idm))
68end do
69end subroutine
70
subroutine bdipole
Definition bdipole.f90:10
subroutine genbs
Definition genbs.f90:7
real(8), parameter gfacte
Definition modmain.f90:1276
real(8), dimension(3) bfieldc
Definition modmain.f90:269
logical tbdip
Definition modmain.f90:643
real(8), dimension(:,:,:), allocatable bxcmt
Definition modmain.f90:636
integer ngtot
Definition modmain.f90:390
logical ncmag
Definition modmain.f90:240
real(8), dimension(3, maxatoms, maxspecies) bfcmt
Definition modmain.f90:273
logical spinpol
Definition modmain.f90:228
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
real(8), dimension(:,:), allocatable bxcir
Definition modmain.f90:636
integer, dimension(maxatoms *maxspecies) idxia
Definition modmain.f90:45
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8), dimension(:,:), allocatable bsir
Definition modmain.f90:658
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:,:), pointer, contiguous bsirc
Definition modmain.f90:660
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition modmain.f90:656
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer ndmag
Definition modmain.f90:238
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rbshtip(nr, nri, rfmt)
Definition rbshtip.f90:7
subroutine rfirftoc(rfir, rfirc)
Definition rfirftoc.f90:7
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7