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
53 use types_oct_m
54 use unit_oct_m
57
58 implicit none
59
60 private
61 public :: &
62 dos_t, &
63 dos_init, &
67
68 type dos_t
69 private
70 real(real64) :: emin
71 real(real64) :: emax
72 integer :: epoints
73 real(real64) :: de
74
75 real(real64) :: gamma
76
77 logical :: computepdos
78
79 integer :: ldos_nenergies = -1
80 real(real64), allocatable :: ldos_energies(:)
81
82 integer(int64) :: method
83 contains
84 final :: dos_end
85 end type dos_t
86
87contains
88
90 subroutine dos_init(this, namespace, st, kpoints)
91 type(dos_t), intent(out) :: this
92 type(namespace_t), intent(in) :: namespace
93 type(states_elec_t), intent(in) :: st
94 type(kpoints_t), intent(in) :: kpoints
95
96 real(real64) :: evalmin, evalmax, eextend
97 integer :: npath, ie
98 type(block_t) :: blk
99
100 push_sub(dos_init)
101
102 !The range of the dos is only calculated for physical points,
103 !without the one from a k-point path
104 npath = kpoints%nkpt_in_path()
105 if (st%nik > npath) then
106 evalmin = minval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
107 evalmax = maxval(st%eigenval(1:st%nst, 1:(st%nik-npath)))
108 else !In case we only have a path, e.g., a bandstructure calculation
109 evalmin = minval(st%eigenval(1:st%nst, 1:st%nik))
110 evalmax = maxval(st%eigenval(1:st%nst, 1:st%nik))
111 end if
112 ! we extend the energy mesh by this amount
113 eextend = (evalmax - evalmin) / m_four
114
115 !%Variable DOSMethod
116 !%Type integer
117 !%Default smear
118 !%Section Output
119 !%Description
120 !% Selects the method that is used to calculate the DOS.
121 !% Currently not used for PDOS and LDOS.
122 !%Option smear 0
123 !% Smearing with a Lorentzian smearing function.
124 !%Option tetrahedra 1
125 !% Linear tetrahedron method as in P. E. Bloechl, et al., <i>Phys. Rev. B</i> <b>49</b>, 16223 (1994).
126 !% Requires a regular Monkhorst-Pack generated by Octopus.
127 !%End
128 call parse_variable(namespace, 'DOSMethod', option__dosmethod__smear, this%method)
129
130 !%Variable DOSEnergyMin
131 !%Type float
132 !%Section Output
133 !%Description
134 !% Lower bound for the energy mesh of the DOS.
135 !% The default is the lowest eigenvalue, minus a quarter of the total range of eigenvalues.
136 !% This is ignored for the joint density of states, and the minimal energy is always set to zero.
137 !%End
138 call parse_variable(namespace, 'DOSEnergyMin', evalmin - eextend, this%emin, units_inp%energy)
139
140 !%Variable DOSEnergyMax
141 !%Type float
142 !%Section Output
143 !%Description
144 !% Upper bound for the energy mesh of the DOS.
145 !% The default is the highest eigenvalue, plus a quarter of the total range of eigenvalues.
146 !%End
147 call parse_variable(namespace, 'DOSEnergyMax', evalmax + eextend, this%emax, units_inp%energy)
148
149 !%Variable DOSEnergyPoints
150 !%Type integer
151 !%Default 500
152 !%Section Output
153 !%Description
154 !% Determines how many energy points <tt>Octopus</tt> should use for
155 !% the DOS energy grid.
156 !%End
157 call parse_variable(namespace, 'DOSEnergyPoints', 500, this%epoints)
158
159 !%Variable DOSGamma
160 !%Type float
161 !%Default 0.008 Ha
162 !%Section Output
163 !%Description
164 !% Determines the width of the Lorentzian which is used for the DOS sum.
165 !%End
166 call parse_variable(namespace, 'DOSGamma', 0.008_real64, this%gamma)
168 !%Variable DOSComputePDOS
169 !%Type logical
170 !%Default false
171 !%Section Output
172 !%Description
173 !% Determines if projected dos are computed or not.
174 !% At the moment, the PDOS is computed from the bare pseudo-atomic orbitals, directly taken from
175 !% the pseudopotentials. The orbitals are not orthonormalized, in order to preserve their
176 !% atomic orbitals character. As a consequence, the sum of the different PDOS does not integrate
177 !% to the total DOS.
178 !%
179 !% The radii of the orbitals are controled by the threshold defined by <tt>AOThreshold</tt>,
180 !% and the fact that they are normalized or not by <tt>AONormalize</tt>.
181 !%End
182 call parse_variable(namespace, 'DOSComputePDOS', .false., this%computepdos)
183
184 ! spacing for energy mesh
185 this%de = (this%emax - this%emin) / (this%epoints - 1)
186
187 !%Variable LDOSEnergies
188 !%Type block
189 !%Section Output
190 !%Description
191 !% Specifies the energies at which the LDOS is computed.
192 !%End
193 if (parse_block(global_namespace, 'LDOSEnergies', blk) == 0) then
194 ! There is one high symmetry k-point per line
195 this%ldos_nenergies = parse_block_cols(blk, 0)
196
197 safe_allocate(this%ldos_energies(1:this%ldos_nenergies))
198 do ie = 1, this%ldos_nenergies
199 call parse_block_float(blk, 0, ie-1, this%ldos_energies(ie))
200 end do
201 call parse_block_end(blk)
202 else
203 this%ldos_nenergies = -1
204 end if
205
206 pop_sub(dos_init)
207 end subroutine dos_init
208
210 subroutine dos_end(this)
211 type(dos_t), intent(inout) :: this
212
213 push_sub(dos_end)
214
215 safe_deallocate_a(this%ldos_energies)
216 this%ldos_nenergies = -1
217
218 pop_sub(dos_end)
219 end subroutine
220
221 ! ---------------------------------------------------------
223 subroutine dos_write_dos(this, dir, st, box, ions, mesh, hm, namespace)
224 type(dos_t), intent(in) :: this
225 character(len=*), intent(in) :: dir
226 type(states_elec_t), target, intent(in) :: st
227 class(box_t), intent(in) :: box
228 type(ions_t), target, intent(in) :: ions
229 class(mesh_t), intent(in) :: mesh
230 type(hamiltonian_elec_t), intent(in) :: hm
231 type(namespace_t), intent(in) :: namespace
232
233 integer :: ie, ik, ist, is, ns, maxdos, ib, ind
234 integer, allocatable :: iunit(:)
235 real(real64) :: energy
236 real(real64), allocatable :: tdos(:)
237 real(real64), allocatable :: dos(:,:,:)
238 character(len=64) :: filename,format_str
239 logical :: normalize
240
241 integer :: ii, ll, mm, nn, work, norb, work2
242 integer :: ia, iorb, idim
243 real(real64) :: threshold
244 real(real64), allocatable :: ddot(:,:,:)
245 complex(real64), allocatable :: zdot(:,:,:)
246 real(real64), allocatable :: weight(:,:,:)
247 type(orbitalset_t) :: os
248 type(wfs_elec_t), pointer :: epsib
249
250 real(real64) :: e_simplex(4)
251 real(real64) :: dos_simplex
252
253 push_sub(dos_write_dos)
254
255 ! shortcuts
256 ns = 1
257 if (st%d%nspin == 2) ns = 2
258
259 if (mpi_world%is_root()) then
260 ! space for state-dependent DOS
261 safe_allocate(dos(1:this%epoints, 1:st%nst, 0:ns-1))
262 safe_allocate(iunit(0:ns-1))
263
264 ! compute band/spin-resolved density of states
265 do ist = 1, st%nst
266
267 do is = 0, ns-1
268 if (ns > 1) then
269 write(filename, '(a,i4.4,a,i1.1,a)') 'dos-', ist, '-', is+1,'.dat'
270 else
271 write(filename, '(a,i4.4,a)') 'dos-', ist, '.dat'
272 end if
273 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
274 ! write header
275 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], band-resolved DOS'
276 end do
277
278 do ie = 1, this%epoints
279 energy = this%emin + (ie - 1) * this%de
280 dos(ie, ist, :) = m_zero
281 select case (this%method)
282 case (option__dosmethod__smear)
283 ! sum up Lorentzians
284 do ik = 1, st%nik, ns
285 do is = 0, ns-1
286 dos(ie, ist, is) = dos(ie, ist, is) + st%kweights(ik+is) * m_one/m_pi * &
287 this%gamma / ((energy - st%eigenval(ist, ik+is))**2 + this%gamma**2)
288 end do
289 end do
290 case (option__dosmethod__tetrahedra)
291 ! sum up tetrahedra
292 assert(associated(hm%kpoints%simplex))
293 call profiling_in("DOS_WRITE_TETRAHEDRA")
294 do ii = 1, hm%kpoints%simplex%n_simplices
295 do is = 0, ns - 1
296 do ll = 1, size(hm%kpoints%simplex%simplices, 2)
297 ik = hm%kpoints%simplex%simplices(ii, ll)
298 ik = ns * (ik - 1) + 1
299 e_simplex(ll) = st%eigenval(ist, ik+is)
300 end do
301 call simplex_dos_3d(e_simplex, energy, dos_simplex)
302 dos(ie, ist, is) = dos(ie, ist, is) + dos_simplex / hm%kpoints%simplex%n_simplices
303 end do
304 end do
305 call profiling_out("DOS_WRITE_TETRAHEDRA")
306 case default
307 assert(.false.)
308 end select
309 do is = 0, ns-1
310 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
311 units_from_atomic(unit_one / units_out%energy, dos(ie, ist, is))
312 call messages_info(1, iunit(is))
313 end do
314 end do
315
316 do is = 0, ns-1
317 call io_close(iunit(is))
318 end do
319 end do
320
321 safe_allocate(tdos(1))
322
323 ! for spin-polarized calculations also output spin-resolved tDOS
324 if (st%d%nspin == spin_polarized) then
325 do is = 0, ns-1
326 write(filename, '(a,i1.1,a)') 'total-dos-', is+1,'.dat'
327 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
328 ! write header
329 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS (spin-resolved)'
330
331 do ie = 1, this%epoints
332 energy = this%emin + (ie - 1) * this%de
333 tdos(1) = m_zero
334 do ist = 1, st%nst
335 tdos(1) = tdos(1) + dos(ie, ist, is)
336 end do
337 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
338 units_from_atomic(unit_one / units_out%energy, tdos(1))
339 call messages_info(1, iunit(is))
340 end do
341
342 call io_close(iunit(is))
343 end do
344 end if
345
346
347 iunit(0) = io_open(trim(dir)//'/'//'total-dos.dat', namespace, action='write')
348 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS'
349
350 ! compute total density of states
351 do ie = 1, this%epoints
352 energy = this%emin + (ie - 1) * this%de
353 tdos(1) = m_zero
354 do ist = 1, st%nst
355 do is = 0, ns-1
356 tdos(1) = tdos(1) + dos(ie, ist, is)
357 end do
358 end do
359 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
360 units_from_atomic(unit_one / units_out%energy, tdos(1))
361 call messages_info(1, iunit(0))
362 end do
363
364 call io_close(iunit(0))
365
366 safe_deallocate_a(tdos)
367
368
369 ! write Fermi file
370 iunit(0) = io_open(trim(dir)//'/'//'total-dos-efermi.dat', namespace, action='write')
371 write(message(1), '(3a)') '# Fermi energy [', trim(units_abbrev(units_out%energy)), &
372 '] in a format compatible with total-dos.dat'
373
374 ! this is the maximum that tdos can reach
375 maxdos = st%smear%el_per_state * st%nst
376
377 write(message(2), '(2f12.6)') units_from_atomic(units_out%energy, st%smear%e_fermi), m_zero
378 write(message(3), '(f12.6,i6)') units_from_atomic(units_out%energy, st%smear%e_fermi), maxdos
379
380 call messages_info(3, iunit(0))
381 call io_close(iunit(0))
382
383 end if
384
385 if (this%computepdos) then
386
387 !These variables are defined in basis_set/orbitalbasis.F90
388 call parse_variable(namespace, 'AOThreshold', 0.01_real64, threshold)
389 call parse_variable(namespace, 'AONormalize', .true., normalize)
390
391 do ia = 1, ions%natoms
392 !We first count how many orbital set we have
393 work = orbitalset_utils_count(ions%atom(ia)%species)
394
395 !We loop over the orbital sets of the atom ia
396 do norb = 1, work
397 call orbitalset_init(os)
398 os%spec => ions%atom(ia)%species
399
400 !We count the orbitals
401 work2 = 0
402 do iorb = 1, os%spec%get_niwfs()
403 call os%spec%get_iwf_ilm(iorb, 1, ii, ll, mm)
404 call os%spec%get_iwf_n(iorb, 1, nn)
405 if (ii == norb) then
406 os%ll = ll
407 os%nn = nn
408 os%ii = ii
409 os%radius = atomic_orbital_get_radius(os%spec, mesh, iorb, 1, &
410 option__aotruncation__ao_full, threshold)
411 work2 = work2 + 1
412 end if
413 end do
414 os%norbs = work2
415 os%ndim = 1
416 os%use_submesh = .false.
417 os%allocated_on_mesh = .true.
418 os%spec => ions%atom(ia)%species
419
420 do work = 1, os%norbs
421 ! We obtain the orbital
422 if (states_are_real(st)) then
423 call dget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
424 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
425 os, work, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
426 normalize = normalize)
427 else
428 call zget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
429 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
430 os, work, os%radius, os%ndim, &
431 use_mesh = .not. hm%phase%is_allocated() .and. .not. os%use_submesh, &
432 normalize = normalize)
433 end if
434 end do
435
436 if (hm%phase%is_allocated()) then
437 ! In case of complex wavefunction, we allocate the array for the phase correction
438 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
439 os%phase(:,:) = m_zero
440
441 if (.not. os%use_submesh) then
442 safe_allocate(os%eorb_mesh(1:mesh%np, 1:os%norbs, 1:os%ndim, st%d%kpt%start:st%d%kpt%end))
443 os%eorb_mesh(:,:,:,:) = m_zero
444 else
445 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
446 os%eorb_submesh(:,:,:,:) = m_zero
447 end if
448
449 if (accel_is_enabled() .and. st%d%dim == 1) then
450 os%ldorbs_eorb = max(pad_pow2(os%sphere%np), 1)
451 if(.not. os%use_submesh) os%ldorbs_eorb = max(pad_pow2(os%sphere%mesh%np), 1)
452
453 safe_allocate(os%buff_eorb(st%d%kpt%start:st%d%kpt%end))
454 do ik= st%d%kpt%start, st%d%kpt%end
455 call accel_create_buffer(os%buff_eorb(ik), accel_mem_read_only, type_cmplx, os%ldorbs_eorb*os%norbs)
456 end do
457 end if
458
459 call orbitalset_update_phase(os, box%dim, st%d%kpt, hm%kpoints, (st%d%ispin==spin_polarized), &
460 vec_pot = hm%hm_base%uniform_vector_potential, &
461 vec_pot_var = hm%hm_base%vector_potential)
462 end if
463
464 if (mpi_world%is_root()) then
465 if (os%nn /= 0) then
466 write(filename,'(a, i3.3, a1, a, i1.1, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
467 os%nn, l_notation(os%ll), '.dat'
468 else
469 write(filename,'(a, i3.3, a1, a, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
470 l_notation(os%ll), '.dat'
471 end if
472
473 iunit(0) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
474 ! write header
475 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), &
476 '], projected DOS (total and orbital resolved)'
477 end if
478
479 if (states_are_real(st)) then
480 safe_allocate(ddot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
481 else
482 safe_allocate(zdot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
483 end if
484
485 safe_allocate(weight(1:os%norbs,1:st%nik,1:st%nst))
486 weight(1:os%norbs,1:st%nik,1:st%nst) = m_zero
487
488 do ik = st%d%kpt%start, st%d%kpt%end
489 do ib = st%group%block_start, st%group%block_end
490
491 if (hm%phase%is_allocated()) then
492 safe_allocate(epsib)
493 call st%group%psib(ib, ik)%copy_to(epsib)
494 call hm%phase%apply_to(mesh, mesh%np, .false., epsib, src = st%group%psib(ib, ik))
495 else
496 epsib => st%group%psib(ib, ik)
497 end if
498
499 if (states_are_real(st)) then
500 call dorbitalset_get_coeff_batch(os, st%d%dim, epsib, ddot(:, :, :))
501 do ist = 1, st%group%psib(ib, ik)%nst
502 ind = st%group%psib(ib, ik)%ist(ist)
503 do iorb = 1, os%norbs
504 do idim = 1, st%d%dim
505 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(ddot(idim, iorb, ist))**2
506 end do
507 end do
508 end do
509 else
510 call zorbitalset_get_coeff_batch(os, st%d%dim, epsib, zdot(:, :, :))
511 do ist = 1, st%group%psib(ib, ik)%nst
512 ind = st%group%psib(ib, ik)%ist(ist)
513 do iorb = 1, os%norbs
514 do idim = 1, st%d%dim
515 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(zdot(idim, iorb, ist))**2
516 end do
517 end do
518 end do
519 end if
520
521 if (hm%phase%is_allocated()) then
522 call epsib%end(copy=.false.)
523 safe_deallocate_p(epsib)
524 end if
525
526 end do
527 end do
528
529 if (st%parallel_in_states .or. st%d%kpt%parallel) then
530 call comm_allreduce(st%st_kpt_mpi_grp, weight)
531 end if
532
533 safe_deallocate_a(ddot)
534 safe_deallocate_a(zdot)
535
536 if (mpi_world%is_root()) then
537 write(format_str,'(a,i2,a)') '(', os%norbs+2, 'f12.6)'
538 safe_allocate(tdos(1:os%norbs))
539 do ie = 1, this%epoints
540 energy = this%emin + (ie - 1) * this%de
541 do iorb = 1, os%norbs
542 tdos(iorb) = m_zero
543 do ist = 1, st%nst
544 do ik = 1, st%nik
545 tdos(iorb) = tdos(iorb) + weight(iorb,ik,ist) * m_one/m_pi * &
546 this%gamma / ((energy - st%eigenval(ist, ik))**2 + this%gamma**2)
547 end do
548 end do
549 end do
550
551 write(iunit(0), trim(format_str)) units_from_atomic(units_out%energy, energy), &
552 units_from_atomic(unit_one / units_out%energy, sum(tdos)), &
553 (units_from_atomic(unit_one / units_out%energy, tdos(iorb)), iorb=1,os%norbs)
554 end do
555 safe_deallocate_a(tdos)
556 call io_close(iunit(0))
557 end if
558
559 call orbitalset_end(os)
560 safe_deallocate_a(weight)
561
562 end do
563
564 end do
565
566 end if
567
568 safe_deallocate_a(iunit)
569 safe_deallocate_a(dos)
570
571 pop_sub(dos_write_dos)
572 end subroutine dos_write_dos
573
574 ! ---------------------------------------------------------
576 subroutine dos_write_jdos(this, dir, st, ions, hm, namespace)
577 type(dos_t), intent(in) :: this
578 character(len=*), intent(in) :: dir
579 type(states_elec_t), intent(in) :: st
580 type(ions_t), target, intent(in) :: ions
581 type(hamiltonian_elec_t), intent(in) :: hm
582 type(namespace_t), intent(in) :: namespace
583
584 integer :: ie, ik, val, cond, is, ns
585 integer, allocatable :: iunit(:)
586 real(real64) :: energy
587 real(real64) :: tjdos(1)
588 real(real64), allocatable :: jdos(:,:)
589 character(len=64) :: filename
590
591 push_sub(dos_write_jdos)
592
593 ! shortcuts
594 ns = 1
595 if (st%d%nspin == 2) ns = 2
596
597 if (mpi_world%is_root()) then
598 ! space for state-dependent DOS
599 safe_allocate(jdos(1:this%epoints, 0:ns-1))
600 safe_allocate(iunit(0:ns-1))
601 jdos = m_zero
602
603 ! compute band/spin-resolved density of states
604 do val = 1, st%nst
605 do cond = val, st%nst
606 do ik = 1, st%nik, ns
607 do is = 0, ns-1
608 if(st%occ(val, ik+is) < m_epsilon) cycle
609 if(st%occ(cond, ik+is) > m_epsilon) cycle
610 do ie = 1, this%epoints
611 energy = (ie - 1) * this%de
612 ! sum up Lorentzians
613 jdos(ie, is) = jdos(ie, is) + st%kweights(ik+is) * m_one/m_pi * &
614 this%gamma / ((energy - (st%eigenval(cond, ik+is)-st%eigenval(val, ik+is)))**2 + this%gamma**2)
615 end do
616 end do
617 end do
618 end do
619 end do
620
621 ! for spin-polarized calculations also output spin-resolved tDOS
622 if (st%d%nspin > 1) then
623 do is = 0, ns-1
624 write(filename, '(a,i1.1,a)') 'total-jdos-', is+1,'.dat'
625 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
626 ! write header
627 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS (spin-resolved)'
628
629 do ie = 1, this%epoints
630 energy = (ie - 1) * this%de
631 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
632 units_from_atomic(unit_one / units_out%energy, jdos(ie, is))
633 call messages_info(1, iunit(is))
634 end do
635
636 call io_close(iunit(is))
637 end do
638 end if
639
640
641 iunit(0) = io_open(trim(dir)//'/'//'total-jdos.dat', namespace, action='write')
642 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS'
643
644 ! compute total joint density of states
645 do ie = 1, this%epoints
646 energy = (ie - 1) * this%de
647 tjdos(1) = m_zero
648 do is = 0, ns-1
649 tjdos(1) = tjdos(1) + jdos(ie, is)
650 end do
651 write(message(1), '(2f12.6)') units_from_atomic(units_out%energy, energy), &
652 units_from_atomic(unit_one / units_out%energy, tjdos(1))
653 call messages_info(1, iunit(0))
654 end do
655
656 call io_close(iunit(0))
657 end if
658
659 safe_deallocate_a(iunit)
660 safe_deallocate_a(jdos)
661
662 pop_sub(dos_write_jdos)
663 end subroutine dos_write_jdos
664
665
666 ! ---------------------------------------------------------
668 subroutine dos_write_ldos(this, dir, st, ions, mesh, how, namespace)
669 type(dos_t), intent(in) :: this
670 character(len=*), intent(in) :: dir
671 type(states_elec_t), intent(in) :: st
672 type(ions_t), target, intent(in) :: ions
673 class(mesh_t), intent(in) :: mesh
674 integer(int64), intent(in) :: how
675 type(namespace_t), intent(in) :: namespace
676
677 integer :: ie, ik, ist, is, ns, ip, ierr
678 character(len=MAX_PATH_LEN) :: fname, name
679 real(real64) :: weight
680 real(real64), allocatable :: ldos(:,:,:), dpsi(:,:), abs_psi2(:)
681 complex(real64), allocatable :: zpsi(:,:)
682 type(unit_t) :: fn_unit
683
684 push_sub(dos_write_ldos)
685
686 if (this%ldos_nenergies < 1) then
687 message(1) = "LDOSEnergies must be defined for Output=ldos"
688 call messages_fatal(1, namespace=namespace)
689 end if
690
691 ! shortcuts
692 ns = 1
693 if (st%d%nspin == 2) ns = 2
694
695 fn_unit = units_out%length**(-ions%space%dim) / units_out%energy
696
697 ! space for state-dependent DOS
698 safe_allocate(ldos(1:mesh%np, 1:this%ldos_nenergies, 1:ns))
699 ldos = m_zero
700
701 safe_allocate(abs_psi2(1:mesh%np))
702 if (states_are_real(st)) then
703 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
704 else
705 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
706 end if
707
708 ! compute band/spin-resolved density of states
709 do ik = st%d%kpt%start, st%d%kpt%end
710 is = st%d%get_spin_index(ik)
711 do ist = st%st_start, st%st_end
712
713 ! TODO: This will be replaced by a poin-wise batch multiplication
714 if (states_are_real(st)) then
715 call states_elec_get_state(st, mesh, ist, ik, dpsi)
716 do ip = 1, mesh%np
717 abs_psi2(ip) = dpsi(ip, 1)**2
718 end do
719 if (st%d%dim > 1) then
720 abs_psi2(ip) = abs_psi2(ip) + dpsi(ip, 2)**2
721 end if
722 else
723 call states_elec_get_state(st, mesh, ist, ik, zpsi)
724 do ip = 1, mesh%np
725 abs_psi2(ip) = real(conjg(zpsi(ip, 1)) * zpsi(ip, 1), real64)
726 end do
727 if (st%d%dim > 1) then
728 abs_psi2(ip) = abs_psi2(ip) + real(conjg(zpsi(ip, 2)) * zpsi(ip, 2), real64)
729 end if
730 end if
731
732 ! sum up Lorentzians
733 do ie = 1, this%ldos_nenergies
734 weight = st%kweights(ik) * m_one/m_pi * &
735 this%gamma / ((this%ldos_energies(ie) - st%eigenval(ist, ik))**2 + this%gamma**2)
736
737 call lalg_axpy(mesh%np, weight, abs_psi2, ldos(:, ie, is))
738 end do
739 end do
740 end do !ist
741
742 safe_deallocate_a(dpsi)
743 safe_deallocate_a(zpsi)
744 safe_deallocate_a(abs_psi2)
745
746 if (st%parallel_in_states .or. st%d%kpt%parallel) then
747 call comm_allreduce(st%st_kpt_mpi_grp, ldos)
748 end if
749
750 do is = 1, ns
751 do ie = 1, this%ldos_nenergies
752 write(name, '(a,i3.3)') 'ldos_en-', ie
753 fname = get_filename_with_spin(name, st%d%nspin, is)
754
755 call dio_function_output(how, dir, fname, namespace, ions%space, mesh, &
756 ldos(:, ie, is), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
757 end do
758 end do
759
760 safe_deallocate_a(ldos)
761
762 pop_sub(dos_write_ldos)
763 end subroutine dos_write_ldos
764
765
766end module dos_oct_m
767
768!! Local Variables:
769!! mode: f90
770!! coding: utf-8
771!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
pure logical function, public accel_is_enabled()
Definition: accel.F90:395
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:118
subroutine dos_end(this)
Finalizer for the dos_t object.
Definition: dos.F90:306
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:319
subroutine, public dos_write_jdos(this, dir, st, ions, hm, namespace)
Computes and output the joint DOS (LDOS)
Definition: dos.F90:672
subroutine, public dos_init(this, namespace, st, kpoints)
Initializes the dot_t object.
Definition: dos.F90:186
subroutine, public dos_write_ldos(this, dir, st, ions, mesh, how, namespace)
Computes and output the local DOS (LDOS)
Definition: dos.F90:764
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_one
Definition: global.F90:192
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:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
Definition: math.F90:892
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
type(namespace_t), public global_namespace
Definition: namespace.F90:134
subroutine, public orbitalset_init(this)
Definition: orbitalset.F90:208
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:234
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:284
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
Definition: orbitalset.F90:588
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
Definition: output_low.F90:117
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:233
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public simplex_dos_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
Definition: simplex.F90:373
pure logical function, public states_are_real(st)
type(type_t), public type_cmplx
Definition: types.F90:136
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
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:143
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)