The Elk Code
nuclei.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine nuclei
7 use modmain
8 implicit none
9 ! local variables
10 integer is,nr
11 real(8) rmn,t1
12 ! external functions
13 real(8), external :: radnucl
14 do is=1,nspecies
15  nr=nrmt(is)
16  rmn=rminsp(is)
17  t1=dble(nr-1)/log(rmt(is)/rmn)
18 ! approximate nuclear radius
19  rnucl(is)=radnucl(spzn(is))
20 ! number of radial mesh points to nuclear radius
21  nrnucl(is)=nint(t1*log(rnucl(is)/rmn))+1
22  nrnucl(is)=max(min(nrnucl(is),nr),1)
23 ! nuclear volume determined from integrating to nuclear radius; this is to
24 ! ensure integrals over the nuclear volume are correctly normalised
25  volnucl(is)=fourpi*sum(wr2mt(1:nrnucl(is),is))
26 ! Thomson radius
27  rtmsn(is)=abs(spzn(is))/solsc**2
28 ! number of radial mesh points to Thomson radius
29  nrtmsn(is)=nint(t1*log(rtmsn(is)/rmn))+1
30  nrtmsn(is)=max(min(nrtmsn(is),nr),1)
31 ! Thomson volume
32  voltmsn(is)=fourpi*sum(wr2mt(1:nrtmsn(is),is))
33 end do
34 end subroutine
35 
real(8), dimension(maxspecies) rnucl
Definition: modmain.f90:85
integer, dimension(maxspecies) nrtmsn
Definition: modmain.f90:95
real(8), dimension(maxspecies) rtmsn
Definition: modmain.f90:91
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
real(8), dimension(maxspecies) voltmsn
Definition: modmain.f90:93
real(8) solsc
Definition: modmain.f90:1253
integer, dimension(maxspecies) nrnucl
Definition: modmain.f90:89
subroutine nuclei
Definition: nuclei.f90:7
real(8), dimension(maxspecies) rminsp
Definition: modmain.f90:103
integer nspecies
Definition: modmain.f90:34
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
real(8), parameter fourpi
Definition: modmain.f90:1234
real(8), dimension(maxspecies) volnucl
Definition: modmain.f90:87
real(8), dimension(:,:), allocatable wr2mt
Definition: modmain.f90:183
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150