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