38 use w90_library,
only: lib_common_type, &
54 wannier_opts, w90main, stdout, stderr, ierr)
55 type(namespace_t),
intent(in) :: namespace
56 class(space_t),
intent(in) :: space
57 class(mesh_t),
intent(in) :: mesh
58 type(ions_t),
intent(in) :: ions
59 type(kpoints_t),
intent(in) :: kpoints
60 type(states_elec_t),
intent(in) :: st
61 type(wan_opts_t),
intent(in) :: wannier_opts
62 type(lib_common_type),
intent(inout) :: w90main
63 integer,
intent(inout) :: stdout, stderr, ierr
66 integer,
allocatable :: distk(:)
77 call w90_set_option(w90main,
'kpoints', -kpoints%reduced%red_point(1:3,:))
78 call w90_set_option(w90main,
'mp_grid', kpoints%nik_axis(1:3))
79 call w90_set_option(w90main,
'num_bands', wannier_opts%num_bands)
80 call w90_set_option(w90main,
'num_wann', wannier_opts%num_wann)
84 nk = product(kpoints%nik_axis(1:3))
85 safe_allocate(distk(1:nk))
86 if (st%d%kpt%parallel)
then
88 call w90_set_comm(w90main, st%d%kpt%mpi_grp%comm%MPI_VAL)
94 ik = st%d%get_kpoint_index(iq)
95 if (distk(ik) == -1)
then
96 distk(ik) = st%d%kpt%node(iq)
97 else if (distk(ik) /= st%d%kpt%node(iq))
then
98 message(1) =
'Error: Wannier90 k-point parallelization: same k-point' // &
99 'assigned to different nodes.'
105 call w90_set_comm(w90main, st%system_grp%comm%MPI_VAL)
109 call w90_set_option(w90main,
'distk', distk)
110 safe_deallocate_a(distk)
120 if (wannier_opts%dump_inputs)
then
121 call w90_set_option(w90main,
"dump_inputs", .
true.)
125 call w90_input_setopt(w90main, wannier_opts%prefix, stdout, stderr, ierr)
130 call w90_input_reader(w90main, stdout, stderr, ierr)
132 write(
message(1),
'(a,i0)')
'Error in w90_input_reader: ', ierr
146 type(lib_common_type),
intent(inout) :: w90main
147 type(ions_t),
intent(in) :: ions
149 integer :: ia, n_atoms
150 real(real64),
allocatable :: atom_pos_frac(:,:)
151 character(len=2),
allocatable :: atom_labels(:)
155 n_atoms = ions%natoms
158 safe_allocate(atom_pos_frac(1:3, 1:n_atoms))
159 safe_allocate(atom_labels(1:n_atoms))
165 atom_pos_frac(1:3, ia) = ions%latt%cart_to_red(ions%pos(1:3, ia))
167 atom_labels(ia) = trim(ions%atom(ia)%label)
171 call w90_set_option(w90main,
'atoms_frac', atom_pos_frac)
173 call w90_set_option(w90main,
'symbols', atom_labels)
175 safe_deallocate_a(atom_pos_frac)
176 safe_deallocate_a(atom_labels)
187 type(lib_common_type),
intent(inout) :: w90main
190 character(len=80) :: filename, line
191 integer :: w90_win, io
197 filename = trim(adjustl(wannier_opts%prefix)) //
'.win'
199 message(1) =
"oct-wannier90: Parsing projections from "//filename
202 inquire(file=filename,exist=exist)
203 if (.not. exist)
then
204 message(1) =
'oct-wannier90: Cannot find specified Wannier90 win file.'
209 w90_win =
io_open(trim(filename), action=
'read')
211 read(w90_win,
'(A)', iostat=io) line
212 if (is_iostat_end(io))
exit
214 if (trim(line) ==
'begin projections')
then
216 read(w90_win,
'(A)') line
217 if (trim(line) ==
'end projections')
exit
218 call w90_set_option(w90main,
'projections', trim(line))
235 character(len=80) :: filename, line, dummy, dummy1, dummy2
236 integer :: w90_win, io
241 filename = trim(adjustl(wannier_opts%prefix)) //
'.win'
243 message(1) =
"oct-wannier90: Parsing "//filename
246 inquire(file=filename,exist=exist)
247 if (.not. exist)
then
248 message(1) =
'oct-wannier90: Cannot find specified Wannier90 win file.'
253 w90_win =
io_open(trim(filename), action=
'read')
255 read(w90_win,
'(A)', iostat=io) line
256 if (is_iostat_end(io))
exit
258 if (trim(line) ==
'begin kpoints')
then
259 message(1) =
'oct-wannier90: Warning: kpoints specified in .win file' // &
260 ' are ignored. Using kpoints from octopus input instead.'
262 else if (trim(line) ==
'begin atoms_frac')
then
263 message(1) =
'oct-wannier90: Warning: atomic positions specified in .win file' // &
264 ' are ignored. Using positions from octopus input instead.'
266 else if (trim(line) ==
'begin unit_cell_cart')
then
267 message(1) =
'oct-wannier90: Warning: unit cell vectors specified in .win file' // &
268 ' are ignored. Using unit cell vectors from octopus input instead.'
271 if (index(line,
'=') > 0)
then
272 read(line, *, iostat=io) dummy, dummy2, dummy1
273 if (dummy ==
'num_bands' .or. dummy ==
'num_wann')
then
274 message(1) =
'oct-wannier90: Warning: Option '//trim(dummy)//
' specified in .win file' // &
275 ' is ignored. Using value from octopus input instead.'
279 read(line, *, iostat=io) dummy, dummy1
280 if(dummy ==
'mp_grid')
then
281 message(1) =
'oct-wannier90: Warning: mp_grid specified in .win file' // &
282 ' is ignored. Using k-point mesh from octopus input instead.'
294 band_index, space, mesh, latt, st, kpoints, projection)
295 type(lib_common_type),
intent(in) :: w90main
298 class(
mesh_t),
intent(in) :: mesh
302 logical,
intent(in) :: exclude_list(:)
303 integer,
intent(in) :: band_index(:)
305 complex(real64),
contiguous,
intent(out) :: projection(:, :, :)
308 type(
wan_proj_t),
allocatable :: proj_input_oct(:)
318 safe_allocate(proj_input_oct(1:w90main%num_proj))
319 do iw = 1, w90main%num_proj
320 proj_input_oct(iw)%l = w90main%proj_input(iw)%l
321 proj_input_oct(iw)%m = w90main%proj_input(iw)%m
322 proj_input_oct(iw)%s = w90main%proj_input(iw)%s
323 proj_input_oct(iw)%radial = w90main%proj_input(iw)%radial
324 proj_input_oct(iw)%site = w90main%proj_input(iw)%site
325 proj_input_oct(iw)%s_qaxis = w90main%proj_input(iw)%s_qaxis
326 proj_input_oct(iw)%z = w90main%proj_input(iw)%z
327 proj_input_oct(iw)%x = w90main%proj_input(iw)%x
328 proj_input_oct(iw)%zona = w90main%proj_input(iw)%zona
332 band_index, space, mesh, latt, st, kpoints, &
333 w90main%wvfn_read%spin_channel, w90main%num_proj, &
334 proj_input_oct, projection)
336 safe_deallocate_a(proj_input_oct)
353 type(lib_common_type),
intent(in) :: w90main
354 logical,
intent(in) :: exclude_list(:)
355 integer,
intent(in) :: band_index(:)
356 class(
mesh_t),
intent(in) :: mesh
360 complex(real64),
contiguous,
intent(out) :: overlap(:, :, :, :)
365 mesh, st, ions, w90main%wvfn_read%spin_channel, &
366 w90main%kmesh_info%nnlist, w90main%kmesh_info%nncell, &
subroutine, public io_close(iunit, grp)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module defines the meshes, which are used in Octopus.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_t), public unit_angstrom
For XYZ files.
Interface module to Wannier 90 library.
subroutine wannier90lib_set_projections(w90main, wannier_opts)
Parse projections.
subroutine, public wannier90lib_init_w90main(namespace, space, mesh, ions, kpoints, st, wannier_opts, w90main, stdout, stderr, ierr)
Initialize wannier90 library data.
subroutine wannier90lib_set_atom_data(w90main, ions)
Set atom data in wannier90 library.
subroutine, public wannier90lib_create_wannier90_mmn(w90main, exclude_list, band_index, mesh, st, ions, overlap)
Kohn-Sham State Overlap Matrix.
subroutine, public wannier90lib_create_wannier90_amn(w90main, wan_opts, exclude_list, band_index, space, mesh, latt, st, kpoints, projection)
Calculate wannier90 Projection Matrix.
subroutine wannier90lib_warn_inputs(wannier_opts)
Check for options in .win file.
subroutine, public wannier_create_amn(wan_opts, exclude_list, band_index, space, mesh, latt, st, kpoints, spin_channel, num_proj, proj_input, projection)
subroutine, public wannier_create_mmn(exclude_list, band_index, mesh, st, ions, spin_channel, nnlist, nncell, overlap)
Kohn-Sham State Overlap Matrix.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.
this structure emulates the proj_type used by wannier90 library it is needed in the calculation of th...
batches of electronic states