corespec_io.f90 9.09 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22
!--------------------------------------------------------------------------------
! Copyright (c) 2017 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
! This file is part of FLEUR and available as free software under the conditions
! of the MIT license as expressed in the LICENSE file in more detail.
!--------------------------------------------------------------------------------

MODULE m_corespec_io

  USE m_corespec
  USE m_types
  USE m_juDFT

  IMPLICIT NONE

  CONTAINS

!===============================================================================
!
!  S U B R O U T I N E   C O R E S P E C _ I N I T
!
!-------------------------------------------------------------------------------
!
23
  SUBROUTINE corespec_init(input,atoms,coreSpecInput)
24 25 26

    IMPLICIT NONE

27 28 29
    TYPE(t_input),INTENT(IN)         :: input
    TYPE(t_atoms),INTENT(IN)         :: atoms
    TYPE(t_coreSpecInput),INTENT(IN) :: coreSpecInput
30 31 32 33 34 35 36 37 38 39 40

    INTEGER                    :: ui,i
    LOGICAL                    :: lexist

    namelist /csinp/ csi

    l_cs = .false.

! initialization of input parameters: type csi

    csi%verb = 1
41
    csi%atomType = 0
42 43 44
    csi%edge = ""
    csi%edgeidx = 0
    csi%lx = -1
Matthias Redies's avatar
Matthias Redies committed
45 46 47 48
    csi%ek0 = 0.0
    csi%emn = -2.0
    csi%emx = 20.0
    csi%ein = 0.1
49

50 51 52 53 54 55
    csi%nqphi = 12
    csi%nqr = 20
    csi%alpha_ex=0.1 
    csi%beta_ex=0.08 
    csi%I0 = 10 

56 57 58 59 60 61 62 63
! reading of input parameters from 'corespec_inp' file

    call iounit(ui)
    inquire(file="corespec_inp",exist=lexist)
    if(lexist) then
      open(ui,file="corespec_inp",status='old')
      read(ui,nml=csinp)
      close(ui)
64 65
    else if(input%l_coreSpec) THEN
      csi = coreSpecInput
66 67 68 69
    else
      return
    endif

70 71 72 73
    smeno = "corespec_init"

    write(*,'(/,a)') trim(smeno)//ssep

74 75 76 77 78 79
    IF(ANY(atoms%nlo(:).NE.0)) THEN
       WRITE(*,*) 'PLEASE NOTE:'
       WRITE(*,*) 'EELS + LOs is only for experienced users.'
       WRITE(*,*) 'I hope you know what you are doing.'
!       CALL juDFT_error("EELS + LOs not available at the moment!" ,calledby ="corespec_io")
    END IF
80 81 82 83 84

! sanity check of the input parameters; if they are not correct, program stops
! unit conversion if necessary
! if csi%verb = 1, detailed information is written to stdout

85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126
!   csv%nqphi
    if(csi%nqphi.lt.0) then
      write(*,csmsgs)  trim(smeno),"found csi%nqphi < 0 !"//csmsgerr ; stop
    endif
    csv%nqphi = csi%nqphi
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"Mesh number of angles: ","csi%nqphi = ",csv%nqphi,"will be used"
 
! csi%nqr
    if(csi%nqr.lt.0) then
      write(*,csmsgs)  trim(smeno),"found csi%nqr < 0 !"//csmsgerr ; stop
    endif
    csv%nqr = csi%nqr
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"Mesh number of radii: ","csi%nqr = ",csv%nqr,"will be used"

! csi%alpha_ex
    if(csi%alpha_ex.lt.0) then
      write(*,csmsgs)  trim(smeno),"found csi%alpha_ex < 0 !"//csmsgerr ; stop
    endif
    csv%alpha_ex = csi%alpha_ex
    if(csi%verb.eq.1) write(*,csmsgsfs)  trim(smeno),&
           &"Experimental convergence angle of incoming e-: ","csi%alpha_ex = ",csv%alpha_ex,"will be used"

