File: genbs.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (80 lines) | stat: -rw-r--r-- 2,036 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80

! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

subroutine genbs
use modmain
use modomp
implicit none
! local variables
integer idm,is,ia,ias
integer nr,nri,nrc,nrci,npc
integer nthd,ithd
real(8) cb,t1
! allocatable arrays
real(8), allocatable :: rfmt(:,:)
if (.not.spinpol) return
! coupling constant of the external field (g_e/4c)
cb=gfacte/(4.d0*solsc)
!------------------------------------!
!     muffin-tin Kohn-Sham field     !
!------------------------------------!
call omp_hold(natmtot,nthd)
allocate(rfmt(npcmtmax,0:nthd-1))
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(ithd,is,ia,nr,nri) &
!$OMP PRIVATE(nrc,nrci,npc,idm,t1) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ias=1,natmtot
  ithd=omp_get_thread_num()
  is=idxis(ias)
  ia=idxia(ias)
  nr=nrmt(is)
  nri=nrmti(is)
  nrc=nrcmt(is)
  nrci=nrcmti(is)
  npc=npcmt(is)
! exchange-correlation magnetic field in spherical coordinates
  do idm=1,ndmag
    call rfmtftoc(nr,nri,bxcmt(:,ias,idm),rfmt(:,ithd))
    call rbsht(nrc,nrci,rfmt(:,ithd),bsmt(:,ias,idm))
  end do
! add the external magnetic field
  t1=cb*(bfcmt(3,ia,is)+bfieldc(3))
  bsmt(1:npc,ias,ndmag)=bsmt(1:npc,ias,ndmag)+t1
  if (ncmag) then
    do idm=1,2
      t1=cb*(bfcmt(idm,ia,is)+bfieldc(idm))
      bsmt(1:npc,ias,idm)=bsmt(1:npc,ias,idm)+t1
    end do
  end if
end do
!$OMP END DO
!$OMP END PARALLEL
deallocate(rfmt)
call omp_free(nthd)
!-----------------------------------------------!
!     interstitial Kohn-Sham magnetic field     !
!-----------------------------------------------!
call omp_hold(ndmag,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(t1) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do idm=1,ndmag
  if (ncmag) then
    t1=cb*bfieldc(idm)
  else
    t1=cb*bfieldc(3)
  end if
! multiply by characteristic function
  bsir(:,idm)=(bxcir(:,idm)+t1)*cfunir(:)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
return
end subroutine