Octopus
dos.F90
Go to the documentation of this file.
1!! Copyright (C) 2017 N. Tancogne-Dejean
2!! Copyright (C) 2025 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23module dos_oct_m
24 use accel_oct_m
26 use box_oct_m
27 use comm_oct_m
28 use debug_oct_m
30 use global_oct_m
32 use io_oct_m
34 use ions_oct_m
35 use, intrinsic :: iso_fortran_env
38 use math_oct_m
39 use mesh_oct_m
41 use mpi_oct_m
46 use parser_oct_m
52 use types_oct_m
53 use unit_oct_m
56
57 implicit none
58
59 private
60 public :: &
61 dos_t, &
62 dos_init, &
66
67 type dos_t
68 private
69 real(real64) :: emin
70 real(real64) :: emax
71 integer :: epoints
72 real(real64) :: de
73
74 real(real64) :: gamma
75
76 logical :: computepdos
77
78 integer :: ldos_nenergies = -1
79 real(real64), allocatable :: ldos_energies(:)
80 contains
81 final :: dos_end
82 end type dos_t
83
84contains
85
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
92
93 real(real64) :: evalmin, evalmax, eextend
94 integer :: npath, ie
95 type(block_t) :: blk
96
97 push_sub(dos_init)
98
99 !The range of the dos is only calculated for physical points,
100 !without the one from a k-point path
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)))
105 else !In case we only have a path, e.g., a bandstructure calculation
106 evalmin = minval(st%eigenval(1:st%nst, 1:st%nik))
107 evalmax = maxval(st%eigenval(1:st%nst, 1:st%nik))
108 end if
109 ! we extend the energy mesh by this amount
110 eextend = (evalmax - evalmin) / m_four
111
112 !%Variable DOSEnergyMin
113 !%Type float
114 !%Section Output
115 !%Description
116 !% Lower bound for the energy mesh of the DOS.
117 !% The default is the lowest eigenvalue, minus a quarter of the total range of eigenvalues.
118 !% This is ignored for the joint density of states, and the minimal energy is always set to zero.
119 !%End
120 call parse_variable(namespace, 'DOSEnergyMin', evalmin - eextend, this%emin, units_inp%energy)
121
122 !%Variable DOSEnergyMax
123 !%Type float
124 !%Section Output
125 !%Description
126 !% Upper bound for the energy mesh of the DOS.
127 !% The default is the highest eigenvalue, plus a quarter of the total range of eigenvalues.
128 !%End
129 call parse_variable(namespace, 'DOSEnergyMax', evalmax + eextend, this%emax, units_inp%energy)
130
131 !%Variable DOSEnergyPoints
132 !%Type integer
133 !%Default 500
134 !%Section Output
135 !%Description
136 !% Determines how many energy points <tt>Octopus</tt> should use for
137 !% the DOS energy grid.
138 !%End
139 call parse_variable(namespace, 'DOSEnergyPoints', 500, this%epoints)
140
141 !%Variable DOSGamma
142 !%Type float
143 !%Default 0.008 Ha
144 !%Section Output
145 !%Description
146 !% Determines the width of the Lorentzian which is used for the DOS sum.
147 !%End
148 call parse_variable(namespace, 'DOSGamma', units_from_atomic(units_inp%energy, 0.008_real64), this%gamma)
149 this%gamma = units_to_atomic(units_inp%energy, this%gamma)
150
151 !%Variable DOSComputePDOS
152 !%Type logical
153 !%Default false
154 !%Section Output
155 !%Description
156 !% Determines if projected dos are computed or not.
157 !% At the moment, the PDOS is computed from the bare pseudo-atomic orbitals, directly taken from
158 !% the pseudopotentials. The orbitals are not orthonormalized, in order to preserve their
159 !% atomic orbitals character. As a consequence, the sum of the different PDOS does not integrate
160 !% to the total DOS.
161 !%
162 !% The radii of the orbitals are controled by the threshold defined by <tt>AOThreshold</tt>,
163 !% and the fact that they are normalized or not by <tt>AONormalize</tt>.
164 !%End
165 call parse_variable(namespace, 'DOSComputePDOS', .false., this%computepdos)
166
167 ! spacing for energy mesh
168 this%de = (this%emax - this%emin) / (this%epoints - 1)
170 !%Variable LDOSEnergies
171 !%Type block
172 !%Section Output
173 !%Description
174 !% Specifies the energies at which the LDOS is computed.
175 !%End
176 if (parse_block(global_namespace, 'LDOSEnergies', blk) == 0) then
177 ! There is one high symmetry k-point per line
178 this%ldos_nenergies = parse_block_cols(blk, 0)
179
180 safe_allocate(this%ldos_energies(1:this%ldos_nenergies))
181 do ie = 1, this%ldos_nenergies
182 call parse_block_float(blk, 0, ie-1, this%ldos_energies(ie))
183 end do
184 call parse_block_end(blk)
185 else
186 this%ldos_nenergies = -1
187 end if
188
189 pop_sub(dos_init)
190 end subroutine dos_init
191
193 subroutine dos_end(this)
194 type(dos_t), intent(inout) :: this
195
196 push_sub(dos_end)
197
198 safe_deallocate_a(this%ldos_energies)
199 this%ldos_nenergies = -1
200
201 pop_sub(dos_end)
202 end subroutine
203
204 ! ---------------------------------------------------------
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
209 type(states_elec_t), target, intent(in) :: st
210 class(box_t), intent(in) :: box
211 type(ions_t), target, intent(in) :: ions
212 class(mesh_t), intent(in) :: mesh
213 type(hamiltonian_elec_t), intent(in) :: hm
214 type(namespace_t), intent(in) :: namespace
215
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
222 logical :: normalize
223
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(:,:,:)
230 type(orbitalset_t) :: os
231 type(wfs_elec_t), pointer :: epsib
232
233 push_sub(dos_write_dos)
234
235 ! shortcuts
236 ns = 1
237 if (st%d%nspin == 2) ns = 2
238
239 if (mpi_grp_is_root(mpi_world)) then
240 ! space for state-dependent DOS
241 safe_allocate(dos(1:this%epoints, 1:st%nst, 0:ns-1))
242 safe_allocate(iunit(0:ns-1))
243
244 ! compute band/spin-resolved density of states
245 do ist = 1, st%nst
246
247 do is = 0, ns-1
248 if (ns > 1) then
249 write(filename, '(a,i5.5,a,i1.1,a)') 'dos-', ist, '-', is+1,'.dat'
250 else
251 write(filename, '(a,i5.5,a)') 'dos-', ist, '.dat'
252 end if
253 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
254 ! write header
255 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], band-resolved DOS'
256 end do
257
258 do ie = 1, this%epoints
259 energy = this%emin + (ie - 1) * this%de
260 dos(ie, ist, :) = m_zero
261 ! sum up Lorentzians
262 do ik = 1, st%nik, ns
263 do is = 0, ns-1
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)
266 end do
267 end do
268 do is = 0, ns-1
269 write(message(1), '(2es19.12)') units_from_atomic(units_out%energy, energy), &
270 units_from_atomic(unit_one / units_out%energy, dos(ie, ist, is))
271 call messages_info(1, iunit(is))
272 end do
273 end do
274
275 do is = 0, ns-1
276 call io_close(iunit(is))
277 end do
278 end do
279
280 safe_allocate(tdos(1))
281
282 ! for spin-polarized calculations also output spin-resolved tDOS
283 if (st%d%nspin == spin_polarized) then
284 do is = 0, ns-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')
287 ! write header
288 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS (spin-resolved)'
289
290 do ie = 1, this%epoints
291 energy = this%emin + (ie - 1) * this%de
292 tdos(1) = m_zero
293 do ist = 1, st%nst
294 tdos(1) = tdos(1) + dos(ie, ist, is)
295 end do
296 write(message(1), '(2es19.12)') units_from_atomic(units_out%energy, energy), &
297 units_from_atomic(unit_one / units_out%energy, tdos(1))
298 call messages_info(1, iunit(is))
299 end do
300
301 call io_close(iunit(is))
302 end do
303 end if
304
305
306 iunit(0) = io_open(trim(dir)//'/'//'total-dos.dat', namespace, action='write')
307 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS'
308
309 ! compute total density of states
310 do ie = 1, this%epoints
311 energy = this%emin + (ie - 1) * this%de
312 tdos(1) = m_zero
313 do ist = 1, st%nst
314 do is = 0, ns-1
315 tdos(1) = tdos(1) + dos(ie, ist, is)
316 end do
317 end do
318 write(message(1), '(2es19.12)') units_from_atomic(units_out%energy, energy), &
319 units_from_atomic(unit_one / units_out%energy, tdos(1))
320 call messages_info(1, iunit(0))
321 end do
322
323 call io_close(iunit(0))
324
325 safe_deallocate_a(tdos)
326
327
328 ! write Fermi file
329 iunit(0) = io_open(trim(dir)//'/'//'total-dos-efermi.dat', namespace, action='write')
330 write(message(1), '(3a)') '# Fermi energy [', trim(units_abbrev(units_out%energy)), &
331 '] in a format compatible with total-dos.dat'
332
333 ! this is the maximum that tdos can reach
334 maxdos = st%smear%el_per_state * st%nst
335
336 write(message(2), '(2es19.12)') units_from_atomic(units_out%energy, st%smear%e_fermi), m_zero
337 write(message(3), '(es19.12,i7)') units_from_atomic(units_out%energy, st%smear%e_fermi), maxdos
338
339 call messages_info(3, iunit(0))
340 call io_close(iunit(0))
341
342 end if
343
344 if (this%computepdos) then
345
346 !These variables are defined in basis_set/orbitalbasis.F90
347 call parse_variable(namespace, 'AOThreshold', 0.01_real64, threshold)
348 call parse_variable(namespace, 'AONormalize', .true., normalize)
349
350 do ia = 1, ions%natoms
351 !We first count how many orbital set we have
352 work = orbitalset_utils_count(ions%atom(ia)%species)
353
354 !We loop over the orbital sets of the atom ia
355 do norb = 1, work
356 call orbitalset_init(os)
357 os%spec => ions%atom(ia)%species
358
359 !We count the orbitals
360 work2 = 0
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)
364 if (ii == norb) then
365 os%ll = ll
366 os%nn = nn
367 os%ii = ii
368 os%radius = atomic_orbital_get_radius(os%spec, mesh, iorb, 1, &
369 option__aotruncation__ao_full, threshold)
370 work2 = work2 + 1
371 end if
372 end do
373 os%norbs = work2
374 os%ndim = 1
375 os%use_submesh = .false.
376 os%allocated_on_mesh = .true.
377 os%spec => ions%atom(ia)%species
378
379 do work = 1, os%norbs
380 ! We obtain the orbital
381 if (states_are_real(st)) then
382 call dget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
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)
386 else
387 call zget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
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)
392 end if
393 end do
394
395 if (hm%phase%is_allocated()) then
396 ! In case of complex wavefunction, we allocate the array for the phase correction
397 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
398 os%phase(:,:) = m_zero
399
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
403 else
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
406 end if
407
408 if (accel_is_enabled() .and. st%d%dim == 1) then
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)
411
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
414 call accel_create_buffer(os%buff_eorb(ik), accel_mem_read_only, type_cmplx, os%ldorbs_eorb*os%norbs)
415 end do
416 end if
417
418 call orbitalset_update_phase(os, box%dim, st%d%kpt, hm%kpoints, (st%d%ispin==spin_polarized), &
419 vec_pot = hm%hm_base%uniform_vector_potential, &
420 vec_pot_var = hm%hm_base%vector_potential)
421 else
422 if (states_are_real(st)) then
423 call dorbitalset_transfer_to_device(os, st%d%kpt, .not. os%use_submesh)
424 else
425 call zorbitalset_transfer_to_device(os, st%d%kpt, .not. os%use_submesh)
426 end if
427 end if
428
429 if (mpi_grp_is_root(mpi_world)) then
430 if (os%nn /= 0) then
431 write(filename,'(a, i4.4, a1, a, i1.1, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
432 os%nn, l_notation(os%ll), '.dat'
433 else
434 write(filename,'(a, i4.4, a1, a, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
435 l_notation(os%ll), '.dat'
436 end if
437
438 iunit(0) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
439 ! write header
440 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), &
441 '], projected DOS (total and orbital resolved)'
442 end if
443
444 if (states_are_real(st)) then
445 safe_allocate(ddot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
446 else
447 safe_allocate(zdot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
448 end if
449
450 safe_allocate(weight(1:os%norbs,1:st%nik,1:st%nst))
451 weight(1:os%norbs,1:st%nik,1:st%nst) = m_zero
452
453 do ik = st%d%kpt%start, st%d%kpt%end
454 do ib = st%group%block_start, st%group%block_end
455
456 if (hm%phase%is_allocated()) then
457 safe_allocate(epsib)
458 call st%group%psib(ib, ik)%copy_to(epsib)
459 call hm%phase%apply_to(mesh, mesh%np, .false., epsib, src = st%group%psib(ib, ik))
460 else
461 epsib => st%group%psib(ib, ik)
462 end if
463
464 if (states_are_real(st)) then
465 call dorbitalset_get_coeff_batch(os, st%d%dim, epsib, ddot(:, :, :))
466 do ist = 1, st%group%psib(ib, ik)%nst
467 ind = st%group%psib(ib, ik)%ist(ist)
468 do iorb = 1, os%norbs
469 do idim = 1, st%d%dim
470 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(ddot(idim, iorb, ist))**2
471 end do
472 end do
473 end do
474 else
475 call zorbitalset_get_coeff_batch(os, st%d%dim, epsib, zdot(:, :, :))
476 do ist = 1, st%group%psib(ib, ik)%nst
477 ind = st%group%psib(ib, ik)%ist(ist)
478 do iorb = 1, os%norbs
479 do idim = 1, st%d%dim
480 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(zdot(idim, iorb, ist))**2
481 end do
482 end do
483 end do
484 end if
485
486 if (hm%phase%is_allocated()) then
487 call epsib%end(copy=.false.)
488 safe_deallocate_p(epsib)
489 end if
490
491 end do
492 end do
493
494 if (st%parallel_in_states .or. st%d%kpt%parallel) then
495 call comm_allreduce(st%st_kpt_mpi_grp, weight)
496 end if
497
498 safe_deallocate_a(ddot)
499 safe_deallocate_a(zdot)
500
501 if (mpi_grp_is_root(mpi_world)) then
502 write(format_str,'(a,i5,a)') '(', os%norbs+2, 'es19.12)'
503 safe_allocate(tdos(1:os%norbs))
504 do ie = 1, this%epoints
505 energy = this%emin + (ie - 1) * this%de
506 do iorb = 1, os%norbs
507 tdos(iorb) = m_zero
508 do ist = 1, st%nst
509 do ik = 1, st%nik
510 tdos(iorb) = tdos(iorb) + weight(iorb,ik,ist) * m_one/m_pi * &
511 this%gamma / ((energy - st%eigenval(ist, ik))**2 + this%gamma**2)
512 end do
513 end do
514 end do
515
516 write(iunit(0), trim(format_str)) units_from_atomic(units_out%energy, energy), &
517 units_from_atomic(unit_one / units_out%energy, sum(tdos)), &
518 (units_from_atomic(unit_one / units_out%energy, tdos(iorb)), iorb=1,os%norbs)
519 end do
520 safe_deallocate_a(tdos)
521 call io_close(iunit(0))
522 end if
523
524 call orbitalset_end(os)
525 safe_deallocate_a(weight)
526
527 end do
528
529 end do
530
531 end if
532
533 safe_deallocate_a(iunit)
534 safe_deallocate_a(dos)
535
536 pop_sub(dos_write_dos)
537 end subroutine dos_write_dos
538
539 ! ---------------------------------------------------------
541 subroutine dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
542 type(dos_t), intent(in) :: this
543 character(len=*), intent(in) :: dir
544 type(states_elec_t), intent(in) :: st
545 class(box_t), intent(in) :: box
546 type(ions_t), target, intent(in) :: ions
547 class(mesh_t), intent(in) :: mesh
548 type(hamiltonian_elec_t), intent(in) :: hm
549 type(namespace_t), intent(in) :: namespace
550
551 integer :: ie, ik, val, cond, is, ns
552 integer, allocatable :: iunit(:)
553 real(real64) :: energy
554 real(real64) :: tjdos(1)
555 real(real64), allocatable :: jdos(:,:)
556 character(len=64) :: filename
557
558 push_sub(dos_write_jdos)
559
560 ! shortcuts
561 ns = 1
562 if (st%d%nspin == 2) ns = 2
563
564 if (mpi_grp_is_root(mpi_world)) then
565 ! space for state-dependent DOS
566 safe_allocate(jdos(1:this%epoints, 0:ns-1))
567 safe_allocate(iunit(0:ns-1))
568 jdos = m_zero
569
570 ! compute band/spin-resolved density of states
571 do val = 1, st%nst
572 do cond = val, st%nst
573 do ik = 1, st%nik, ns
574 do is = 0, ns-1
575 if(st%occ(val, ik+is) < m_epsilon) cycle
576 if(st%occ(cond, ik+is) > m_epsilon) cycle
577 do ie = 1, this%epoints
578 energy = (ie - 1) * this%de
579 ! sum up Lorentzians
580 jdos(ie, is) = jdos(ie, is) + st%kweights(ik+is) * m_one/m_pi * &
581 this%gamma / ((energy - (st%eigenval(cond, ik+is)-st%eigenval(val, ik+is)))**2 + this%gamma**2)
582 end do
583 end do
584 end do
585 end do
586 end do
587
588 ! for spin-polarized calculations also output spin-resolved tDOS
589 if (st%d%nspin > 1) then
590 do is = 0, ns-1
591 write(filename, '(a,i1.1,a)') 'total-jdos-', is+1,'.dat'
592 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
593 ! write header
594 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS (spin-resolved)'
595
596 do ie = 1, this%epoints
597 energy = (ie - 1) * this%de
598 write(message(1), '(2es19.12)') units_from_atomic(units_out%energy, energy), &
599 units_from_atomic(unit_one / units_out%energy, jdos(ie, is))
600 call messages_info(1, iunit(is))
601 end do
602
603 call io_close(iunit(is))
604 end do
605 end if
606
607
608 iunit(0) = io_open(trim(dir)//'/'//'total-jdos.dat', namespace, action='write')
609 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS'
610
611 ! compute total joint density of states
612 do ie = 1, this%epoints
613 energy = (ie - 1) * this%de
614 tjdos(1) = m_zero
615 do is = 0, ns-1
616 tjdos(1) = tjdos(1) + jdos(ie, is)
617 end do
618 write(message(1), '(2es19.12)') units_from_atomic(units_out%energy, energy), &
619 units_from_atomic(unit_one / units_out%energy, tjdos(1))
620 call messages_info(1, iunit(0))
621 end do
622
623 call io_close(iunit(0))
624 end if
625
626 safe_deallocate_a(iunit)
627 safe_deallocate_a(jdos)
628
629 pop_sub(dos_write_jdos)
630 end subroutine dos_write_jdos
631
632
633 ! ---------------------------------------------------------
635 subroutine dos_write_ldos(this, dir, st, box, ions, mesh, how, namespace)
636 type(dos_t), intent(in) :: this
637 character(len=*), intent(in) :: dir
638 type(states_elec_t), intent(in) :: st
639 class(box_t), intent(in) :: box
640 type(ions_t), target, intent(in) :: ions
641 class(mesh_t), intent(in) :: mesh
642 integer(int64), intent(in) :: how
643 type(namespace_t), intent(in) :: namespace
644
645 integer :: ie, ik, ist, is, ns, ip, ierr
646 character(len=MAX_PATH_LEN) :: fname, name
647 real(real64) :: weight
648 real(real64), allocatable :: ldos(:,:,:), dpsi(:,:), abs_psi2(:)
649 complex(real64), allocatable :: zpsi(:,:)
650 type(unit_t) :: fn_unit
651
652 push_sub(dos_write_ldos)
653
654 if (this%ldos_nenergies < 1) then
655 message(1) = "LDOSEnergies must be defined for Output=ldos"
656 call messages_fatal(1, namespace=namespace)
657 end if
658
659 ! shortcuts
660 ns = 1
661 if (st%d%nspin == 2) ns = 2
662
663 fn_unit = units_out%length**(-ions%space%dim) / units_out%energy
664
665 ! space for state-dependent DOS
666 safe_allocate(ldos(1:mesh%np, 1:this%ldos_nenergies, 1:ns))
667 ldos = m_zero
668
669 safe_allocate(abs_psi2(1:mesh%np))
670 if (states_are_real(st)) then
671 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
672 else
673 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
674 end if
675
676 ! compute band/spin-resolved density of states
677 do ik = st%d%kpt%start, st%d%kpt%end
678 is = st%d%get_spin_index(ik)
679 do ist = st%st_start, st%st_end
680
681 ! TODO: This will be replaced by a poin-wise batch multiplication
682 if (states_are_real(st)) then
683 call states_elec_get_state(st, mesh, ist, ik, dpsi)
684 do ip = 1, mesh%np
685 abs_psi2(ip) = dpsi(ip, 1)**2
686 if (st%d%dim > 1) then
687 abs_psi2(ip) = abs_psi2(ip) + dpsi(ip, 2)**2
688 end if
689 end do
690 else
691 call states_elec_get_state(st, mesh, ist, ik, zpsi)
692 do ip = 1, mesh%np
693 abs_psi2(ip) = real(conjg(zpsi(ip, 1)) * zpsi(ip, 1), real64)
694 if (st%d%dim > 1) then
695 abs_psi2(ip) = abs_psi2(ip) + real(conjg(zpsi(ip, 2)) * zpsi(ip, 2), real64)
696 end if
697 end do
698 end if
699
700 ! sum up Lorentzians
701 do ie = 1, this%ldos_nenergies
702 weight = st%kweights(ik) * m_one/m_pi * &
703 this%gamma / ((this%ldos_energies(ie) - st%eigenval(ist, ik))**2 + this%gamma**2)
704
705 call lalg_axpy(mesh%np, weight, abs_psi2, ldos(:, ie, is))
706 end do
707 end do
708 end do !ist
709
710 safe_deallocate_a(dpsi)
711 safe_deallocate_a(zpsi)
712 safe_deallocate_a(abs_psi2)
713
714 if (st%parallel_in_states .or. st%d%kpt%parallel) then
715 call comm_allreduce(st%st_kpt_mpi_grp, ldos)
716 end if
717
718 do is = 1, ns
719 do ie = 1, this%ldos_nenergies
720 write(name, '(a,i5.5)') 'ldos_en-', ie
721 fname = get_filename_with_spin(name, st%d%nspin, is)
722
723 call dio_function_output(how, dir, fname, namespace, ions%space, mesh, &
724 ldos(:, ie, is), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
725 end do
726 end do
727
728 safe_deallocate_a(ldos)
729
730 pop_sub(dos_write_ldos)
731 end subroutine dos_write_ldos
732
733
734end module dos_oct_m
735
736!! Local Variables:
737!! mode: f90
738!! coding: utf-8
739!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
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.
Definition: dos.F90:116
subroutine dos_end(this)
Finalizer for the dos_t object.
Definition: dos.F90:287
subroutine, public dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
Computes and output the DOS and the projected DOS (PDOS)
Definition: dos.F90:300
subroutine, public dos_write_ldos(this, dir, st, box, ions, mesh, how, namespace)
Computes and output the local DOS (LDOS)
Definition: dos.F90:729
subroutine, public dos_init(this, namespace, st, kpoints)
Initializes the dot_t object.
Definition: dos.F90:181
subroutine, public dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
Computes and output the joint DOS (LDOS)
Definition: dos.F90:635
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_one
Definition: global.F90:189
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.
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
Definition: math.F90:890
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
type(namespace_t), public global_namespace
Definition: namespace.F90:132
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:206
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:232
subroutine, public dorbitalset_transfer_to_device(os, kpt, use_mesh)
Allocate and transfer the orbitals to the device.
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...
Definition: orbitalset.F90:282
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
Definition: orbitalset.F90:586
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
subroutine, public zorbitalset_transfer_to_device(os, kpt, use_mesh)
Allocate and transfer the orbitals to the device.
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
Definition: output_low.F90:115
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.
Definition: output_low.F90:232
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
pure logical function, public states_are_real(st)
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
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
Definition: box.F90:141
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)