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