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 !%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%dim
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 select case (hm%kpoints%simplex%dim)
343 case (3, 10) ! Triangle with 3 corners (with optionally 7 neighbors for the improved simplex method)
344 call simplex_dos_2d(e_simplex(1:hm%kpoints%simplex%dim), energy, dos_simplex)
345 case (4, 20) ! Tetrahedron with 4 corners (with optionally 16 neighbors for the improved simplex method)
346 call simplex_dos_3d(e_simplex(1:hm%kpoints%simplex%dim), energy, dos_simplex)
347 case default
348 assert(.false.)
349 end select
350 dos(ie, ist, is) = dos(ie, ist, is) + dos_simplex / hm%kpoints%simplex%n_simplices
351 end do
352 end do
353 call profiling_out("DOS_WRITE_TETRAHEDRA")
354 case default
355 assert(.false.)
356 end select
357 do is = 0, ns-1
358 write(message(1), '(2es25.16E3)') units_from_atomic(units_out%energy, energy), &
359 units_from_atomic(unit_one / units_out%energy, dos(ie, ist, is))
360 call messages_info(1, iunit(is))
361 end do
362 end do
363
364 do is = 0, ns-1
365 call io_close(iunit(is))
366 end do
367 end do
368
369 safe_allocate(tdos(1))
370
371 ! for spin-polarized calculations also output spin-resolved tDOS
372 if (st%d%nspin == spin_polarized) then
373 do is = 0, ns-1
374 write(filename, '(a,i1.1,a)') 'total-dos-', is+1,'.dat'
375 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
376 ! write header
377 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS (spin-resolved)'
378
379 do ie = 1, this%epoints
380 energy = this%emin + (ie - 1) * this%de
381 tdos(1) = m_zero
382 do ist = 1, st%nst
383 tdos(1) = tdos(1) + dos(ie, ist, is)
384 end do
385 write(message(1), '(2es25.16E3)') units_from_atomic(units_out%energy, energy), &
386 units_from_atomic(unit_one / units_out%energy, tdos(1))
387 call messages_info(1, iunit(is))
388 end do
389
390 call io_close(iunit(is))
391 end do
392 end if
393
394
395 iunit(0) = io_open(trim(dir)//'/'//'total-dos.dat', namespace, action='write')
396 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total DOS'
397
398 ! compute total density of states
399 do ie = 1, this%epoints
400 energy = this%emin + (ie - 1) * this%de
401 tdos(1) = m_zero
402 do ist = 1, st%nst
403 do is = 0, ns-1
404 tdos(1) = tdos(1) + dos(ie, ist, is)
405 end do
406 end do
407 write(message(1), '(2es25.16E3)') units_from_atomic(units_out%energy, energy), &
408 units_from_atomic(unit_one / units_out%energy, tdos(1))
409 call messages_info(1, iunit(0))
410 end do
411
412 call io_close(iunit(0))
413
414 safe_deallocate_a(tdos)
415
416
417 ! write Fermi file
418 iunit(0) = io_open(trim(dir)//'/'//'total-dos-efermi.dat', namespace, action='write')
419 write(message(1), '(3a)') '# Fermi energy [', trim(units_abbrev(units_out%energy)), &
420 '] in a format compatible with total-dos.dat'
421
422 ! this is the maximum that tdos can reach
423 maxdos = st%smear%el_per_state * st%nst
424
425 write(message(2), '(2es25.16E3)') units_from_atomic(units_out%energy, st%smear%e_fermi), m_zero
426 write(message(3), '(es25.16E3,i7)') units_from_atomic(units_out%energy, st%smear%e_fermi), maxdos
427
428 call messages_info(3, iunit(0))
429 call io_close(iunit(0))
430
431 end if
432
433 if (this%computepdos) then
434
435 !These variables are defined in basis_set/orbitalbasis.F90
436 call parse_variable(namespace, 'AOThreshold', 0.01_real64, threshold)
437 call parse_variable(namespace, 'AONormalize', .true., normalize)
438
439 do ia = 1, ions%natoms
440 !We first count how many orbital set we have
441 work = orbitalset_utils_count(ions%atom(ia)%species)
442
443 !We loop over the orbital sets of the atom ia
444 do norb = 1, work
445 call orbitalset_init(os)
446 os%spec => ions%atom(ia)%species
447
448 !We count the orbitals
449 work2 = 0
450 do iorb = 1, os%spec%get_niwfs()
451 call os%spec%get_iwf_ilm(iorb, 1, ii, ll, mm)
452 call os%spec%get_iwf_n(iorb, 1, nn)
453 if (ii == norb) then
454 os%ll = ll
455 os%nn = nn
456 os%ii = ii
457 os%radius = atomic_orbital_get_radius(os%spec, mesh, iorb, 1, &
458 option__aotruncation__ao_full, threshold)
459 work2 = work2 + 1
460 end if
461 end do
462 os%norbs = work2
463 os%ndim = 1
464 os%use_submesh = .false.
465 os%allocated_on_mesh = .true.
466 os%spec => ions%atom(ia)%species
467
468 do work = 1, os%norbs
469 ! We obtain the orbital
470 if (states_are_real(st)) then
471 call dget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
472 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
473 os, work, os%radius, os%ndim, use_mesh=.not.os%use_submesh, &
474 normalize = normalize)
475 else
476 call zget_atomic_orbital(namespace, ions%space, ions%latt, ions%pos(:,ia), &
477 ions%atom(ia)%species, mesh, os%sphere, os%ii, os%ll, os%jj, &
478 os, work, os%radius, os%ndim, &
479 use_mesh = .not. hm%phase%is_allocated() .and. .not. os%use_submesh, &
480 normalize = normalize)
481 end if
482 end do
483
484 if (hm%phase%is_allocated()) then
485 ! In case of complex wavefunction, we allocate the array for the phase correction
486 safe_allocate(os%phase(1:os%sphere%np, st%d%kpt%start:st%d%kpt%end))
487 os%phase(:,:) = m_zero
488
489 if (.not. os%use_submesh) then
490 safe_allocate(os%eorb_mesh(1:mesh%np, 1:os%norbs, 1:os%ndim, st%d%kpt%start:st%d%kpt%end))
491 os%eorb_mesh(:,:,:,:) = m_zero
492 else
493 safe_allocate(os%eorb_submesh(1:os%sphere%np, 1:os%ndim, 1:os%norbs, st%d%kpt%start:st%d%kpt%end))
494 os%eorb_submesh(:,:,:,:) = m_zero
495 end if
496
497 if (accel_is_enabled() .and. st%d%dim == 1) then
498 os%ldorbs_eorb = max(pad_pow2(os%sphere%np), 1)
499 if(.not. os%use_submesh) os%ldorbs_eorb = max(pad_pow2(os%sphere%mesh%np), 1)
500
501 safe_allocate(os%buff_eorb(st%d%kpt%start:st%d%kpt%end))
502 do ik= st%d%kpt%start, st%d%kpt%end
503 call accel_create_buffer(os%buff_eorb(ik), accel_mem_read_only, type_cmplx, os%ldorbs_eorb*os%norbs)
504 end do
505 end if
506
507 call orbitalset_update_phase(os, box%dim, st%d%kpt, hm%kpoints, (st%d%ispin==spin_polarized), &
508 vec_pot = hm%hm_base%uniform_vector_potential, &
509 vec_pot_var = hm%hm_base%vector_potential)
510 else
511 if (states_are_real(st)) then
512 call dorbitalset_transfer_to_device(os, st%d%kpt, .not. os%use_submesh)
513 else
514 call zorbitalset_transfer_to_device(os, st%d%kpt, .not. os%use_submesh)
515 end if
516 end if
517
518 if (st%system_grp%is_root()) then
519 if (os%nn /= 0) then
520 write(filename,'(a, i4.4, a1, a, i1.1, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
521 os%nn, l_notation(os%ll), '.dat'
522 else
523 write(filename,'(a, i4.4, a1, a, a1,a)') 'pdos-at', ia, '-', trim(os%spec%get_label()), &
524 l_notation(os%ll), '.dat'
525 end if
526
527 iunit(0) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
528 ! write header
529 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), &
530 '], projected DOS (total and orbital resolved)'
531 end if
532
533 if (states_are_real(st)) then
534 safe_allocate(ddot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
535 else
536 safe_allocate(zdot(1:st%d%dim, 1:os%norbs, 1:st%block_size))
537 end if
538
539 safe_allocate(weight(1:os%norbs,1:st%nik,1:st%nst))
540 weight(1:os%norbs,1:st%nik,1:st%nst) = m_zero
541
542 do ik = st%d%kpt%start, st%d%kpt%end
543 do ib = st%group%block_start, st%group%block_end
544
545 if (hm%phase%is_allocated()) then
546 safe_allocate(epsib)
547 call st%group%psib(ib, ik)%copy_to(epsib)
548 call hm%phase%apply_to(mesh, mesh%np, .false., epsib, src = st%group%psib(ib, ik))
549 else
550 epsib => st%group%psib(ib, ik)
551 end if
552
553 if (states_are_real(st)) then
554 call dorbitalset_get_coeff_batch(os, st%d%dim, epsib, ddot(:, :, :))
555 do ist = 1, st%group%psib(ib, ik)%nst
556 ind = st%group%psib(ib, ik)%ist(ist)
557 do iorb = 1, os%norbs
558 do idim = 1, st%d%dim
559 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(ddot(idim, iorb, ist))**2
560 end do
561 end do
562 end do
563 else
564 call zorbitalset_get_coeff_batch(os, st%d%dim, epsib, zdot(:, :, :))
565 do ist = 1, st%group%psib(ib, ik)%nst
566 ind = st%group%psib(ib, ik)%ist(ist)
567 do iorb = 1, os%norbs
568 do idim = 1, st%d%dim
569 weight(iorb, ik, ind) = weight(iorb, ik, ind) + st%kweights(ik) * abs(zdot(idim, iorb, ist))**2
570 end do
571 end do
572 end do
573 end if
574
575 if (hm%phase%is_allocated()) then
576 call epsib%end(copy=.false.)
577 safe_deallocate_p(epsib)
578 end if
579
580 end do
581 end do
582
583 if (st%parallel_in_states .or. st%d%kpt%parallel) then
584 call comm_allreduce(st%st_kpt_mpi_grp, weight)
585 end if
586
587 safe_deallocate_a(ddot)
588 safe_deallocate_a(zdot)
589
590 if (st%system_grp%is_root()) then
591 write(format_str,'(a,i5,a)') '(', os%norbs+2, 'es25.16E3)'
592 safe_allocate(tdos(1:os%norbs))
593 do ie = 1, this%epoints
594 energy = this%emin + (ie - 1) * this%de
595 do iorb = 1, os%norbs
596 tdos(iorb) = m_zero
597 do ist = 1, st%nst
598 do ik = 1, st%nik
599 tdos(iorb) = tdos(iorb) + weight(iorb,ik,ist) / this%gamma * &
600 smear_delta_function(smear, (energy - st%eigenval(ist, ik))/this%gamma)
601 end do
602 end do
603 end do
604
605 write(iunit(0), trim(format_str)) units_from_atomic(units_out%energy, energy), &
606 units_from_atomic(unit_one / units_out%energy, sum(tdos)), &
607 (units_from_atomic(unit_one / units_out%energy, tdos(iorb)), iorb=1,os%norbs)
608 end do
609 safe_deallocate_a(tdos)
610 call io_close(iunit(0))
611 end if
612
613 call orbitalset_end(os)
614 safe_deallocate_a(weight)
615
616 end do
617
618 end do
619
620 end if
621
622 safe_deallocate_a(iunit)
623 safe_deallocate_a(dos)
624
625 pop_sub(dos_write_dos)
626 end subroutine dos_write_dos
627
628 ! ---------------------------------------------------------
630 subroutine dos_write_jdos(this, dir, st, ions, hm, namespace)
631 type(dos_t), intent(in) :: this
632 character(len=*), intent(in) :: dir
633 type(states_elec_t), intent(in) :: st
634 type(ions_t), target, intent(in) :: ions
635 type(hamiltonian_elec_t), intent(in) :: hm
636 type(namespace_t), intent(in) :: namespace
637
638 integer :: ie, ik, val, cond, is, ns
639 integer, allocatable :: iunit(:)
640 real(real64) :: energy
641 real(real64) :: tjdos(1)
642 real(real64), allocatable :: jdos(:,:)
643 character(len=64) :: filename
644
645 type(smear_t) :: smear
646
647 push_sub(dos_write_jdos)
648
649 ! shortcuts
650 ns = 1
651 if (st%d%nspin == 2) ns = 2
652
653 ! set up smearing parameters
654 smear%method = this%smear_func
655 smear%MP_n = 1
656
657 if (st%system_grp%is_root()) then
658 ! space for state-dependent DOS
659 safe_allocate(jdos(1:this%epoints, 0:ns-1))
660 safe_allocate(iunit(0:ns-1))
661 jdos = m_zero
662
663 ! compute band/spin-resolved density of states
664 do val = 1, st%nst
665 do cond = val, st%nst
666 do ik = 1, st%nik, ns
667 do is = 0, ns-1
668 if(st%occ(val, ik+is) < m_epsilon) cycle
669 if(st%occ(cond, ik+is) > m_epsilon) cycle
670 do ie = 1, this%epoints
671 energy = (ie - 1) * this%de
672 ! sum up Lorentzians
673 jdos(ie, is) = jdos(ie, is) + st%kweights(ik+is) / this%gamma * &
674 smear_delta_function(smear, (energy - (st%eigenval(cond, ik+is)-st%eigenval(val, ik+is)))/this%gamma)
675 end do
676 end do
677 end do
678 end do
679 end do
680
681 ! for spin-polarized calculations also output spin-resolved tDOS
682 if (st%d%nspin > 1) then
683 do is = 0, ns-1
684 write(filename, '(a,i1.1,a)') 'total-jdos-', is+1,'.dat'
685 iunit(is) = io_open(trim(dir)//'/'//trim(filename), namespace, action='write')
686 ! write header
687 write(iunit(is), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS (spin-resolved)'
688
689 do ie = 1, this%epoints
690 energy = (ie - 1) * this%de
691 write(message(1), '(2es25.16E3)') units_from_atomic(units_out%energy, energy), &
692 units_from_atomic(unit_one / units_out%energy, jdos(ie, is))
693 call messages_info(1, iunit(is))
694 end do
695
696 call io_close(iunit(is))
697 end do
698 end if
699
700
701 iunit(0) = io_open(trim(dir)//'/'//'total-jdos.dat', namespace, action='write')
702 write(iunit(0), '(3a)') '# energy [', trim(units_abbrev(units_out%energy)), '], total JDOS'
703
704 ! compute total joint density of states
705 do ie = 1, this%epoints
706 energy = (ie - 1) * this%de
707 tjdos(1) = m_zero
708 do is = 0, ns-1
709 tjdos(1) = tjdos(1) + jdos(ie, is)
710 end do
711 write(message(1), '(2es25.16E3)') units_from_atomic(units_out%energy, energy), &
712 units_from_atomic(unit_one / units_out%energy, tjdos(1))
713 call messages_info(1, iunit(0))
714 end do
715
716 call io_close(iunit(0))
717 end if
718
719 safe_deallocate_a(iunit)
720 safe_deallocate_a(jdos)
721
722 pop_sub(dos_write_jdos)
723 end subroutine dos_write_jdos
724
726 ! ---------------------------------------------------------
728 subroutine dos_write_ldos(this, dir, st, ions, mesh, how, namespace)
729 type(dos_t), intent(in) :: this
730 character(len=*), intent(in) :: dir
731 type(states_elec_t), intent(in) :: st
732 type(ions_t), target, intent(in) :: ions
733 class(mesh_t), intent(in) :: mesh
734 integer(int64), intent(in) :: how
735 type(namespace_t), intent(in) :: namespace
736
737 integer :: ie, ik, ist, is, ns, ip, ierr
738 character(len=MAX_PATH_LEN) :: fname, name
739 real(real64) :: weight
740 real(real64), allocatable :: ldos(:,:,:), dpsi(:,:), abs_psi2(:)
741 complex(real64), allocatable :: zpsi(:,:)
742 type(unit_t) :: fn_unit
743
744 type(smear_t) :: smear
745
746 push_sub(dos_write_ldos)
747
748 if (this%ldos_nenergies < 1) then
749 message(1) = "LDOSEnergies must be defined for Output=ldos"
750 call messages_fatal(1, namespace=namespace)
751 end if
752
753 ! shortcuts
754 ns = 1
755 if (st%d%nspin == 2) ns = 2
756
757 ! set up smearing parameters
758 smear%method = this%smear_func
759 smear%MP_n = 1
760
761 fn_unit = units_out%length**(-ions%space%dim) / units_out%energy
762
763 ! space for state-dependent DOS
764 safe_allocate(ldos(1:mesh%np, 1:this%ldos_nenergies, 1:ns))
765 ldos = m_zero
766
767 safe_allocate(abs_psi2(1:mesh%np))
768 if (states_are_real(st)) then
769 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
770 else
771 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
772 end if
773
774 ! compute band/spin-resolved density of states
775 do ik = st%d%kpt%start, st%d%kpt%end
776 is = st%d%get_spin_index(ik)
777 do ist = st%st_start, st%st_end
778
779 ! TODO: This will be replaced by a poin-wise batch multiplication
780 if (states_are_real(st)) then
781 call states_elec_get_state(st, mesh, ist, ik, dpsi)
782 do ip = 1, mesh%np
783 abs_psi2(ip) = dpsi(ip, 1)**2
784 end do
785 if (st%d%dim > 1) then
786 abs_psi2(ip) = abs_psi2(ip) + dpsi(ip, 2)**2
787 end if
788 else
789 call states_elec_get_state(st, mesh, ist, ik, zpsi)
790 do ip = 1, mesh%np
791 abs_psi2(ip) = real(conjg(zpsi(ip, 1)) * zpsi(ip, 1), real64)
792 end do
793 if (st%d%dim > 1) then
794 abs_psi2(ip) = abs_psi2(ip) + real(conjg(zpsi(ip, 2)) * zpsi(ip, 2), real64)
795 end if
796 end if
797
798 ! sum up Lorentzians
799 do ie = 1, this%ldos_nenergies
800 weight = st%kweights(ik) / this%gamma * &
801 smear_delta_function(smear, (this%ldos_energies(ie) - st%eigenval(ist, ik))/this%gamma)
802 call lalg_axpy(mesh%np, weight, abs_psi2, ldos(:, ie, is))
803 end do
804 end do
805 end do !ist
806
807 safe_deallocate_a(dpsi)
808 safe_deallocate_a(zpsi)
809 safe_deallocate_a(abs_psi2)
810
811 if (st%parallel_in_states .or. st%d%kpt%parallel) then
812 call comm_allreduce(st%st_kpt_mpi_grp, ldos)
813 end if
814
815 do is = 1, ns
816 do ie = 1, this%ldos_nenergies
817 write(name, '(a,i5.5)') 'ldos_en-', ie
818 fname = get_filename_with_spin(name, st%d%nspin, is)
819
820 call dio_function_output(how, dir, fname, namespace, ions%space, mesh, &
821 ldos(:, ie, is), fn_unit, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
822 end do
823 end do
824
825 safe_deallocate_a(ldos)
826
827 pop_sub(dos_write_ldos)
828 end subroutine dos_write_ldos
829
830
831end module dos_oct_m
832
833!! Local Variables:
834!! mode: f90
835!! coding: utf-8
836!! 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:726
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:824
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_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
Definition: simplex.F90:670
subroutine, public simplex_dos_2d(etriangle, eF, dos)
Get only the DOS contribution of a single triangle.
Definition: simplex.F90:472
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)