! csi%beta_ex
    if(csi%beta_ex.lt.0) then
      write(*,csmsgs)  trim(smeno),"found csi%beta_ex < 0 !"//csmsgerr ; stop
    endif
    csv%beta_ex = csi%beta_ex
    if(csi%verb.eq.1) write(*,csmsgsfs)  trim(smeno),&
           &"Experimental convergence angle of outgoing e-: ","csi%beta_ex = ",csv%beta_ex,"will be used"

! csi%I0
    if(csi%I0.le.0) then
      write(*,csmsgs)  trim(smeno),"found csi%I0 <= 0 !"//csmsgerr ; stop
    endif
    csv%I0 = csi%I0
    if(csi%verb.eq.1) write(*,csmsgsfs)  trim(smeno),&
           &"Intensity of incoming electrons: ","csi%I0 = ",csv%I0,"will be used"



127 128 129
! csi%atomType
    if(csi%atomType.le.0) then
      write(*,csmsgs)  trim(smeno),"found csi%atomType <= 0 !"//csmsgerr ; stop
130
    endif
131 132
    if(csi%atomType.gt.atoms%ntype) then
      write(*,csmsgs)  trim(smeno),"found csi%atomType > atoms%ntype!"//csmsgerr ; stop
133 134
    endif
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
135
           &"atomic type: ","csi%atomType = ",csi%atomType,"will be used"
136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160

! csi%edge -> csv%nc
    csv%nc = 0
    if(len(trim(csi%edge)).ne.0) csv%nc = ichar(csi%edge)-ichar('K')+1
    if(csi%edge.eq."") then
      write(*,csmsgs) trim(smeno),"found empty csi%edge !"//csmsgerr ; stop
    endif
    if(csv%nc.lt.1.or.csv%nc.gt.6) then
      write(*,csmsgs) trim(smeno),&
           &"specify csi%edge as one of {K,L,M,N,O,P} !"//csmsgerr ; stop
    endif
    if(csi%verb.eq.1) write(*,csmsgsss)  trim(smeno),&
           &"edge: ","csi%edge = ",csi%edge,"will be used"
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"main quantum no.: ","csv%nc = ",csv%nc,"will be used"

! csi%edgeidx(:) -> csv%nljc
    if(maxval(csi%edgeidx).gt.(2*csv%nc-1)) then
      write(*,csmsgs) trim(smeno),&
         &"found csi%edgeidx > 2*csv%nc-1 !"//csmsgerr ; stop
    endif
    if(count(csi%edgeidx.gt.0).gt.(2*csv%nc-1)) then
      write(*,csmsgs) trim(smeno),&
         &"found more than 2*csv%nc-1 of csi%edgeidx > 0 !"//csmsgerr ; stop
    endif
161
    if((csv%nc-1)**2+maxval(csi%edgeidx).gt.atoms%ncst(csi%atomType)) then
162
      write(*,csmsgs) trim(smeno),&
163
         &"found (csv%nc-1)^2+maxval(csi%edgeidx) > atoms%ncst(csi%atomType)!"//csmsgerr
164 165 166 167 168 169 170 171 172 173 174 175 176
      stop
    endif
    csv%nljc = count(csi%edgeidx.gt.0)
    if(.not.allocated(csv%lc)) allocate(csv%lc(csv%nljc))
    csv%lc = (/(edgel(csi%edgeidx(i)), i = 1,csv%nljc)/)
!!$    print*,csv%lc
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"nljc edge lines: ","csv%nljc = ",csv%nljc,"will be used"

! csi%lx
    if(csi%lx.lt.0) then
      write(*,csmsgs)  trim(smeno),"found csi%lx < 0 !"//csmsgerr ; stop
    endif
177
    if(csi%lx.gt.atoms%lmax(csi%atomType)) then
178
      write(*,csmsgs)  trim(smeno),&
179
           &"found csi%lx > atoms%lmax(csi%atomType)!"//csmsgerr ; stop
180 181 182 183 184
    endif
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"maximum l: ","csi%lx = ",csi%lx,"will be used"

! csi%ek0
Matthias Redies's avatar
Matthias Redies committed
185
    if(csi%ek0.le.0.0) then
186 187
      write(*,csmsgs)  trim(smeno),"found csi%ek0 <= 0.0 !"//csmsgerr ; stop
    endif
