35 use,
intrinsic :: iso_fortran_env
76 logical :: computepdos
78 integer :: ldos_nenergies = -1
79 real(real64),
allocatable :: ldos_energies(:)
87 subroutine dos_init(this, namespace, st, kpoints)
88 type(dos_t),
intent(out) :: this
89 type(namespace_t),
intent(in) :: namespace
90 type(states_elec_t),
intent(in) :: st
91 type(kpoints_t),
intent(in) :: kpoints
93 real(real64) :: evalmin, evalmax, eextend
101 npath = kpoints%nkpt_in_path()
102 if (st%nik > npath)
then
103 evalmin = minval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
104 evalmax = maxval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
106 evalmin = minval(st%eigenval(1:st%nst, 1:st%nik))
107 evalmax = maxval(st%eigenval(1:st%nst, 1:st%nik))
110 eextend = (evalmax - evalmin) /
m_four
139 call parse_variable(namespace,
'DOSEnergyPoints', 500, this%epoints)
168 this%de = (this%emax - this%emin) / (this%epoints - 1)
180 safe_allocate(this%ldos_energies(1:this%ldos_nenergies))
181 do ie = 1, this%ldos_nenergies
186 this%ldos_nenergies = -1
194 type(
dos_t),
intent(inout) :: this
198 safe_deallocate_a(this%ldos_energies)
199 this%ldos_nenergies = -1
206 subroutine dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
207 type(
dos_t),
intent(in) :: this
208 character(len=*),
intent(in) :: dir
210 class(
box_t),
intent(in) :: box
211 type(
ions_t),
target,
intent(in) :: ions
212 class(
mesh_t),
intent(in) :: mesh
216 integer :: ie, ik, ist, is, ns, maxdos, ib, ind
217 integer,
allocatable :: iunit(:)
218 real(real64) :: energy
219 real(real64),
allocatable :: tdos(:)
220 real(real64),
allocatable :: dos(:,:,:)
221 character(len=64) :: filename,format_str
224 integer :: ii, ll, mm, nn, work, norb, work2
225 integer :: ia, iorb, idim
226 real(real64) :: threshold
227 real(real64),
allocatable :: ddot(:,:,:)
228 complex(real64),
allocatable :: zdot(:,:,:)
229 real(real64),
allocatable :: weight(:,:,:)
237 if (st%d%nspin == 2) ns = 2
241 safe_allocate(dos(1:this%epoints, 1:st%nst, 0:ns-1))
242 safe_allocate(iunit(0:ns-1))
249 write(filename,
'(a,i4.4,a,i1.1,a)')
'dos-', ist,
'-', is+1,
'.dat'
251 write(filename,
'(a,i4.4,a)')
'dos-', ist,
'.dat'
253 iunit(is) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
258 do ie = 1, this%epoints
259 energy = this%emin + (ie - 1) * this%de
262 do ik = 1, st%nik, ns
264 dos(ie, ist, is) = dos(ie, ist, is) + st%kweights(ik+is) *
m_one/
m_pi * &
265 this%gamma / ((energy - st%eigenval(ist, ik+is))**2 + this%gamma**2)
280 safe_allocate(tdos(1))
285 write(filename,
'(a,i1.1,a)')
'total-dos-', is+1,
'.dat'
286 iunit(is) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
288 write(iunit(is),
'(3a)')
'# energy [', trim(
units_abbrev(
units_out%energy)),
'], total DOS (spin-resolved)'
290 do ie = 1, this%epoints
291 energy = this%emin + (ie - 1) * this%de
294 tdos(1) = tdos(1) + dos(ie, ist, is)
306 iunit(0) =
io_open(trim(dir)//
'/'//
'total-dos.dat', namespace, action=
'write')
310 do ie = 1, this%epoints
311 energy = this%emin + (ie - 1) * this%de
315 tdos(1) = tdos(1) + dos(ie, ist, is)
325 safe_deallocate_a(tdos)
329 iunit(0) =
io_open(trim(dir)//
'/'//
'total-dos-efermi.dat', namespace, action=
'write')
331 '] in a format compatible with total-dos.dat'
334 maxdos = st%smear%el_per_state * st%nst
344 if (this%computepdos)
then
347 call parse_variable(namespace,
'AOThreshold', 0.01_real64, threshold)
350 do ia = 1, ions%natoms
357 os%spec => ions%atom(ia)%species
361 do iorb = 1, os%spec%get_niwfs()
362 call os%spec%get_iwf_ilm(iorb, 1, ii, ll, mm)
363 call os%spec%get_iwf_n(iorb, 1, nn)
369 option__aotruncation__ao_full, threshold)
375 os%use_submesh = .false.
376 os%allocated_on_mesh = .
true.
377 os%spec => ions%atom(ia)%species
379 do work = 1, os%norbs
383 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
384 os, work, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
385 normalize = normalize)
388 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
389 os, work, os%radius, os%ndim, &
390 use_mesh = .not. hm%phase%is_allocated() .and. .not. os%use_submesh, &
391 normalize = normalize)
395 if (hm%phase%is_allocated())
then
397 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
400 if (.not. os%use_submesh)
then
401 safe_allocate(os%eorb_mesh(1:mesh%np, 1:os%norbs, 1:os%ndim, st%d%kpt%start:st%d%kpt%end))
402 os%eorb_mesh(:,:,:,:) =
m_zero
404 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
405 os%eorb_submesh(:,:,:,:) =
m_zero
409 os%ldorbs_eorb = max(
pad_pow2(os%sphere%np), 1)
410 if(.not. os%use_submesh) os%ldorbs_eorb = max(
pad_pow2(os%sphere%mesh%np), 1)
412 safe_allocate(os%buff_eorb(st%d%kpt%start:st%d%kpt%end))
413 do ik= st%d%kpt%start, st%d%kpt%end
419 vec_pot = hm%hm_base%uniform_vector_potential, &
420 vec_pot_var = hm%hm_base%vector_potential)
425 write(filename,
'(a, i3.3, a1, a, i1.1, a1,a)')
'pdos-at', ia,
'-', trim(os%spec%get_label()), &
428 write(filename,
'(a, i3.3, a1, a, a1,a)')
'pdos-at', ia,
'-', trim(os%spec%get_label()), &
432 iunit(0) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
435 '], projected DOS (total and orbital resolved)'
439 safe_allocate(ddot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
441 safe_allocate(zdot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
444 safe_allocate(weight(1:os%norbs,1:st%nik,1:st%nst))
445 weight(1:os%norbs,1:st%nik,1:st%nst) =
m_zero
447 do ik = st%d%kpt%start, st%d%kpt%end
448 do ib = st%group%block_start, st%group%block_end
450 if (hm%phase%is_allocated())
then
452 call st%group%psib(ib, ik)%copy_to(epsib)
453 call hm%phase%apply_to(mesh, mesh%np, .false., epsib, src = st%group%psib(ib, ik))
455 epsib => st%group%psib(ib, ik)
460 do ist = 1, st%group%psib(ib, ik)%nst
461 ind = st%group%psib(ib, ik)%ist(ist)
462 do iorb = 1, os%norbs
463 do idim = 1, st%d%dim
464 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(ddot(idim, iorb, ist))**2
470 do ist = 1, st%group%psib(ib, ik)%nst
471 ind = st%group%psib(ib, ik)%ist(ist)
472 do iorb = 1, os%norbs
473 do idim = 1, st%d%dim
474 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(zdot(idim, iorb, ist))**2
480 if (hm%phase%is_allocated())
then
481 call epsib%end(copy=.false.)
482 safe_deallocate_p(epsib)
488 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
492 safe_deallocate_a(ddot)
493 safe_deallocate_a(zdot)
496 write(format_str,
'(a,i2,a)')
'(', os%norbs+2,
'f12.6)'
497 safe_allocate(tdos(1:os%norbs))
498 do ie = 1, this%epoints
499 energy = this%emin + (ie - 1) * this%de
500 do iorb = 1, os%norbs
504 tdos(iorb) = tdos(iorb) + weight(iorb,ik,ist) *
m_one/
m_pi * &
505 this%gamma / ((energy - st%eigenval(ist, ik))**2 + this%gamma**2)
514 safe_deallocate_a(tdos)
519 safe_deallocate_a(weight)
527 safe_deallocate_a(iunit)
528 safe_deallocate_a(dos)
535 subroutine dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
536 type(
dos_t),
intent(in) :: this
537 character(len=*),
intent(in) :: dir
539 class(
box_t),
intent(in) :: box
540 type(
ions_t),
target,
intent(in) :: ions
541 class(
mesh_t),
intent(in) :: mesh
545 integer :: ie, ik, val, cond, is, ns
546 integer,
allocatable :: iunit(:)
547 real(real64) :: energy
548 real(real64) :: tjdos(1)
549 real(real64),
allocatable :: jdos(:,:)
550 character(len=64) :: filename
556 if (st%d%nspin == 2) ns = 2
560 safe_allocate(jdos(1:this%epoints, 0:ns-1))
561 safe_allocate(iunit(0:ns-1))
566 do cond = val, st%nst
567 do ik = 1, st%nik, ns
570 if(st%occ(cond, ik+is) >
m_epsilon) cycle
571 do ie = 1, this%epoints
572 energy = (ie - 1) * this%de
574 jdos(ie, is) = jdos(ie, is) + st%kweights(ik+is) *
m_one/
m_pi * &
575 this%gamma / ((energy - (st%eigenval(cond, ik+is)-st%eigenval(val, ik+is)))**2 + this%gamma**2)
583 if (st%d%nspin > 1)
then
585 write(filename,
'(a,i1.1,a)')
'total-jdos-', is+1,
'.dat'
586 iunit(is) =
io_open(trim(dir)//
'/'//trim(filename), namespace, action=
'write')
588 write(iunit(is),
'(3a)')
'# energy [', trim(
units_abbrev(
units_out%energy)),
'], total JDOS (spin-resolved)'
590 do ie = 1, this%epoints
591 energy = (ie - 1) * this%de
602 iunit(0) =
io_open(trim(dir)//
'/'//
'total-jdos.dat', namespace, action=
'write')
606 do ie = 1, this%epoints
607 energy = (ie - 1) * this%de
610 tjdos(1) = tjdos(1) + jdos(ie, is)
620 safe_deallocate_a(iunit)
621 safe_deallocate_a(jdos)
629 subroutine dos_write_ldos(this, dir, st, box, ions, mesh, how, namespace)
630 type(
dos_t),
intent(in) :: this
631 character(len=*),
intent(in) :: dir
633 class(
box_t),
intent(in) :: box
634 type(
ions_t),
target,
intent(in) :: ions
635 class(
mesh_t),
intent(in) :: mesh
636 integer(int64),
intent(in) :: how
639 integer :: ie, ik, ist, is, ns, ip, ierr
640 character(len=MAX_PATH_LEN) :: fname, name
641 real(real64) :: weight
642 real(real64),
allocatable :: ldos(:,:,:), dpsi(:,:), abs_psi2(:)
643 complex(real64),
allocatable :: zpsi(:,:)
648 if (this%ldos_nenergies < 1)
then
649 message(1) =
"LDOSEnergies must be defined for Output=ldos"
655 if (st%d%nspin == 2) ns = 2
660 safe_allocate(ldos(1:mesh%np, 1:this%ldos_nenergies, 1:ns))
663 safe_allocate(abs_psi2(1:mesh%np))
665 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
667 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
671 do ik = st%d%kpt%start, st%d%kpt%end
672 is = st%d%get_spin_index(ik)
673 do ist = st%st_start, st%st_end
679 abs_psi2(ip) = dpsi(ip, 1)**2
681 if (st%d%dim > 1)
then
682 abs_psi2(ip) = abs_psi2(ip) + dpsi(ip, 2)**2
687 abs_psi2(ip) = real(conjg(zpsi(ip, 1)) * zpsi(ip, 1), real64)
689 if (st%d%dim > 1)
then
690 abs_psi2(ip) = abs_psi2(ip) + real(conjg(zpsi(ip, 2)) * zpsi(ip, 2), real64)
695 do ie = 1, this%ldos_nenergies
696 weight = st%kweights(ik) *
m_one/
m_pi * &
697 this%gamma / ((this%ldos_energies(ie) - st%eigenval(ist, ik))**2 + this%gamma**2)
699 call lalg_axpy(mesh%np, weight, abs_psi2, ldos(:, ie, is))
704 safe_deallocate_a(dpsi)
705 safe_deallocate_a(zpsi)
706 safe_deallocate_a(abs_psi2)
708 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
713 do ie = 1, this%ldos_nenergies
714 write(name,
'(a,i3.3)')
'ldos_en-', ie
718 ldos(:, ie, is), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
722 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_coeff_batch(os, ndim, psib, dot, reduce)
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, 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.
batches of electronic states