File: potksulr.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (72 lines) | stat: -rw-r--r-- 1,929 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

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

subroutine potksulr
use modmain
use modulr
use modomp
implicit none
! local variables
integer iq,idm,is,ias
integer npc,nthd
! allocatable arrays
real(8), allocatable :: rfmt1(:),rfmt2(:)
! compute the ultra long-range Coulomb potential
call potclulr
! compute the ultra long-range exchange-correlation potential and fields
call potxculr
! subtract the normal Kohn-Sham potential for Q=0
allocate(rfmt1(npcmtmax),rfmt2(npcmtmax))
do ias=1,natmtot
  is=idxis(ias)
  npc=npcmt(is)
  call rfmtftoc(nrmt(is),nrmti(is),vsmt(:,ias),rfmt1)
  call rbsht(nrcmt(is),nrcmti(is),rfmt1,rfmt2)
  vsqmt(1:npc,ias,1)=vsqmt(1:npc,ias,1)-rfmt2(1:npc)
end do
deallocate(rfmt1,rfmt2)
vsqir(:,1)=vsqir(:,1)-vsir(:)
! multiply vsqir by the characteristic function
call omp_hold(nqpt,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do iq=1,nqpt
  vsqir(:,iq)=vsqir(:,iq)*cfunir(:)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
if (.not.spinpol) return
! subtract the normal Kohn-Sham magnetic field for Q=0 in the muffin-tins
do idm=1,ndmag
  do ias=1,natmtot
    is=idxis(ias)
    npc=npcmt(is)
    bsqmt(1:npc,ias,idm,1)=bsqmt(1:npc,ias,idm,1)-bsmt(1:npc,ias,idm)
  end do
end do
! multiply bsqir by the characteristic function
call omp_hold(nqpt,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(idm) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do iq=1,nqpt
  do idm=1,ndmag
    bsqir(:,idm,iq)=bsqir(:,idm,iq)*cfunir(:)
  end do
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
! subtract the normal Kohn-Sham magnetic field for Q=0 in the interstitial
! (this is already multiplied by the characteristic function)
do idm=1,ndmag
  bsqir(:,idm,1)=bsqir(:,idm,1)-bsir(:,idm)
end do
return
end subroutine