Octopus
isdf_utils.F90
Go to the documentation of this file.
1!! Copyright (C) 2025. Alexander Buccheri
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, U
17
18#include "global.h"
19
21 use, intrinsic :: iso_fortran_env, only: real64
22 use debug_oct_m
23 use global_oct_m
25 use mpi_oct_m, only: mpi_world
30 implicit none
31
32 private
33 public :: &
36
37contains
38
40 subroutine output_matrix(namespace, fname, matrix)
41 type(namespace_t), intent(in) :: namespace
42 character(len=*), intent(in) :: fname
43 real(real64), contiguous, intent(in) :: matrix(:, :)
44
45 integer :: i, j, unit
46
47 push_sub(output_matrix)
48
49 write(message(1),'(a)') "ISDF. Outputting: " // trim(adjustl(fname))
50 call messages_info(1, namespace=namespace, debug_only=.true.)
51
52 if (mpi_world%is_root()) then
53
54 open(newunit=unit, file=trim(adjustl(fname)))
55 do j = 1, size(matrix, 2)
56 do i = 1, size(matrix, 1)
57 write(unit, *) matrix(i, j)
58 enddo
59 enddo
60 close(unit)
61
62 endif
63
64 pop_sub(output_matrix)
65
66 end subroutine output_matrix
67
77 function highest_occupied_index(st, ik_index) result(i_max_occ)
78 type(states_elec_t), intent(in) :: st
79 integer, intent(in) :: ik_index
80 integer :: i_max_occ
81
82 integer :: ist
83
85
86 if (st%smear%method /= smear_semiconductor) then
87 i_max_occ = st%nst
88 endif
89
90 i_max_occ = 0
91 do ist = 1, st%nst
92 if (abs(st%occ(ist, ik_index)) < m_min_occ) exit
93 i_max_occ = ist
94 enddo
95 assert(i_max_occ > 0)
96
98
99 end function highest_occupied_index
100
101end module isdf_utils_oct_m
102
103!! Local Variables:
104!! mode: f90
105!! coding: utf-8
106!! End:
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:211
subroutine, public output_matrix(namespace, fname, matrix)
Helper routine to output a 2D matrix.
Definition: isdf_utils.F90:134
integer function, public highest_occupied_index(st, ik_index)
Return the index of highest occupied Kohn-Sham state for k-point ik.
Definition: isdf_utils.F90:171
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:270
integer, parameter, public smear_semiconductor
Definition: smear.F90:171
int true(void)