The Elk Code
 
Loading...
Searching...
No Matches
init2.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
6subroutine init2
7use modmain
8use modrdm
9use modphonon
10use modgw
11use modvars
12use modmpi
13implicit none
14! local variables
15logical lsym(48)
16integer isym,iv(3)
17real(8) boxl(3,0:3)
18real(8) ts0,ts1
19
20call timesec(ts0)
21
22!---------------------!
23! q-point set !
24!---------------------!
25! check if the system is an isolated molecule
26if (molecule) ngridq(:)=1
27! store the point group symmetries for reducing the q-point set
28if (reduceq == 0) then
29 nsymqpt=1
30 symqpt(:,:,1)=symlat(:,:,1)
31else
32 lsym(:)=.false.
33 do isym=1,nsymcrys
34 lsym(lsplsymc(isym))=.true.
35 end do
36 nsymqpt=0
37 do isym=1,nsymlat
38 if (lsym(isym)) then
40 symqpt(:,:,nsymqpt)=symlat(:,:,isym)
41 end if
42 end do
43end if
44if (any(task == [105,180,185,320,330,331])) then
45! equal k- and q-point grids for nesting function, BSE and linear-reposnse TDDFT
46 ngridq(:)=ngridk(:)
47else if (any(task == [5,300,600,601,620]).or.(xctype(1) < 0).or.ksgwrho) then
48! allow the q-point grid to be smaller than the k-point grid for OEP,
49! Hartree-Fock, RDMFT and GW
50 if (any(ngridq(:) < 1)) ngridq(:)=ngridk(:)
51else
52 ngridq(:)=abs(ngridq(:))
53end if
54! check that the q-point and k-point grids are commensurate for some tasks
55if (any(task == [5,205,240,241,300,600,601,620]).or.(xctype(1) < 0).or. &
56 ksgwrho) then
57 iv(:)=mod(ngridk(:),ngridq(:))
58 if (any(iv(:) /= 0)) then
59 write(*,*)
60 write(*,'("Error(init2): k-point grid incommensurate with q-point grid")')
61 write(*,'(" ngridk : ",3I6)') ngridk
62 write(*,'(" ngridq : ",3I6)') ngridq
63 write(*,*)
64 stop
65 end if
66end if
67! allocate the q-point arrays
68if (allocated(ivqiq)) deallocate(ivqiq)
69allocate(ivqiq(0:ngridq(1)-1,0:ngridq(2)-1,0:ngridq(3)-1))
70if (allocated(ivqiqnr)) deallocate(ivqiqnr)
71allocate(ivqiqnr(0:ngridq(1)-1,0:ngridq(2)-1,0:ngridq(3)-1))
73if (allocated(ivq)) deallocate(ivq)
74allocate(ivq(3,nqptnr))
75if (allocated(vql)) deallocate(vql)
76allocate(vql(3,nqptnr))
77if (allocated(vqc)) deallocate(vqc)
78allocate(vqc(3,nqptnr))
79if (allocated(wqpt)) deallocate(wqpt)
80allocate(wqpt(nqptnr))
81! set up the q-point box (offset should always be zero)
82boxl(:,:)=0.d0
83boxl(1,1)=1.d0; boxl(2,2)=1.d0; boxl(3,3)=1.d0
84! generate the q-point set
85! (note that the vectors vql and vqc are in the first Brillouin zone)
88! write the q-points to QPOINTS.OUT
89if (mp_mpi) call writeqpts
90! write to VARIABLES.OUT
91if (wrtvars) then
92 call writevars('nsymqpt',iv=nsymqpt)
93 call writevars('symqpt',nv=9*nsymqpt,iva=symqpt)
94 call writevars('ngridq',nv=3,iva=ngridq)
95 call writevars('nqpt',iv=nqpt)
96 call writevars('ivqiq',nv=nqptnr,iva=ivqiq)
97 call writevars('ivq',nv=3*nqptnr,iva=ivq)
98 call writevars('vql',nv=3*nqptnr,rva=vql)
99 call writevars('wqpt',nv=nqpt,rva=wqpt)
100end if
101
102!--------------------------------------------------------!
103! OEP, Hartree-Fock, RDMFT, BSE and GW variables !
104!--------------------------------------------------------!
105if (any(task == [5,180,185,300,320,330,331,600,601,620]).or. &
106 (xctype(1) < 0).or.ksgwrho) then
107! determine the regularised Coulomb Green's function for small q
108 call gengclq
109! output the Coulomb Green's function to GCLQ.OUT
110 if (mp_mpi) call writegclq
111! initialise OEP variables
112 if (xctype(1) < 0) call initoep
113end if
114if (task == 300) then
115 if (allocated(vclmat)) deallocate(vclmat)
116 allocate(vclmat(nstsv,nstsv,nkpt))
117 if (allocated(dkdc)) deallocate(dkdc)
118 allocate(dkdc(nstsv,nstsv,nkpt))
119end if
120
121call timesec(ts1)
122timeinit=timeinit+ts1-ts0
123
124end subroutine
125
subroutine gengclq
Definition gengclq.f90:10
subroutine genppts(tfbz, nsym, sym, ngridp, npptnr, epslat, bvec, boxl, nppt, ipvip, ipvipnr, ivp, vpl, vpc, wppt, wpptnr)
Definition genppts.f90:11
subroutine init2
Definition init2.f90:7
subroutine initoep
Definition initoep.f90:7
Definition modgw.f90:6
logical ksgwrho
Definition modgw.f90:38
integer reduceq
Definition modmain.f90:519
integer nsymlat
Definition modmain.f90:342
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer, dimension(:,:), allocatable ivq
Definition modmain.f90:529
integer nqptnr
Definition modmain.f90:527
integer, dimension(3) xctype
Definition modmain.f90:588
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
logical molecule
Definition modmain.f90:47
integer nsymqpt
Definition modmain.f90:521
integer, dimension(3, 3, 48) symqpt
Definition modmain.f90:523
real(8), dimension(:,:), allocatable vqc
Definition modmain.f90:547
integer nqpt
Definition modmain.f90:525
real(8) epslat
Definition modmain.f90:24
integer nkpt
Definition modmain.f90:461
real(8), dimension(:,:), allocatable vql
Definition modmain.f90:545
real(8), dimension(:), allocatable wqpt
Definition modmain.f90:549
integer, dimension(3) ngridk
Definition modmain.f90:448
integer task
Definition modmain.f90:1298
integer nstsv
Definition modmain.f90:886
integer, dimension(:,:,:), allocatable ivqiqnr
Definition modmain.f90:533
real(8) wqptnr
Definition modmain.f90:551
real(8) timeinit
Definition modmain.f90:1212
integer nsymcrys
Definition modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition modmain.f90:344
integer, dimension(3) ngridq
Definition modmain.f90:515
integer, dimension(maxsymcrys) lsplsymc
Definition modmain.f90:364
logical mp_mpi
Definition modmpi.f90:17
complex(8), dimension(:,:,:), allocatable vclmat
Definition modrdm.f90:13
complex(8), dimension(:,:,:), allocatable dkdc
Definition modrdm.f90:15
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 timesec(ts)
Definition timesec.f90:10
subroutine writegclq
Definition writegclq.f90:10
subroutine writeqpts
Definition writeqpts.f90:7