The Elk Code
 
Loading...
Searching...
No Matches
writelsj.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2007 J. K. Dewhurst, S. Sharma, C. Ambrosch-Draxl and
3! F. Cricchio. This file is distributed under the terms of the GNU General
4! Public License. See the file COPYING for license details.
5
6subroutine writelsj
7use modmain
8use moddftu
9use modmpi
10use modtest
11implicit none
12! local variables
13integer kst,ik,ist
14integer is,ia,ias
15real(8) xl(3),xs(3)
16! allocatable arrays
17real(8), allocatable :: xj(:,:)
18complex(8), allocatable :: dmat(:,:,:,:,:)
19! initialise universal variables
20call init0
21call init1
22! read density and potentials from file
23call readstate
24! find the new linearisation energies
25call linengy
26! generate the APW radial functions
27call genapwfr
28! generate the local-orbital radial functions
29call genlofr
30! allocate local arrays
31allocate(xj(3,natmtot))
32allocate(dmat(lmmaxo,nspinor,lmmaxo,nspinor,natmtot))
33if (task == 15) then
34! compute total L, S and J
35! read in the occupation numbers
36 do ik=1,nkpt
37 call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
38 end do
39! generate the density matrix in each muffin-tin
40 call gendmat(.false.,.false.,0,lmaxo,lmmaxo,dmat)
41 if (mp_mpi) then
42 open(50,file='LSJ.OUT',form='FORMATTED',action='WRITE')
43 write(50,*)
44 write(50,'("Expectation values are computed only over the muffin-tin")')
45! loop over species and atoms
46 do is=1,nspecies
47 write(50,*)
48 write(50,'("Species : ",I4," (",A,")")') is,trim(spsymb(is))
49 do ia=1,natoms(is)
50 ias=idxas(ia,is)
51! calculate the expectation value of L and S
52 call dmatls(dmat(:,:,:,:,ias),xl,xs)
53! J = L + S
54 xj(:,ias)=xl(:)+xs(:)
55 write(50,'(" atom : ",I4)') ia
56 write(50,'(" L : ",3G18.10)') xl(:)
57 write(50,'(" S : ",3G18.10)') xs(:)
58 write(50,'(" J : ",3G18.10)') xj(:,ias)
59! end loop over atoms and species
60 end do
61 end do
62 close(50)
63 write(*,*)
64 write(*,'("Info(writelsj):")')
65 write(*,'(" muffin-tin L, S and J expectation values written to LSJ.OUT")')
66 end if
67! write J to test file
68 call writetest(15,'total muffin-tin angular momentum',nv=3*natmtot,tol=1.d-3,&
69 rva=xj)
70else
71! compute L, S and J for all states in kstlist
72 if (mp_mpi) then
73 open(50,file='LSJ_KST.OUT',form='FORMATTED',action='WRITE')
74 write(50,*)
75 write(50,'("Expectation values are computed only over the muffin-tin")')
76 end if
77 do kst=1,nkstlist
78 ik=kstlist(1,kst)
79 ist=kstlist(2,kst)
80 if ((ik < 1).or.(ik > nkpt)) then
81 write(*,*)
82 write(*,'("Error(writelsj): k-point out of range : ",I8)') ik
83 write(*,*)
84 stop
85 end if
86 if ((ist < 1).or.(ist > nstsv)) then
87 write(*,*)
88 write(*,'("Error(writelsj): state out of range : ",I8)') ist
89 write(*,*)
90 stop
91 end if
92! select a particular wavefunction using its occupancy
93 occsv(:,:)=0.d0
94 occsv(ist,ik)=1.d0/wkpt(ik)
95! no symmetrisation required
96 nsymcrys=1
97 eqatoms(:,:,:)=.false.
98! generate the density matrix in each muffin-tin
99 call gendmat(.false.,.false.,0,lmaxo,lmmaxo,dmat)
100 if (mp_mpi) then
101! loop over species and atoms
102 do is=1,nspecies
103 do ia=1,natoms(is)
104 ias=idxas(ia,is)
105! calculate the expectation value of L and S
106 call dmatls(dmat(:,:,:,:,ias),xl,xs)
107! J = L + S
108 xj(:,ias)=xl(:)+xs(:)
109 write(50,*)
110 write(50,'("k-point : ",I6,3G18.10)') ik,vkl(:,ik)
111 write(50,'("state : ",I6)') ist
112 write(50,'("species : ",I4," (",A,"), atom : ",I4)') is, &
113 trim(spsymb(is)),ia
114 write(50,'(" L : ",3G18.10)') xl(:)
115 write(50,'(" S : ",3G18.10)') xs(:)
116 write(50,'(" J : ",3G18.10)') xj(:,ias)
117 end do
118 end do
119 end if
120 end do
121 if (mp_mpi) then
122 close(50)
123 write(*,*)
124 write(*,'("Info(writelsj):")')
125 write(*,'(" muffin-tin L, S and J expectation values for each k-point and &
126 &state")')
127 write(*,'(" in kstlist written to LSJ_KST.OUT")')
128 end if
129! write J to test file
130 call writetest(16,'muffin-tin angular momentum for one state',nv=3*natmtot, &
131 tol=1.d-3,rva=xj)
132end if
133deallocate(xj,dmat)
134end subroutine
135
subroutine dmatls(dmat, xl, xs)
Definition dmatls.f90:7
subroutine genapwfr
Definition genapwfr.f90:10
subroutine gendmat(tspndg, tlmdg, lmin, lmax, ld, dmat)
Definition gendmat.f90:7
subroutine genlofr
Definition genlofr.f90:10
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition getoccsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine linengy
Definition linengy.f90:10
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
integer nspinor
Definition modmain.f90:267
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
character(256) filext
Definition modmain.f90:1300
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer nkstlist
Definition modmain.f90:924
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
integer task
Definition modmain.f90:1298
integer nstsv
Definition modmain.f90:886
integer lmaxo
Definition modmain.f90:201
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
logical, dimension(:,:,:), allocatable eqatoms
Definition modmain.f90:370
integer nsymcrys
Definition modmain.f90:358
integer lmmaxo
Definition modmain.f90:203
integer nspecies
Definition modmain.f90:34
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer, dimension(2, maxkst) kstlist
Definition modmain.f90:926
logical mp_mpi
Definition modmpi.f90:17
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine readstate
Definition readstate.f90:10
subroutine writelsj
Definition writelsj.f90:7