The Elk Code
 
Loading...
Searching...
No Matches
init3.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 S. Sharma, J. K. Dewhurst 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
6subroutine init3
7use modmain
8use modgw
9use modtddft
10use modvars
11use modmpi
12implicit none
13! local variables
14integer ig,iw
15real(8) w1,w2,t1,t2
16
17!-------------------------------------------------------------!
18! response function and perturbation theory variables !
19!-------------------------------------------------------------!
20! G-vectors for response functions
21ngrf=1
22do ig=2,ngvec
23 if (gc(ig) > gmaxrf) then
24 ngrf=ig-1
25 exit
26 end if
27end do
28ngrf=min(ngrf,ngvc)
29! write the G-vectors to file
30if (mp_mpi) call writegvecrf
31! frequencies for reponse functions
32nwrf=1
33if (allocated(wrf)) deallocate(wrf)
34if (any(task == [320,330,331])) then
36 allocate(wrf(nwrf))
37 w1=wplot(1)
38 w2=max(wplot(2),w1)
39 t1=(w2-w1)/dble(nwplot)
40 do iw=1,nwplot
41 t2=w1+t1*dble(iw-1)
42 wrf(iw)=cmplx(t2,swidth,8)
43 end do
44! set the first frequency to zero for the bootstrap functional
45 if ((fxctype(1) == 210).or.(fxctype(1) == 211)) then
46 wrf(1)=cmplx(0.d0,swidth,8)
47 end if
48else if (any(task == [600,601,610,620,640]).or.ksgwrho) then
49! GW Matsubara frequencies
50 call genwgw
51else
52 nwrf=1
53 allocate(wrf(nwrf))
54 wrf(1)=cmplx(0.d0,swidth,8)
55end if
56! write to VARIABLES.OUT
57if (wrtvars) then
58 call writevars('gmaxrf',rv=gmaxrf)
59 call writevars('ngrf',iv=ngrf)
60 call writevars('nwrf',iv=nwrf)
61 call writevars('wrf',nv=nwrf,zva=wrf)
62end if
63
64end subroutine
65
subroutine genwgw
Definition genwgw.f90:7
subroutine init3
Definition init3.f90:7
Definition modgw.f90:6
logical ksgwrho
Definition modgw.f90:38
integer ngvec
Definition modmain.f90:396
integer nwplot
Definition modmain.f90:1070
integer nwrf
Definition modmain.f90:1165
real(8) swidth
Definition modmain.f90:892
real(8), dimension(2) wplot
Definition modmain.f90:1076
integer task
Definition modmain.f90:1298
integer ngrf
Definition modmain.f90:1161
complex(8), dimension(:), allocatable wrf
Definition modmain.f90:1167
real(8) gmaxrf
Definition modmain.f90:1157
integer ngvc
Definition modmain.f90:398
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
logical mp_mpi
Definition modmpi.f90:17
integer, dimension(3) fxctype
Definition modtddft.f90:12
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition modvars.f90:16
logical wrtvars
Definition modvars.f90:9
subroutine writegvecrf