35 use,
intrinsic :: iso_fortran_env
75 logical :: computepdos
77 integer :: ldos_nenergies = -1
78 real(real64),
allocatable :: ldos_energies(:)
86 subroutine dos_init(this, namespace, st, kpoints)
87 type(dos_t),
intent(out) :: this
88 type(namespace_t),
intent(in) :: namespace
89 type(states_elec_t),
intent(in) :: st
90 type(kpoints_t),
intent(in) :: kpoints
92 real(real64) :: evalmin, evalmax, eextend
100 npath = kpoints%nkpt_in_path()
101 if (st%nik > npath)
then
102 evalmin = minval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
103 evalmax = maxval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
105 evalmin = minval(st%eigenval(1:st%nst, 1:st%nik))
106 evalmax = maxval(st%eigenval(1:st%nst, 1:st%nik))
109 eextend = (evalmax - evalmin) /
m_four
138 call parse_variable(namespace,
'DOSEnergyPoints', 500, this%epoints)
167 this%de = (this%emax - this%emin) / (this%epoints - 1)
179 safe_allocate(this%ldos_energies(1:this%ldos_nenergies))
180 do ie = 1, this%ldos_nenergies
184 this%ldos_nenergies = -1
192 type(
dos_t),
intent(inout) :: this
196 safe_deallocate_a(this%ldos_energies)
197 this%ldos_nenergies = -1
204 subroutine dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
205 type(
dos_t),
intent(in) :: this
206 character(len=*),
intent(in) :: dir
208 class(
box_t),
intent(in) :: box
209 type(
ions_t),
target,
intent(in) :: ions
210 class(
mesh_t),
intent(in) :: mesh
214 integer :: ie, ik, ist, is, ns, maxdos
215 integer,
allocatable :: iunit(:)
216 real(real64) :: energy
217 real(real64),
allocatable :: tdos(:)
218 real(real64),
allocatable :: dos(:,:,:)
219 character(len=64) :: filename,format_str
222 integer :: ii, ll, mm, nn, work, norb, work2
223 integer :: ia, iorb, idim
224 real(real64) :: threshold
225 real(real64),
allocatable :: dpsi(:,:), ddot(:,:)
226 complex(real64),
allocatable :: zpsi(:,:), zdot(:,:)
227 real(real64),
allocatable :: weight(:,:,:)
234 if (st%d%nspin == 2) ns = 2
238 safe_allocate(dos(1:this%epoints, 1:st%nst, 0:ns-1))
239 safe_allocate(iunit(0:ns-1))
246 write(filename,
'(a,i4.4,a,i1.1,a)')
'dos-', ist,
'-', is+1,
'.dat'
248 write(filename,
'(a,i4.4,a)')
'dos-', ist,
'.dat'
250 iunit(is) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
255 do ie = 1, this%epoints
256 energy = this%emin + (ie - 1) * this%de
259 do ik = 1, st%nik, ns
261 dos(ie, ist, is) = dos(ie, ist, is) + st%kweights(ik+is) *
m_one/
m_pi * &
262 this%gamma / ((energy - st%eigenval(ist, ik+is))**2 + this%gamma**2)
277 safe_allocate(tdos(1))
282 write(filename,
'(a,i1.1,a)')
'total-dos-', is+1,
'.dat'
283 iunit(is) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
285 write(iunit(is),
'(3a)')
'# energy [', trim(
units_abbrev(
units_out%energy)),
'], total DOS (spin-resolved)'
287 do ie = 1, this%epoints
288 energy = this%emin + (ie - 1) * this%de
291 tdos(1) = tdos(1) + dos(ie, ist, is)
303 iunit(0) =
io_open(trim(dir)//
'/'//
'total-dos.dat', namespace, action=
'write')
307 do ie = 1, this%epoints
308 energy = this%emin + (ie - 1) * this%de
312 tdos(1) = tdos(1) + dos(ie, ist, is)
322 safe_deallocate_a(tdos)
326 iunit(0) =
io_open(trim(dir)//
'/'//
'total-dos-efermi.dat', namespace, action=
'write')
328 '] in a format compatible with total-dos.dat'
331 maxdos = st%smear%el_per_state * st%nst
341 if (this%computepdos)
then
344 call parse_variable(namespace,
'AOThreshold', 0.01_real64, threshold)
348 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
350 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
354 do ia = 1, ions%natoms
361 os%spec => ions%atom(ia)%species
365 do iorb = 1, os%spec%get_niwfs()
366 call os%spec%get_iwf_ilm(iorb, 1, ii, ll, mm)
367 call os%spec%get_iwf_n(iorb, 1, nn)
373 option__aotruncation__ao_full, threshold)
379 os%use_submesh = .false.
381 do work = 1, os%norbs
385 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
386 os, work, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
387 normalize = normalize)
390 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
391 os, work, os%radius, os%ndim, &
392 use_mesh = .not. hm%phase%is_allocated() .and. .not. os%use_submesh, &
393 normalize = normalize)
397 if (hm%phase%is_allocated())
then
399 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
402 if (.not. os%use_submesh)
then
403 safe_allocate(os%eorb_mesh(1:mesh%np, 1:os%norbs, 1:os%ndim, st%d%kpt%start:st%d%kpt%end))
404 os%eorb_mesh(:,:,:,:) =
m_zero
406 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
407 os%eorb_submesh(:,:,:,:) =
m_zero
411 os%ldorbs = max(
pad_pow2(os%sphere%np), 1)
412 if(.not. os%use_submesh) os%ldorbs = max(
pad_pow2(os%sphere%mesh%np), 1)
414 safe_allocate(os%buff_eorb(st%d%kpt%start:st%d%kpt%end))
415 do ik= st%d%kpt%start, st%d%kpt%end
421 vec_pot = hm%hm_base%uniform_vector_potential, &
422 vec_pot_var = hm%hm_base%vector_potential)
427 write(filename,
'(a, i3.3, a1, a, i1.1, a1,a)')
'pdos-at', ia,
'-', trim(os%spec%get_label()), &
430 write(filename,
'(a, i3.3, a1, a, a1,a)')
'pdos-at', ia,
'-', trim(os%spec%get_label()), &
434 iunit(0) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
437 '], projected DOS (total and orbital resolved)'
441 safe_allocate(ddot(1:st%d%dim,1:os%norbs))
443 safe_allocate(zdot(1:st%d%dim,1:os%norbs))
446 safe_allocate(weight(1:os%norbs,1:st%nik,1:st%nst))
447 weight(1:os%norbs,1:st%nik,1:st%nst) =
m_zero
449 do ist = st%st_start, st%st_end
450 do ik = st%d%kpt%start, st%d%kpt%end
451 if (abs(st%kweights(ik)) <=
m_epsilon) cycle
455 do iorb = 1, os%norbs
456 do idim = 1, st%d%dim
457 weight(iorb,ik,ist) = weight(iorb,ik,ist) + st%kweights(ik)*abs(ddot(idim,iorb))**2
462 if (hm%phase%is_allocated())
then
464 call hm%phase%apply_to_single(zpsi, mesh%np, st%d%dim, ik, .false.)
469 do iorb = 1, os%norbs
470 do idim = 1, st%d%dim
471 weight(iorb,ik,ist) = weight(iorb,ik,ist) + st%kweights(ik)*abs(zdot(idim,iorb))**2
478 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
482 safe_deallocate_a(ddot)
483 safe_deallocate_a(zdot)
486 write(format_str,
'(a,i2,a)')
'(', os%norbs+2,
'f12.6)'
487 safe_allocate(tdos(1:os%norbs))
488 do ie = 1, this%epoints
489 energy = this%emin + (ie - 1) * this%de
490 do iorb = 1, os%norbs
494 tdos(iorb) = tdos(iorb) + weight(iorb,ik,ist) *
m_one/
m_pi * &
495 this%gamma / ((energy - st%eigenval(ist, ik))**2 + this%gamma**2)
504 safe_deallocate_a(tdos)
509 safe_deallocate_a(weight)
515 safe_deallocate_a(dpsi)
516 safe_deallocate_a(zpsi)
519 safe_deallocate_a(iunit)
520 safe_deallocate_a(dos)
527 subroutine dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
528 type(
dos_t),
intent(in) :: this
529 character(len=*),
intent(in) :: dir
531 class(
box_t),
intent(in) :: box
532 type(
ions_t),
target,
intent(in) :: ions
533 class(
mesh_t),
intent(in) :: mesh
537 integer :: ie, ik, val, cond, is, ns
538 integer,
allocatable :: iunit(:)
539 real(real64) :: energy
540 real(real64) :: tjdos(1)
541 real(real64),
allocatable :: jdos(:,:)
542 character(len=64) :: filename
548 if (st%d%nspin == 2) ns = 2
552 safe_allocate(jdos(1:this%epoints, 0:ns-1))
553 safe_allocate(iunit(0:ns-1))
558 do cond = val, st%nst
559 do ik = 1, st%nik, ns
562 if(st%occ(cond, ik+is) >
m_epsilon) cycle
563 do ie = 1, this%epoints
564 energy = (ie - 1) * this%de
566 jdos(ie, is) = jdos(ie, is) + st%kweights(ik+is) *
m_one/
m_pi * &
567 this%gamma / ((energy - (st%eigenval(cond, ik+is)-st%eigenval(val, ik+is)))**2 + this%gamma**2)
575 if (st%d%nspin > 1)
then
577 write(filename,
'(a,i1.1,a)')
'total-jdos-', is+1,
'.dat'
578 iunit(is) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
580 write(iunit(is),
'(3a)')
'# energy [', trim(
units_abbrev(
units_out%energy)),
'], total JDOS (spin-resolved)'
582 do ie = 1, this%epoints
583 energy = (ie - 1) * this%de
594 iunit(0) =
io_open(trim(dir)//
'/'//
'total-jdos.dat', namespace, action=
'write')
598 do ie = 1, this%epoints
599 energy = (ie - 1) * this%de
602 tjdos(1) = tjdos(1) + jdos(ie, is)
612 safe_deallocate_a(iunit)
613 safe_deallocate_a(jdos)
621 subroutine dos_write_ldos(this, dir, st, box, ions, mesh, how, namespace)
622 type(
dos_t),
intent(in) :: this
623 character(len=*),
intent(in) :: dir
625 class(
box_t),
intent(in) :: box
626 type(
ions_t),
target,
intent(in) :: ions
627 class(
mesh_t),
intent(in) :: mesh
628 integer(int64),
intent(in) :: how
631 integer :: ie, ik, ist, is, ns, ip, ierr
632 character(len=MAX_PATH_LEN) :: fname, name
633 real(real64) :: weight
634 real(real64),
allocatable :: ldos(:,:,:), dpsi(:,:), abs_psi2(:)
635 complex(real64),
allocatable :: zpsi(:,:)
640 if (this%ldos_nenergies < 1)
then
641 message(1) =
"LDOSEnergies must be defined for Output=ldos"
647 if (st%d%nspin == 2) ns = 2
652 safe_allocate(ldos(1:mesh%np, 1:this%ldos_nenergies, 1:ns))
655 safe_allocate(abs_psi2(1:mesh%np))
657 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
659 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
663 do ik = st%d%kpt%start, st%d%kpt%end
664 is = st%d%get_spin_index(ik)
665 do ist = st%st_start, st%st_end
671 abs_psi2(ip) = dpsi(ip, 1)**2
673 if (st%d%dim > 1)
then
674 abs_psi2(ip) = abs_psi2(ip) + dpsi(ip, 2)**2
679 abs_psi2(ip) = real(conjg(zpsi(ip, 1)) * zpsi(ip, 1), real64)
681 if (st%d%dim > 1)
then
682 abs_psi2(ip) = abs_psi2(ip) + real(conjg(zpsi(ip, 2)) * zpsi(ip, 2), real64)
687 do ie = 1, this%ldos_nenergies
688 weight = st%kweights(ik) *
m_one/
m_pi * &
689 this%gamma / ((this%ldos_energies(ie) - st%eigenval(ist, ik))**2 + this%gamma**2)
691 call lalg_axpy(mesh%np, weight, abs_psi2, ldos(:, ie, is))
696 safe_deallocate_a(dpsi)
697 safe_deallocate_a(zpsi)
698 safe_deallocate_a(abs_psi2)
700 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
705 do ie = 1, this%ldos_nenergies
706 write(name,
'(a,i3.3)')
'ldos_en-', ie
710 ldos(:, ie, is), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
714 safe_deallocate_a(ldos)
constant times a vector plus a vector
pure logical function, public accel_is_enabled()
integer, parameter, public accel_mem_read_only
subroutine, public zget_atomic_orbital(namespace, space, latt, pos, species, mesh, sm, ii, ll, jj, os, orbind, radius, d_dim, use_mesh, normalize, index_shift)
This routine returns the atomic orbital basis – provided by the pseudopotential structure in geo.
character(len=1), dimension(0:3), parameter, public l_notation
real(real64) function, public atomic_orbital_get_radius(species, mesh, iorb, ispin, truncation, threshold)
subroutine, public dget_atomic_orbital(namespace, space, latt, pos, species, mesh, sm, ii, ll, jj, os, orbind, radius, d_dim, use_mesh, normalize, index_shift)
This routine returns the atomic orbital basis – provided by the pseudopotential structure in geo.
Module that handles computing and output of various density of states.
subroutine dos_end(this)
Finalizer for the dos_t object.
subroutine, public dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
Computes and output the DOS and the projected DOS (PDOS)
subroutine, public dos_write_ldos(this, dir, st, box, ions, mesh, how, namespace)
Computes and output the local DOS (LDOS)
subroutine, public dos_init(this, namespace, st, kpoints)
Initializes the dot_t object.
subroutine, public dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
Computes and output the joint DOS (LDOS)
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_one
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine, public io_close(iunit, grp)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module is intended to contain "only mathematical" functions and procedures.
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
This module defines the meshes, which are used in Octopus.
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)
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
type(namespace_t), public global_namespace
subroutine, public orbitalset_init(this)
subroutine, public orbitalset_end(this)
subroutine, public orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase correction to the global phase in case the orbital crosses the border of the simulato...
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
subroutine, public zorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
integer function, public orbitalset_utils_count(species, iselect)
Count the number of orbital sets we have for a given atom.
this module contains the low-level part of the output system
character(len=max_path_len) function, public get_filename_with_spin(output, nspin, spin_index)
Returns the filame as output, or output-spX is spin polarized.
integer function, public parse_block(namespace, name, blk, check_varinfo_)
pure logical function, public states_are_real(st)
type(type_t), public type_cmplx
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
class to tell whether a point is inside or outside
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.