The Elk Code
setupw90.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 Manh Duc Le, 2017-18 Arsenii Gerasimov, Yaroslav Kvashnin
3 ! and Lars Nordstrom. This file is distributed under the terms of the GNU
4 ! General Public License. See the file COPYING for license details.
5 
6 subroutine setupw90
7 use modmain
8 use modw90
9 implicit none
10 ! local variables
11 integer is,ia,ias
12 real(8) real_lattice(3,3),recip_lattice(3,3)
13 ! allocatable arrays
14 integer, allocatable :: proj_l(:),proj_m(:),proj_radial(:),proj_s(:)
15 integer, allocatable :: exclude_bands(:)
16 real(8), allocatable :: atoms_cart(:,:),proj_site(:,:)
17 real(8), allocatable :: proj_z(:,:),proj_x(:,:)
18 real(8), allocatable :: proj_zona(:),proj_s_qaxis(:,:)
19 character(256), allocatable :: atom_symbols(:)
20 ! allocate global arrays
21 if (allocated(nnlist)) deallocate(nnlist)
22 allocate(nnlist(nkptnr,num_nnmax))
23 if (allocated(nncell)) deallocate(nncell)
24 allocate(nncell(3,nkptnr,num_nnmax))
25 ! allocate local arrays
26 allocate(proj_l(num_bands),proj_m(num_bands))
27 allocate(proj_radial(num_bands),proj_s(num_bands))
28 allocate(exclude_bands(num_bands))
29 allocate(atoms_cart(3,natmtot),proj_site(3,num_bands))
30 allocate(proj_z(3,num_bands),proj_x(3,num_bands))
31 allocate(proj_zona(num_bands),proj_s_qaxis(3,num_bands))
32 allocate(atom_symbols(natmtot))
33 real_lattice=br_ang*transpose(avec)
34 recip_lattice=(1.d0/br_ang)*transpose(bvec)
35 do ias=1,natmtot
36  is=idxis(ias)
37  ia=idxia(ias)
38  atom_symbols(ias)=trim(spsymb(is))
39  atoms_cart(:,ias)=br_ang*atposc(:,ia,is)
40 end do
41 call wannier_setup(seedname,ngridk,nkptnr,real_lattice,recip_lattice,vkl, &
42  num_bands,natmtot,atom_symbols,atoms_cart,.false.,spinpol,nntot,nnlist, &
43  nncell,num_bands,num_wann,proj_site,proj_l,proj_m,proj_radial,proj_z,proj_x, &
44  proj_zona,exclude_bands,proj_s,proj_s_qaxis)
45 deallocate(proj_l,proj_m)
46 deallocate(proj_radial,proj_s,exclude_bands)
47 deallocate(atoms_cart,proj_site,proj_z,proj_x)
48 deallocate(proj_zona,proj_s_qaxis,atom_symbols)
49 end subroutine
50 
integer num_bands
Definition: modw90.f90:20
logical spinpol
Definition: modmain.f90:228
integer num_wann
Definition: modw90.f90:18
integer nkptnr
Definition: modmain.f90:463
subroutine wannier_setup(seed__name, mp_grid_loc, num_kpts_loc, real_lattice_loc, recip_lattice_loc, kpt_latt_loc, num_bands_tot, num_atoms_loc, atom_symbols_loc, atoms_cart_loc, gamma_only_loc, spinors_loc, nntot_loc, nnlist_loc, nncell_loc, num_bands_loc, num_wann_loc, proj_site_loc, proj_l_loc, proj_m_loc, proj_radial_loc, proj_z_loc, proj_x_loc, proj_zona_loc, exclude_bands_loc, proj_s_loc, proj_s_qaxis_loc)
Definition: w90_stub.f90:14
subroutine setupw90
Definition: setupw90.f90:7
character(256) seedname
Definition: modw90.f90:12
integer, parameter num_nnmax
Definition: modw90.f90:30
integer, dimension(:,:,:), allocatable nncell
Definition: modw90.f90:36
Definition: modw90.f90:6
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
integer, dimension(3) ngridk
Definition: modmain.f90:448
integer nntot
Definition: modw90.f90:32
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), parameter br_ang
Definition: modmain.f90:1269
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
integer natmtot
Definition: modmain.f90:40
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
character(64), dimension(maxspecies) spsymb
Definition: modmain.f90:78
integer, dimension(:,:), allocatable nnlist
Definition: modw90.f90:34