Matthias Redies's avatar
Matthias Redies committed
188 189 190
    csi%ek0 = csi%ek0*1000.0 ! conversion from keV to eV
    csv%gamma = 1.0+csi%ek0/mec2
    csv%beta = sqrt(1.0-1.0/(csv%gamma**2))
191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206
    if(csi%verb.eq.1) then
      write(*,csmsgses)  trim(smeno),&
           &"kinetic energy of incoming electrons: ","csi%ek0 = ",csi%ek0,&
           &"eV will be used"
      write(*,csmsgses)  trim(smeno),&
           &"Lorentz factor (gamma): ","csv%gamma = ",csv%gamma,&
           &"will be used"
      write(*,csmsgses)  trim(smeno),&
           &"v/c factor (beta): ","csv%beta = ",csv%beta,&
           &"will be used"
    endif

! csi%emn csi%emx csi%ein -> csv%nex csv%egrid(0:csv%nex)
    if(csi%emn.gt.csi%emx) then
      write(*,csmsgs)  trim(smeno),"found csi%emn > csi%emx !"//csmsgerr ; stop
    endif
Matthias Redies's avatar
Matthias Redies committed
207 208
    if(csi%ein.le.0.0) then
      write(*,csmsgs)  trim(smeno),"found csi%ein <= 0.0 !"//csmsgerr ; stop
209 210 211 212 213 214 215 216 217 218
    endif
    if(((csi%emx-csi%emn)/csi%ein)-int((csi%emx-csi%emn)/csi%ein).ne.0) then
      write(*,csmsgs)  trim(smeno),&
           &"found non-integer (csi%emx-csi%emn)/csi%ein !"//csmsgerr ; stop
    endif
    csv%nex = int((csi%emx-csi%emn)/csi%ein)
    if(.not.allocated(csv%egrid)) allocate(csv%egrid(0:csv%nex))
    csv%egrid = (/(csi%emn+csi%ein*dble(i), i = 0,csv%nex)/)
    csv%nen = 0
!!$    do i = 0,csv%nex
Matthias Redies's avatar
Matthias Redies committed
219
!!$      if(csv%egrid(i).ge.0.0) then
220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242
!!$        csv%nen = i
!!$        exit
!!$      endif
!!$    enddo

!!$    print*,csv%egrid
!!$    print*,csv%nen
    if(csi%verb.eq.1) write(*,csmsgsfs)  trim(smeno),&
           &"energy spectrum lower bound: ","csi%emn = ",csi%emn,&
           &"eV will be used"
    if(csi%verb.eq.1) write(*,csmsgsfs)  trim(smeno),&
           &"energy spectrum upper bound: ","csi%emx = ",csi%emx,&
           &"eV will be used"
    if(csi%verb.eq.1) write(*,csmsgsfs)  trim(smeno),&
           &"energy spectrum increment: ","csi%ein = ",csi%ein,"eV will be used"
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"no. of energy spectrum grid points: ","csv%nex = 0 : ",&
           &csv%nex,"will be used"
    if(csi%verb.eq.1) write(*,csmsgsis)  trim(smeno),&
           &"minimum index for which egrid >= 0: ","csv%nen = ",&
           &csv%nen,"will be used"

    if(.not.allocated(csv%eedge)) allocate(csv%eedge(csv%nljc))
Matthias Redies's avatar
Matthias Redies committed
243
    csv%eedge = 0.0
244
    if(.not.allocated(csv%occ)) allocate(csv%occ(csv%nljc))
Matthias Redies's avatar
Matthias Redies committed
245
    csv%occ = 0.0
246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283

    l_cs = .true.

    if(csi%verb.eq.1) write(*,*) ""

  end subroutine corespec_init
!
!===============================================================================
!===============================================================================
!
!  S U B R O U T I N E   I O U N I T
!
!-------------------------------------------------------------------------------
!
  subroutine iounit(unit)
!
! Assignes unused integer to I/O unit
!
    implicit none
!
    integer, intent(out) :: unit
!
    integer :: i
    logical :: fopen
!
    i = 9
    fopen = .true.
    do while(fopen)
      i = i+1
      inquire(unit=i,opened=fopen)
    enddo
!
    unit = i
!
  end subroutine iounit
!
!===============================================================================

284 285


286
end module m_corespec_io