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
55
56 implicit none
57
58 private
59 public :: &
60 dos_t, &
61 dos_init, &
65
66 type dos_t
67 private
68 real(real64) :: emin
69 real(real64) :: emax
70 integer :: epoints
71 real(real64) :: de
72
73 real(real64) :: gamma
74
75 logical :: computepdos
76
77 integer :: ldos_nenergies = -1
78 real(real64), allocatable :: ldos_energies(:)
79 contains
80 final :: dos_end
81 end type dos_t
82
83contains
84
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
91
92 real(real64) :: evalmin, evalmax, eextend
93 integer :: npath, ie
94 type(block_t) :: blk
95
96 push_sub(dos_init)
97
98 !The range of the dos is only calculated for physical points,
99 !without the one from a k-point path
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)))
104 else !In case we only have a path, e.g., a bandstructure calculation
105 evalmin = minval(st%eigenval(1:st%nst, 1:st%nik))
106 evalmax = maxval(st%eigenval(1:st%nst, 1:st%nik))
107 end if
108 ! we extend the energy mesh by this amount
109 eextend = (evalmax - evalmin) / m_four
110
111 !%Variable DOSEnergyMin
112 !%Type float
113 !%Section Output
114 !%Description
115 !% Lower bound for the energy mesh of the DOS.
116 !% The default is the lowest eigenvalue, minus a quarter of the total range of eigenvalues.
117 !% This is ignored for the joint density of states, and the minimal energy is always set to zero.
118 !%End
119 call parse_variable(namespace, 'DOSEnergyMin', evalmin - eextend, this%emin, units_inp%energy)
120
121 !%Variable DOSEnergyMax
122 !%Type float
123 !%Section Output
124 !%Description
125 !% Upper bound for the energy mesh of the DOS.
126 !% The default is the highest eigenvalue, plus a quarter of the total range of eigenvalues.
127 !%End
128 call parse_variable(namespace, 'DOSEnergyMax', evalmax + eextend, this%emax, units_inp%energy)
129
130 !%Variable DOSEnergyPoints
131 !%Type integer
132 !%Default 500
133 !%Section Output
134 !%Description
135 !% Determines how many energy points <tt>Octopus</tt> should use for
136 !% the DOS energy grid.
137 !%End
138 call parse_variable(namespace, 'DOSEnergyPoints', 500, this%epoints)
139
140 !%Variable DOSGamma
141 !%Type float
142 !%Default 0.008 Ha
143 !%Section Output
144 !%Description
145 !% Determines the width of the Lorentzian which is used for the DOS sum.
146 !%End
147 call parse_variable(namespace, 'DOSGamma', units_from_atomic(units_inp%energy, 0.008_real64), this%gamma)
148 this%gamma = units_to_atomic(units_inp%energy, this%gamma)
149
150 !%Variable DOSComputePDOS
151 !%Type logical
152 !%Default false
153 !%Section Output
154 !%Description
155 !% Determines if projected dos are computed or not.
156 !% At the moment, the PDOS is computed from the bare pseudo-atomic orbitals, directly taken from
157 !% the pseudopotentials. The orbitals are not orthonormalized, in order to preserve their
158 !% atomic orbitals character. As a consequence, the sum of the different PDOS does not integrate
159 !% to the total DOS.
160 !%
161 !% The radii of the orbitals are controled by the threshold defined by <tt>AOThreshold</tt>,
162 !% and the fact that they are normalized or not by <tt>AONormalize</tt>.
163 !%End
164 call parse_variable(namespace, 'DOSComputePDOS', .false., this%computepdos)
165
166 ! spacing for energy mesh
167 this%de = (this%emax - this%emin) / (this%epoints - 1)
169 !%Variable LDOSEnergies
170 !%Type block
171 !%Section Output
172 !%Description
173 !% Specifies the energies at which the LDOS is computed.
174 !%End
175 if (parse_block(global_namespace, 'LDOSEnergies', blk) == 0) then
176 ! There is one high symmetry k-point per line
177 this%ldos_nenergies = parse_block_cols(blk, 0)
178
179 safe_allocate(this%ldos_energies(1:this%ldos_nenergies))
180 do ie = 1, this%ldos_nenergies
181 call parse_block_float(blk, 0, ie-1, this%ldos_energies(ie))
182 end do
183 else
184 this%ldos_nenergies = -1
185 end if
186
187 pop_sub(dos_init)
188 end subroutine dos_init
189
191 subroutine dos_end(this)
192 type(dos_t), intent(inout) :: this
193
194 push_sub(dos_end)
195
196 safe_deallocate_a(this%ldos_energies)
197 this%ldos_nenergies = -1
198
199 pop_sub(dos_end)
200 end subroutine
201
202 ! ---------------------------------------------------------
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
207 type(states_elec_t), intent(in) :: st
208 class(box_t), intent(in) :: box
209 type(ions_t), target, intent(in) :: ions
210 class(mesh_t), intent(in) :: mesh
211 type(hamiltonian_elec_t), intent(in) :: hm
212 type(namespace_t), intent(in) :: namespace
213
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
220 logical :: normalize
221
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(:,:,:)
228 type(orbitalset_t) :: os
229
230 push_sub(dos_write_dos)
231
232 ! shortcuts
233 ns = 1
234 if (st%d%nspin == 2) ns = 2
235
236 if (mpi_grp_is_root(mpi_world)) then
237 ! space for state-dependent DOS
238 safe_allocate(dos(1:this%epoints, 1:st%nst, 0:ns-1))
239 safe_allocate(iunit(0:ns-1))
240
241 ! compute band/spin-resolved density of states
242 do ist = 1, st%nst
243
244 do is = 0, ns-1
245 if (ns > 1) then
246 write(filename, '(a,i4.4,a,i1.1,a)') 'dos-', ist, '-', is+1,'.dat'
247 else
248 write(filename, '(a,i4.4,a)') 'dos-', ist, '.dat'
249 end if
250 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
251 ! write header
252 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], band-resolved DOS'
253 end do
254
255 do ie = 1, this%epoints
256 energy = this%emin + (ie - 1) * this%de
257 dos(ie, ist, :) = m_zero
258 ! sum up Lorentzians
259 do ik = 1, st%nik, ns
260 do is = 0, ns-1
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)
263 end do
264 end do
265 do is = 0, ns-1
266 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
267 units_from_atomic(unit_one / units_out%energy, dos(ie, ist, is))
268 call messages_info(1, iunit(is))
269 end do
270 end do
271
272 do is = 0, ns-1
273 call io_close(iunit(is))
274 end do
275 end do
276
277 safe_allocate(tdos(1))
278
279 ! for spin-polarized calculations also output spin-resolved tDOS
280 if (st%d%nspin == spin_polarized) then
281 do is = 0, ns-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')
284 ! write header
285 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS (spin-resolved)'
286
287 do ie = 1, this%epoints
288 energy = this%emin + (ie - 1) * this%de
289 tdos(1) = m_zero
290 do ist = 1, st%nst
291 tdos(1) = tdos(1) + dos(ie, ist, is)
292 end do
293 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
294 units_from_atomic(unit_one / units_out%energy, tdos(1))
295 call messages_info(1, iunit(is))
296 end do
298 call io_close(iunit(is))
299 end do
300 end if
301
302
303 iunit(0) = io_open(trim(dir)//'/'//'total-dos.dat', namespace, action='write')
304 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS'
305
306 ! compute total density of states
307 do ie = 1, this%epoints
308 energy = this%emin + (ie - 1) * this%de
309 tdos(1) = m_zero
310 do ist = 1, st%nst
311 do is = 0, ns-1
312 tdos(1) = tdos(1) + dos(ie, ist, is)
313 end do
314 end do
315 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
316 units_from_atomic(unit_one / units_out%energy, tdos(1))
317 call messages_info(1, iunit(0))
318 end do
319
320 call io_close(iunit(0))
321
322 safe_deallocate_a(tdos)
323
324
325 ! write Fermi file
326 iunit(0) = io_open(trim(dir)//'/'//'total-dos-efermi.dat', namespace, action='write')
327 write(message(1), '(3a)') '# Fermi energy [', trim(units_abbrev(units_out%energy)), &
328 '] in a format compatible with total-dos.dat'
329
330 ! this is the maximum that tdos can reach
331 maxdos = st%smear%el_per_state * st%nst
332
333 write(message(2), '(2f12.6)') units_from_atomic(units_out%energy, st%smear%e_fermi), m_zero
334 write(message(3), '(f12.6,i6)') units_from_atomic(units_out%energy, st%smear%e_fermi), maxdos
335
336 call messages_info(3, iunit(0))
337 call io_close(iunit(0))
338
339 end if
340
341 if (this%computepdos) then
342
343 !These variables are defined in basis_set/orbitalbasis.F90
344 call parse_variable(namespace, 'AOThreshold', 0.01_real64, threshold)
345 call parse_variable(namespace, 'AONormalize', .true., normalize)
346
347 if (states_are_real(st)) then
348 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
349 else
350 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
351 end if
352
353
354 do ia = 1, ions%natoms
355 !We first count how many orbital set we have
356 work = orbitalset_utils_count(ions%atom(ia)%species)
357
358 !We loop over the orbital sets of the atom ia
359 do norb = 1, work
360 call orbitalset_init(os)
361 os%spec => ions%atom(ia)%species
362
363 !We count the orbitals
364 work2 = 0
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)
368 if (ii == norb) then
369 os%ll = ll
370 os%nn = nn
371 os%ii = ii
372 os%radius = atomic_orbital_get_radius(os%spec, mesh, iorb, 1, &
373 option__aotruncation__ao_full, threshold)
374 work2 = work2 + 1
375 end if
376 end do
377 os%norbs = work2
378 os%ndim = 1
379 os%use_submesh = .false.
380
381 do work = 1, os%norbs
382 ! We obtain the orbital
383 if (states_are_real(st)) then
384 call dget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
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)
388 else
389 call zget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
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)
394 end if
395 end do
396
397 if (hm%phase%is_allocated()) then
398 ! In case of complex wavefunction, we allocate the array for the phase correction
399 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
400 os%phase(:,:) = m_zero
401
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
405 else
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
408 end if
409
410 if (accel_is_enabled() .and. st%d%dim == 1) then
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)
413
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
416 call accel_create_buffer(os%buff_eorb(ik), accel_mem_read_only, type_cmplx, os%ldorbs*os%norbs)
417 end do
418 end if
419
420 call orbitalset_update_phase(os, box%dim, st%d%kpt, hm%kpoints, (st%d%ispin==spin_polarized), &
421 vec_pot = hm%hm_base%uniform_vector_potential, &
422 vec_pot_var = hm%hm_base%vector_potential)
423 end if
424
425 if (mpi_grp_is_root(mpi_world)) then
426 if (os%nn /= 0) then
427 write(filename,'(a, i3.3, a1, a, i1.1, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
428 os%nn, l_notation(os%ll), '.dat'
429 else
430 write(filename,'(a, i3.3, a1, a, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
431 l_notation(os%ll), '.dat'
432 end if
433
434 iunit(0) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
435 ! write header
436 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), &
437 '], projected DOS (total and orbital resolved)'
438 end if
439
440 if (states_are_real(st)) then
441 safe_allocate(ddot(1:st%d%dim,1:os%norbs))
442 else
443 safe_allocate(zdot(1:st%d%dim,1:os%norbs))
444 end if
445
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
448
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
452 if (states_are_real(st)) then
453 call states_elec_get_state(st, mesh, ist, ik, dpsi)
454 call dorbitalset_get_coefficients(os, st%d%dim, dpsi, ik, .false., ddot(:,1:os%norbs))
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
458 end do
459 end do
460 else
461 call states_elec_get_state(st, mesh, ist, ik, zpsi)
462 if (hm%phase%is_allocated()) then
463 ! Apply the phase that contains both the k-point and vector-potential terms.
464 call hm%phase%apply_to_single(zpsi, mesh%np, st%d%dim, ik, .false.)
465 end if
466 call zorbitalset_get_coefficients(os, st%d%dim, zpsi, ik, hm%phase%is_allocated(), &
467 zdot(:,1:os%norbs))
468
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
472 end do
473 end do
474 end if
475 end do
476 end do
477
478 if (st%parallel_in_states .or. st%d%kpt%parallel) then
479 call comm_allreduce(st%st_kpt_mpi_grp, weight)
480 end if
481
482 safe_deallocate_a(ddot)
483 safe_deallocate_a(zdot)
484
485 if (mpi_grp_is_root(mpi_world)) then
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
491 tdos(iorb) = m_zero
492 do ist = 1, st%nst
493 do ik = 1, st%nik
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)
496 end do
497 end do
498 end do
499
500 write(iunit(0), trim(format_str)) units_from_atomic(units_out%energy, energy), &
501 units_from_atomic(unit_one / units_out%energy, sum(tdos)), &
502 (units_from_atomic(unit_one / units_out%energy, tdos(iorb)), iorb=1,os%norbs)
503 end do
504 safe_deallocate_a(tdos)
505 call io_close(iunit(0))
506 end if
507
508 call orbitalset_end(os)
509 safe_deallocate_a(weight)
510
511 end do
512
513 end do
514
515 safe_deallocate_a(dpsi)
516 safe_deallocate_a(zpsi)
517 end if
518
519 safe_deallocate_a(iunit)
520 safe_deallocate_a(dos)
521
522 pop_sub(dos_write_dos)
523 end subroutine dos_write_dos
524
525 ! ---------------------------------------------------------
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
530 type(states_elec_t), intent(in) :: st
531 class(box_t), intent(in) :: box
532 type(ions_t), target, intent(in) :: ions
533 class(mesh_t), intent(in) :: mesh
534 type(hamiltonian_elec_t), intent(in) :: hm
535 type(namespace_t), intent(in) :: namespace
536
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
543
544 push_sub(dos_write_jdos)
545
546 ! shortcuts
547 ns = 1
548 if (st%d%nspin == 2) ns = 2
549
550 if (mpi_grp_is_root(mpi_world)) then
551 ! space for state-dependent DOS
552 safe_allocate(jdos(1:this%epoints, 0:ns-1))
553 safe_allocate(iunit(0:ns-1))
554 jdos = m_zero
555
556 ! compute band/spin-resolved density of states
557 do val = 1, st%nst
558 do cond = val, st%nst
559 do ik = 1, st%nik, ns
560 do is = 0, ns-1
561 if(st%occ(val, ik+is) < m_epsilon) cycle
562 if(st%occ(cond, ik+is) > m_epsilon) cycle
563 do ie = 1, this%epoints
564 energy = (ie - 1) * this%de
565 ! sum up Lorentzians
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)
568 end do
569 end do
570 end do
571 end do
572 end do
573
574 ! for spin-polarized calculations also output spin-resolved tDOS
575 if (st%d%nspin > 1) then
576 do is = 0, ns-1
577 write(filename, '(a,i1.1,a)') 'total-jdos-', is+1,'.dat'
578 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
579 ! write header
580 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS (spin-resolved)'
581
582 do ie = 1, this%epoints
583 energy = (ie - 1) * this%de
584 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
585 units_from_atomic(unit_one / units_out%energy, jdos(ie, is))
586 call messages_info(1, iunit(is))
587 end do
588
589 call io_close(iunit(is))
590 end do
591 end if
592
593
594 iunit(0) = io_open(trim(dir)//'/'//'total-jdos.dat', namespace, action='write')
595 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS'
596
597 ! compute total joint density of states
598 do ie = 1, this%epoints
599 energy = (ie - 1) * this%de
600 tjdos(1) = m_zero
601 do is = 0, ns-1
602 tjdos(1) = tjdos(1) + jdos(ie, is)
603 end do
604 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
605 units_from_atomic(unit_one / units_out%energy, tjdos(1))
606 call messages_info(1, iunit(0))
607 end do
608
609 call io_close(iunit(0))
610 end if
611
612 safe_deallocate_a(iunit)
613 safe_deallocate_a(jdos)
614
615 pop_sub(dos_write_jdos)
616 end subroutine dos_write_jdos
617
618
619 ! ---------------------------------------------------------
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
624 type(states_elec_t), intent(in) :: st
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
629 type(namespace_t), intent(in) :: namespace
630
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(:,:)
636 type(unit_t) :: fn_unit
637
638 push_sub(dos_write_ldos)
639
640 if (this%ldos_nenergies < 1) then
641 message(1) = "LDOSEnergies must be defined for Output=ldos"
642 call messages_fatal(1, namespace=namespace)
643 end if
644
645 ! shortcuts
646 ns = 1
647 if (st%d%nspin == 2) ns = 2
648
649 fn_unit = units_out%length**(-ions%space%dim) / units_out%energy
650
651 ! space for state-dependent DOS
652 safe_allocate(ldos(1:mesh%np, 1:this%ldos_nenergies, 1:ns))
653 ldos = m_zero
654
655 safe_allocate(abs_psi2(1:mesh%np))
656 if (states_are_real(st)) then
657 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
658 else
659 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
660 end if
661
662 ! compute band/spin-resolved density of states
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
666
667 ! TODO: This will be replaced by a poin-wise batch multiplication
668 if (states_are_real(st)) then
669 call states_elec_get_state(st, mesh, ist, ik, dpsi)
670 do ip = 1, mesh%np
671 abs_psi2(ip) = dpsi(ip, 1)**2
672 end do
673 if (st%d%dim > 1) then
674 abs_psi2(ip) = abs_psi2(ip) + dpsi(ip, 2)**2
675 end if
676 else
677 call states_elec_get_state(st, mesh, ist, ik, zpsi)
678 do ip = 1, mesh%np
679 abs_psi2(ip) = real(conjg(zpsi(ip, 1)) * zpsi(ip, 1), real64)
680 end do
681 if (st%d%dim > 1) then
682 abs_psi2(ip) = abs_psi2(ip) + real(conjg(zpsi(ip, 2)) * zpsi(ip, 2), real64)
683 end if
684 end if
685
686 ! sum up Lorentzians
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)
690
691 call lalg_axpy(mesh%np, weight, abs_psi2, ldos(:, ie, is))
692 end do
693 end do
694 end do !ist
695
696 safe_deallocate_a(dpsi)
697 safe_deallocate_a(zpsi)
698 safe_deallocate_a(abs_psi2)
699
700 if (st%parallel_in_states .or. st%d%kpt%parallel) then
701 call comm_allreduce(st%st_kpt_mpi_grp, ldos)
702 end if
703
704 do is = 1, ns
705 do ie = 1, this%ldos_nenergies
706 write(name, '(a,i3.3)') 'ldos_en-', ie
707 fname = get_filename_with_spin(name, st%d%nspin, is)
708
709 call dio_function_output(how, dir, fname, namespace, ions%space, mesh, &
710 ldos(:, ie, is), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
711 end do
712 end do
713
714 safe_deallocate_a(ldos)
715
716 pop_sub(dos_write_ldos)
717 end subroutine dos_write_ldos
718
719
720end module dos_oct_m
721
722!! Local Variables:
723!! mode: f90
724!! coding: utf-8
725!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
pure logical function, public accel_is_enabled()
Definition: accel.F90:400
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:285
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:298
subroutine, public dos_write_ldos(this, dir, st, box, ions, mesh, how, namespace)
Computes and output the local DOS (LDOS)
Definition: dos.F90:715
subroutine, public dos_init(this, namespace, st, kpoints)
Initializes the dot_t object.
Definition: dos.F90:180
subroutine, public dos_write_jdos(this, dir, st, box, ions, mesh, hm, namespace)
Computes and output the joint DOS (LDOS)
Definition: dos.F90:621
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:887
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:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
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:203
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:229
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:279
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
Definition: orbitalset.F90:521
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
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), 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.
int true(void)