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