The Elk Code
potks.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: potks
8 ! !INTERFACE:
9 subroutine potks(txc)
10 ! !USES:
11 use modmain
12 use modxcifc
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! txc : .true. if the exchange-correlation energy density and potentials
15 ! should be calculated (in,logical)
16 ! !DESCRIPTION:
17 ! Computes the Kohn-Sham effective potential by adding together the Coulomb
18 ! and exchange-correlation potentials. Also computes the effective magnetic
19 ! field. See routines {\tt potcoul} and {\tt potxc}.
20 !
21 ! !REVISION HISTORY:
22 ! Created April 2003 (JKD)
23 !EOP
24 !BOC
25 implicit none
26 ! arguments
27 logical, intent(in) :: txc
28 ! local variables
29 integer is,ias,np
30 real(8) ts0,ts1
31 call timesec(ts0)
32 ! compute the Coulomb potential
33 call potcoul
34 ! generate the kinetic energy density for meta-GGA
35 if (any(xcgrad == [3,4,5,6])) call gentau
36 ! compute the Tran-Blaha '09 constant if required
37 if ((xctype(2) == xc_mgga_x_tb09).and.(.not.tc_tb09)) call xc_c_tb09
38 ! compute the exchange-correlation potential and fields
39 if (txc) call potxc(.true.,xctype,rhomt,rhoir,magmt,magir,taumt,tauir,exmt, &
41 ! optimised effective potential exchange potential
42 if (xctype(1) < 0) call oepmain
43 ! remove the source term of the exchange-correlation magnetic field if required
44 if (spinpol.and.nosource) call projsbf
45 ! trim the exchange-correlation potential for |G| > 2 gkmax
46 call trimrfg(vxcir)
47 ! effective potential from sum of Coulomb and exchange-correlation potentials
48 do ias=1,natmtot
49  is=idxis(ias)
50  np=npmt(is)
51  vsmt(1:np,ias)=vclmt(1:np,ias)+vxcmt(1:np,ias)
52 end do
54 ! multiply interstitial part by characteristic function and store on coarse grid
55 call rfirftoc(vsir,vsirc)
56 ! add the diamagnetic term if required
57 if (bfdmag) call potdmag
58 ! generate the effective magnetic fields if required
59 if (spinpol) call genbs
60 call timesec(ts1)
61 timepot=timepot+ts1-ts0
62 end subroutine
63 !EOC
64 
real(8), dimension(:), allocatable wxcir
Definition: modmain.f90:676
integer ngtot
Definition: modmain.f90:390
real(8), dimension(:,:), allocatable vxcmt
Definition: modmain.f90:634
integer, dimension(3) xctype
Definition: modmain.f90:588
real(8), dimension(:), allocatable ecir
Definition: modmain.f90:632
logical spinpol
Definition: modmain.f90:228
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
subroutine xc_c_tb09
Definition: xc_c_tb09.f90:7
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
real(8), dimension(:,:), allocatable ecmt
Definition: modmain.f90:632
logical nosource
Definition: modmain.f90:664
logical tc_tb09
Definition: modmain.f90:680
real(8), dimension(:), allocatable vsir
Definition: modmain.f90:651
real(8), dimension(:,:), allocatable exmt
Definition: modmain.f90:630
subroutine potxc(tsh, xctype_, rhomt_, rhoir_, magmt_, magir_, taumt_, tauir_, exmt_, exir_, ecmt_, ecir_, vxcmt_, vxcir_, bxcmt_, bxcir_, wxcmt_, wxcir_)
Definition: potxc.f90:11
real(8), dimension(:), allocatable vxcir
Definition: modmain.f90:634
subroutine potdmag
Definition: potdmag.f90:10
real(8), dimension(:,:), allocatable wxcmt
Definition: modmain.f90:676
real(8) timepot
Definition: modmain.f90:1225
subroutine potks(txc)
Definition: potks.f90:10
real(8), dimension(:,:,:), allocatable bxcmt
Definition: modmain.f90:636
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
subroutine gentau
Definition: gentau.f90:7
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
real(8), dimension(:,:), allocatable bxcir
Definition: modmain.f90:636
real(8), dimension(:), pointer, contiguous vsirc
Definition: modmain.f90:653
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine oepmain
Definition: oepmain.f90:7
subroutine timesec(ts)
Definition: timesec.f90:10
subroutine genbs
Definition: genbs.f90:7
real(8), dimension(:,:), allocatable tauir
Definition: modmain.f90:672
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
integer natmtot
Definition: modmain.f90:40
subroutine projsbf
Definition: projsbf.f90:7
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
subroutine trimrfg(rfir)
Definition: trimrfg.f90:7
subroutine potcoul
Definition: potcoul.f90:10
real(8), dimension(:,:,:), allocatable taumt
Definition: modmain.f90:672
subroutine rfirftoc(rfir, rfirc)
Definition: rfirftoc.f90:7
integer xcgrad
Definition: modmain.f90:602
real(8), dimension(:), allocatable exir
Definition: modmain.f90:630
logical bfdmag
Definition: modmain.f90:236