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, int64
22 use batch_oct_m, only: batch_packed
23 use comm_oct_m, only: comm_allreduce
24 use debug_oct_m
26 use global_oct_m
27 use io_oct_m, only: io_open, io_close
30 use mpi_oct_m, only: mpi_world, mpi_double_precision
37
38 implicit none
39
40 private
41 public :: &
47
48 ! TODO(Alex) Issue #1195 Extend ISDF to spin-polarised systems
50 integer, private, parameter :: ik = 1
51
52contains
53
55 subroutine output_matrix(namespace, fname, matrix)
56 type(namespace_t), intent(in) :: namespace
57 character(len=*), intent(in) :: fname
58 real(real64), contiguous, intent(in) :: matrix(:, :)
59
60 integer :: i, j, unit
61
62 push_sub(output_matrix)
63
64 write(message(1),'(a)') "ISDF. Outputting: " // trim(adjustl(fname))
65 call messages_info(1, namespace=namespace, debug_only=.true.)
66
67 if (mpi_world%is_root()) then
68
69 open(newunit=unit, file=trim(adjustl(fname)))
70 do j = 1, size(matrix, 2)
71 do i = 1, size(matrix, 1)
72 write(unit, *) matrix(i, j)
73 enddo
74 enddo
75 close(unit)
76
77 endif
78
79 pop_sub(output_matrix)
80
81 end subroutine output_matrix
82
89 subroutine isdf_potential(namespace, poisson_solver, isdf_vectors, isdf_dist, V_r_nu)
90 type(namespace_t), intent(in ) :: namespace
91 type(poisson_t), intent(in ) :: poisson_solver
92 real(real64), intent(in ), contiguous :: isdf_vectors(:, :)
93 type(distributed_t), intent(in) :: isdf_dist
94
95 real(real64), intent(out), contiguous :: V_r_nu(:, :)
96
97 integer :: np, i_nu
98 integer, allocatable :: recv_counts(:), displs(:)
99 real(real64), allocatable :: V_r_nu_local(:, :)
100
101 push_sub_with_profile(isdf_potential)
102
103 ! Local buffer
104 np = size(v_r_nu, 1)
105 safe_allocate(v_r_nu_local(1:np, isdf_dist%start:isdf_dist%end))
106
107 do i_nu = isdf_dist%start, isdf_dist%end
108 call dpoisson_solve(poisson_solver, namespace, v_r_nu_local(:, i_nu), isdf_vectors(:, i_nu), all_nodes=.false.)
109 enddo
110
111 ! As an initial step, gather V_r_nu rather than keeping it distributed
112 safe_allocate(recv_counts(1:isdf_dist%mpi_grp%size))
113 safe_allocate(displs(1:isdf_dist%mpi_grp%size))
114
115 ! isdf_dist%num corresponds to number of columns per process
116 recv_counts(:) = np * isdf_dist%num(:)
117 call mpi_displacements(recv_counts, displs)
118
119 call isdf_dist%mpi_grp%allgatherv(v_r_nu_local, size(v_r_nu_local), mpi_double_precision, &
120 v_r_nu, recv_counts, displs, mpi_double_precision)
121
122 safe_deallocate_a(v_r_nu_local)
123
124 pop_sub_with_profile(isdf_potential)
125
126 end subroutine isdf_potential
127
128
134 integer function local_number_of_states(st, max_state)
135 type(states_elec_t), intent(in) :: st
136 integer, intent(in) :: max_state
137
138 integer :: minst, maxst
139
140 push_sub(local_number_of_states)
141
142 minst = states_elec_block_min(st, st%group%block_start)
143 maxst = min(states_elec_block_max(st, st%group%block_end), max_state)
144
145 ! Handles when max_state is part of a block with index < (block_start of current process)
146 ! resulting in no states on this process being used in the ISDF expansion
147 local_number_of_states = max(maxst - minst + 1, 0)
148
151 end function local_number_of_states
152
153
158 subroutine gather_psi_mu_over_states(st, psi_mu, psi_global)
159 type(states_elec_t), intent(in ) :: st
160 real(real64), contiguous, intent(in ) :: psi_mu(:, :)
161 real(real64), contiguous, intent(out) :: psi_global(:, :)
162
163 integer :: max_state, n_int_g, ic, ist, st_end, ist_local
164
166
167 if (st%group%psib(st%group%block_start, ik)%status() /= batch_packed) then
168 message(1) = "Developer Error: Trying to output psi_mu when not BATCH_PACKED"
169 call messages_fatal(1)
170 endif
171
172 ! Total number of interpolation points
173 n_int_g = size(psi_global, 2)
174 assert(size(psi_mu, 2) == size(psi_global, 2))
175
176 ! Total number of states used in ISDF
177 max_state = size(psi_global, 1)
178
179 ! All elements must be zeroed, such that allreduce does not
180 ! sum contributions from uninitialised elements
181 do ic = 1, n_int_g
182 do ist = 1, max_state
183 psi_global(ist, ic) = 0.0_real64
184 enddo
185 enddo
186
187 ! Ensure states does not iterate beyond the max state used in
188 ! ISDF expansion
189 st_end = min(st%st_end, max_state)
190
191 ! Sanity check on number of local states held by psi_mu
192 if (size(psi_mu, 1) > 0) then
193 assert(st_end - st%st_start + 1 == size(psi_mu, 1))
194 endif
195
196 ! Fill a section of psi_global using psi_mu
197 do ic = 1, n_int_g
198 do ist = st%st_start, st_end
199 ist_local = ist - st%st_start + 1
200 psi_global(ist, ic) = psi_mu(ist_local, ic)
201 enddo
202 enddo
203
204 call comm_allreduce(st%mpi_grp, psi_global)
205
207
208 end subroutine gather_psi_mu_over_states
209
210
213 subroutine output_psi_mu_for_all_states(namespace, st, max_state, psi_mu)
214 type(namespace_t), intent(in) :: namespace
215 type(states_elec_t), intent(in) :: st
216 integer, intent(in) :: max_state
217 real(real64), contiguous, intent(in) :: psi_mu(:, :)
218
219 real(real64), allocatable :: psi_global(:, :)
220 integer :: n_int_g, ic, ist, unit
221 character(len=2) :: np_char
222
224
225 if (st%group%psib(st%group%block_start, ik)%status() /= batch_packed) then
226 message(1) = "Trying to output psi_mu when not BATCH_PACKED"
227 call messages_fatal(1, namespace=namespace)
228 endif
230 write(message(1),'(a)') "ISDF: Writing psi_mu (all states/DD)"
231 call messages_info(1)
232
233 n_int_g = size(psi_mu, 2)
234 safe_allocate(psi_global(1:max_state, 1:n_int_g))
235
236 call gather_psi_mu_over_states(st, psi_mu, psi_global)
237
238 write(np_char, '(I2)') st%mpi_grp%size
239 unit = io_open("global_psi_mu_np"//trim(adjustl(np_char))//".txt", namespace, action='write')
240
241 do ic = 1, n_int_g
242 do ist = 1, max_state
243 write(unit, *) psi_global(ist, ic)
244 enddo
245 enddo
246
247 call io_close(unit)
248
249 safe_deallocate_a(psi_global)
250
252
254
255
256end module isdf_utils_oct_m
257
258!! Local Variables:
259!! mode: f90
260!! coding: utf-8
261!! End:
This module implements batches of mesh functions.
Definition: batch.F90:135
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:286
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
subroutine, public output_matrix(namespace, fname, matrix)
Helper routine to output a 2D matrix.
Definition: isdf_utils.F90:151
subroutine, public isdf_potential(namespace, poisson_solver, isdf_vectors, isdf_dist, V_r_nu)
Compute the effective potential in the ISDF vector basis.
Definition: isdf_utils.F90:185
integer function, public local_number_of_states(st, max_state)
Number of states contributing to the expansion, local to current process.
Definition: isdf_utils.F90:230
subroutine, public output_psi_mu_for_all_states(namespace, st, max_state, psi_mu)
Output the gathered psi_mu for all states, such that the matrix is the same irregardless of state par...
Definition: isdf_utils.F90:309
subroutine, public gather_psi_mu_over_states(st, psi_mu, psi_global)
Gather state-distributed psi from multiple processes.
Definition: isdf_utils.F90:254
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
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:875
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
The states_elec_t class contains all electronic wave functions.
int true(void)