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