The Elk Code
genrmesh.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: genrmesh
8 ! !INTERFACE:
9 subroutine genrmesh
10 ! !USES:
11 use modmain
12 use modvars
13 ! !DESCRIPTION:
14 ! Generates the coarse and fine radial meshes for each atomic species in the
15 ! crystal. Also determines which points are in the inner part of the
16 ! muffin-tin using the value of {\tt fracinr}.
17 !
18 ! !REVISION HISTORY:
19 ! Created September 2002 (JKD)
20 !EOP
21 !BOC
22 implicit none
23 ! local variables
24 integer is,nr,nrc,ir,irc,l
25 real(8) t1
26 ! estimate the number of radial mesh points to infinity
27 nrspmax=1
28 do is=1,nspecies
29 ! logarithmic mesh
30  t1=log(rmaxsp(is)/rminsp(is))/log(rmt(is)/rminsp(is))
31  t1=dble(nrmt(is)-1)*t1
32  nrsp(is)=nint(t1)+1
33  nrspmax=max(nrspmax,nrsp(is))
34 end do
35 ! compute and store (R_mt)ˡ
36 if (allocated(rmtl)) deallocate(rmtl)
37 allocate(rmtl(0:lmaxo+3,nspecies))
38 do is=1,nspecies
39  do l=0,lmaxo+3
40  rmtl(l,is)=rmt(is)**l
41  end do
42 end do
43 ! generate the radial meshes
44 if (allocated(rsp)) deallocate(rsp)
45 allocate(rsp(nrspmax,nspecies))
46 if (allocated(rlmt)) deallocate(rlmt)
47 allocate(rlmt(nrmtmax,-lmaxo-1:lmaxo+2,nspecies))
48 if (allocated(wr2mt)) deallocate(wr2mt)
49 allocate(wr2mt(nrmtmax,nspecies))
50 if (allocated(wprmt)) deallocate(wprmt)
51 allocate(wprmt(4,nrmtmax,nspecies))
52 if (allocated(wcrmt)) deallocate(wcrmt)
53 allocate(wcrmt(12,nrmtmax,nspecies))
54 do is=1,nspecies
55 ! generate logarithmic radial mesh
56  t1=log(rmt(is)/rminsp(is))/dble(nrmt(is)-1)
57  do ir=1,nrsp(is)
58  rsp(ir,is)=rminsp(is)*exp(dble(ir-1)*t1)
59  end do
60 ! calculate rˡ on the fine radial mesh
61  nr=nrmt(is)
62  do l=-lmaxo-1,lmaxo+2
63  rlmt(1:nr,l,is)=rsp(1:nr,is)**l
64  end do
65 ! determine the weights for spline integration on the fine radial mesh
66  call wsplint(nr,rsp(:,is),wr2mt(:,is))
67 ! multiply by r²
68  wr2mt(1:nr,is)=wr2mt(1:nr,is)*rlmt(1:nr,2,is)
69 ! determine the weights for partial integration on fine radial mesh
70  call wsplintp(nr,rsp(:,is),wprmt(:,:,is))
71 ! determine the weights for the spline coefficients
72  call wspline(nr,rsp(:,is),wcrmt(:,:,is))
73 end do
74 ! determine the fraction of the muffin-tin radius which defines the inner part
75 if (fracinr < 0.d0) fracinr=sqrt(dble(lmmaxi)/dble(lmmaxo))
76 ! set up the coarse radial meshes and find the inner part of the muffin-tin
77 ! where the maximum angular momentum is lmaxi
78 if (allocated(rcmt)) deallocate(rcmt)
79 allocate(rcmt(nrcmtmax,nspecies))
80 if (allocated(rlcmt)) deallocate(rlcmt)
81 allocate(rlcmt(nrcmtmax,-lmaxo-1:lmaxo+2,nspecies))
82 if (allocated(wr2cmt)) deallocate(wr2cmt)
83 allocate(wr2cmt(nrcmtmax,nspecies))
84 if (allocated(wprcmt)) deallocate(wprcmt)
85 allocate(wprcmt(4,nrcmtmax,nspecies))
86 if (allocated(wcrcmt)) deallocate(wcrcmt)
87 allocate(wcrcmt(12,nrcmtmax,nspecies))
88 do is=1,nspecies
89  t1=fracinr*rmt(is)
90  nrmti(is)=1
91  nrcmti(is)=1
92  do ir=1,nrmt(is),lradstp
93  irc=(ir-1)/lradstp+1
94  rcmt(irc,is)=rsp(ir,is)
95  if (rsp(ir,is) < t1) then
96  nrmti(is)=ir
97  nrcmti(is)=irc
98  end if
99  end do
100 ! store rˡ on the coarse radial mesh
101  do l=-lmaxo-1,lmaxo+2
102  do ir=1,nrmt(is),lradstp
103  irc=(ir-1)/lradstp+1
104  rlcmt(irc,l,is)=rlmt(ir,l,is)
105  end do
106  end do
107 ! determine the weights for spline integration on the coarse radial mesh
108  nrc=nrcmt(is)
109  call wsplint(nrc,rcmt(:,is),wr2cmt(:,is))
110 ! multiply by r²
111  wr2cmt(1:nrc,is)=wr2cmt(1:nrc,is)*rlcmt(1:nrc,2,is)
112 ! determine the weights for partial integration on coarse radial mesh
113  call wsplintp(nrc,rcmt(:,is),wprcmt(:,:,is))
114 ! determine the weights for the spline coefficients
115  call wspline(nrc,rcmt(:,is),wcrcmt(:,:,is))
116 end do
117 ! write to VARIABLES.OUT
118 if (wrtvars) then
119  call writevars('nrsp',nv=nspecies,iva=nrsp)
120  call writevars('nrmt',nv=nspecies,iva=nrmt)
121  call writevars('nrmti',nv=nspecies,iva=nrmti)
122  call writevars('lradstp',iv=lradstp)
123  call writevars('nrcmt',nv=nspecies,iva=nrcmt)
124  call writevars('nrcmti',nv=nspecies,iva=nrcmti)
125  do is=1,nspecies
126  call writevars('rsp',nv=nrmt(is),rva=rsp(:,is))
127  end do
128 end if
129 end subroutine
130 !EOC
131 
real(8), dimension(:,:), allocatable rcmt
Definition: modmain.f90:177
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(maxspecies) rmaxsp
Definition: modmain.f90:105
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
integer lmaxo
Definition: modmain.f90:201
integer nrspmax
Definition: modmain.f90:109
integer nrcmtmax
Definition: modmain.f90:175
subroutine wspline(n, x, wc)
Definition: wspline.f90:7
real(8), dimension(:,:), allocatable rmtl
Definition: modmain.f90:167
subroutine genrmesh
Definition: genrmesh.f90:10
real(8), dimension(:,:,:), allocatable rlcmt
Definition: modmain.f90:181
integer lradstp
Definition: modmain.f90:171
real(8), dimension(:,:,:), allocatable wprcmt
Definition: modmain.f90:191
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
subroutine wsplintp(n, x, wp)
Definition: wsplintp.f90:7
integer, dimension(maxspecies) nrsp
Definition: modmain.f90:107
logical wrtvars
Definition: modvars.f90:9
real(8), dimension(:,:), allocatable wr2cmt
Definition: modmain.f90:189
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(:,:,:), allocatable wcrcmt
Definition: modmain.f90:193
real(8), dimension(maxspecies) rminsp
Definition: modmain.f90:103
integer nspecies
Definition: modmain.f90:34
subroutine wsplint(n, x, w)
Definition: wsplint.f90:7
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
real(8), dimension(:,:,:), allocatable wprmt
Definition: modmain.f90:185
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
real(8), dimension(:,:,:), allocatable wcrmt
Definition: modmain.f90:187
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition: modvars.f90:16
real(8) fracinr
Definition: modmain.f90:209
real(8), dimension(:,:), allocatable wr2mt
Definition: modmain.f90:183
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150