The Elk Code
checkmt.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: checkmt
8 ! !INTERFACE:
9 subroutine checkmt
10 ! !USES:
11 use modmain
12 use modmpi
13 use modvars
14 ! !DESCRIPTION:
15 ! Checks for muffin-tins which are too close together or intersecting. If any
16 ! such muffin-tins are found then the radii of their associated atomic species
17 ! are adjusted so that the minimum distance between their surfaces is
18 ! {\tt rmtdelta}.
19 !
20 ! !REVISION HISTORY:
21 ! Created May 2003 (JKD)
22 ! Modified, October 2011 (JKD)
23 !EOP
24 !BOC
25 implicit none
26 ! local variables
27 integer is,js
28 real(8) dmin,t1,t2
29 real(8) rmtp(nspecies)
30 if (nspecies < 1) return
31 ! store previous muffin-tin radii
32 rmtp(1:nspecies)=rmt(1:nspecies)
33 ! restore original muffin-tin radii read from species files if required
34 if (trmt0) rmt(1:nspecies)=rmt0(1:nspecies)
35 do
36 ! find the minimum distance between muffin-tin surfaces
37  call mtdmin(is,js,dmin)
38  if (dmin > rmtdelta-1.d-4) exit
39 ! adjust muffin-tin radii
40  t1=rmt(is)+rmt(js)
41  t2=(t1+dmin-rmtdelta)/t1
42  rmt(is)=rmt(is)*t2
43  if (is /= js) rmt(js)=rmt(js)*t2
44 end do
45 do is=1,nspecies
46  if (rmt(is) < 0.1d0) then
47  write(*,*)
48  write(*,'("Error(checkmt): muffin-tin radius too small for species ",I4,&
49  &" (",A,")")') is,trim(spsymb(is))
50  write(*,'(" Radius : ",G18.10)') rmt(is)
51  write(*,*)
52  stop
53  end if
54 ! report changed muffin-tin radii
55  t1=abs(rmt(is)-rmtp(is))
56  if (t1 > 1.d-4) then
57  if (mp_mpi) then
58  write(*,'("Info(checkmt): changed muffin-tin radius of species ",I3,&
59  &" (",A,") from ",F8.4," to ",F8.4)') is,trim(spsymb(is)),rmtp(is), &
60  rmt(is)
61  end if
62  end if
63 end do
64 ! write to VARIABLES.OUT
65 if (wrtvars) call writevars('rmt',nv=nspecies,rva=rmt)
66 end subroutine
67 !EOC
68 
logical mp_mpi
Definition: modmpi.f90:17
pure subroutine mtdmin(is, js, dmin)
Definition: mtdmin.f90:10
real(8) rmtdelta
Definition: modmain.f90:160
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
subroutine checkmt
Definition: checkmt.f90:10
real(8), dimension(maxspecies) rmt0
Definition: modmain.f90:162
logical wrtvars
Definition: modvars.f90:9
Definition: modmpi.f90:6
logical trmt0
Definition: modmain.f90:165
character(64), dimension(maxspecies) spsymb
Definition: modmain.f90:78
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition: modvars.f90:16