Octopus
test_isdf.F90
Go to the documentation of this file.
1!! Copyright (C) 2024 - 2025 A 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, USA.
17
18#include "global.h"
19
21module test_isdf_oct_m
22 use, intrinsic :: iso_fortran_env
23
26 use comm_oct_m
27 use debug_oct_m
30 use global_oct_m
34 use lcao_oct_m, only: lcao_run
35 use mesh_oct_m
37 use mpi_oct_m, only: mpi_world, mpi_double_precision
40 use scf_oct_m
42 use types_oct_m
45
46 ! Modules under test
49 use isdf_oct_m
51
52 implicit none
53 private
54 public :: test_isdf
55
56contains
57
60 subroutine test_isdf(namespace, serial)
61 type(namespace_t), intent(in) :: namespace
62 logical, optional, intent(in) :: serial
63
64 type(electrons_t), pointer :: sys
65 real(real64), parameter :: LMM_R_SINGLE_ATOM = 100.0_real64
66
67 ! Centroids
68 type(centroids_t) :: centroids
69 integer(int64), allocatable :: indices(:)
70 integer :: n_int
71
72 ! ISDF
73 type(isdf_options_t) :: isdf_options
74 logical :: output_cubes, run_serial_isdf
75 real(real64), allocatable :: phi(:, :)
76 real(real64), allocatable :: isdf_vectors(:, :), gram_matrix(:, :)
77
78 push_sub(test_isdf)
79
80 run_serial_isdf = optional_default(serial, .false.)
81
82 call messages_write('Info: Testing ISDF')
85 call messages_info(namespace=namespace)
86
87 ! Set up a system of electrons
88 call calc_mode_par%set_parallelization(p_strategy_states, default=.false.)
89 ! If generate_epot is false, the code will crash in v_ks_calc
90 sys => electrons_t(namespace, generate_epot=.true.)
91 call sys%init_parallelization(mpi_world)
92
93 ! Initialise a view of the wave functions
94 ! Note, this does not change the underlying data layout of the states
95 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_float, packed=.true.)
96
97 ! Guess at the states using lcao
98 call lcao_run(namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, &
99 sys%st, sys%ks, sys%hm, lmm_r = lmm_r_single_atom)
100
101 call isdf_options%init(namespace, sys%st)
102
103 call test_isdf_get_interpolation_points(namespace, isdf_options, sys, centroids)
104
105 ! Visualisation exact and approx product functions
106 output_cubes = .false.
107
108 if (run_serial_isdf) then
109 if (sys%st%parallel_in_states .or. sys%gr%parallel_in_domains) then
110 write(message(1), '(A)') 'Reference serial ISDF test requested, but state or domain parallel'
111 write(message(2), '(A)') ' specified in the input.'
112 call messages_fatal(2, namespace=namespace)
113 endif
115 ! Returns isdf_vectors, and phi for reuse
116 indices = centroids%global_mesh_indices()
117 n_int = size(indices)
118
120 sys%gr, sys%st, indices, phi, isdf_vectors)
121
122 safe_allocate(gram_matrix(1:n_int, 1:n_int))
123 call interpolation_vector_gram_matrix(sys%gr, isdf_vectors, gram_matrix)
124 call output_matrix(namespace, "isdf_gram.dat", gram_matrix)
125
126 ! This deallocates phi and isdf_vectors to save memory
127 call quantify_error_and_visualise(namespace, sys%st, sys%space, sys%gr, &
128 sys%ions, phi, indices, isdf_vectors, output_cubes)
129
130 deallocate(indices)
131
132 else
133 call messages_write('Info: Testing ISDF Parallel Implementation')
134 call messages_new_line()
135 call messages_info(namespace=namespace)
136
137 n_int = centroids%npoints_global()
138 if (mpi_world%is_root()) write(*, *) 'ISDF: Total number of centroids: ', n_int
139 write(*, *) "ISDF: Number of centroids local to domain process ", sys%gr%mpi_grp%rank, " is ", centroids%npoints()
140
141 call interpolative_separable_density_fitting_vectors(namespace, sys%gr, sys%st, centroids, isdf_vectors)
142
143 ! Gram matrix of the interpolation vectors offers a quantity reduced over states and domains
144 ! so good for testing parallelisation on
145 safe_allocate(gram_matrix(1:n_int, 1:n_int))
146 call interpolation_vector_gram_matrix(sys%gr, isdf_vectors, gram_matrix)
147 call output_matrix(namespace, "isdf_gram.dat", gram_matrix)
148 endif
149
150 ! Tear down
151 safe_deallocate_a(gram_matrix)
152 safe_deallocate_a(isdf_vectors)
154 safe_deallocate_p(sys)
155
156 pop_sub(test_isdf)
157
158 end subroutine test_isdf
159
160
162 subroutine test_isdf_get_interpolation_points(namespace, isdf_options, sys, centroids)
163 type(namespace_t), intent(in ) :: namespace
164 type(isdf_options_t), intent(in ) :: isdf_options
165 type(electrons_t), pointer, intent(in ) :: sys
166
167 type(centroids_t), intent(out) :: centroids
168
169 integer(int64) :: rng_seed
170 character(len=1), allocatable :: dummy_species(:)
171 real(real64), allocatable :: centroid_positions(:, :)
172
174
175 ! Perform kmeans - Random seed fixed for test reproducibility
176 rng_seed = 101_int64
177 safe_allocate(centroid_positions(1:sys%space%dim, 1:isdf_options%n_interp))
178
179 call sample_initial_centroids(sys%gr, centroid_positions, seed_value=rng_seed)
180
181 if (debug%info) then
182 allocate(dummy_species(isdf_options%n_interp), source='A')
183 call io_mkdir(static_dir, namespace)
184 call write_standard_xyz_file(namespace, static_dir // "/initial_centroids", centroid_positions, dummy_species)
185 endif
186
187 call weighted_kmeans(sys%space, sys%gr, sys%st%rho(1:sys%gr%np, 1), centroid_positions)
188
189 if (debug%info) then
190 dummy_species = 'B'
191 call write_standard_xyz_file(namespace, static_dir // "/final_centroids", centroid_positions, dummy_species)
192 endif
193
194 ! Convert positions to indices
195 call centroids%init(sys%gr, centroid_positions, check_duplicates=.true.)
196 safe_deallocate_a(centroid_positions)
197
198 ! Output the centroid indices so they can be asserted on
199 call centroids%output_all_indices("centroid_indices.dat")
200
202
204
205
206end module test_isdf_oct_m
207
208!! Local Variables:
209!! mode: f90
210!! coding: utf-8
211!! End:
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
type(debug_t), save, public debug
Definition: debug.F90:156
character(len= *), parameter, public static_dir
Definition: global.F90:252
subroutine, public write_standard_xyz_file(namespace, fname, pos, species, header)
Write a standard xyz file with atom labels and positions (in Angstrom).
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
Interoperable Separable Density Fitting (ISDF) molecular implementation.
Definition: isdf.F90:114
subroutine, public interpolative_separable_density_fitting_vectors(namespace, mesh, st, centroids, isdf_vectors)
Top-level routine for computing ISDF vectors.
Definition: isdf.F90:151
subroutine, public interpolation_vector_gram_matrix(mesh, isdf_vectors, gram_matrix)
Compute the Gram matrix for the ISDF interpolation vectors.
Definition: isdf.F90:374
Serial prototype for benchmarking and validating ISDF implementation.
subroutine, public serial_interpolative_separable_density_fitting_vectors(namespace, mesh, st, int_indices, phi, isdf_vectors)
Compute a set of ISDF interpolation vectors in serial, for code validation.
subroutine, public quantify_error_and_visualise(namespace, st, space, mesh, ions, phi, indices, isdf_vectors, output_cubes)
Wrapper for quantifying the error in the expansion of the product basis.
subroutine, public output_matrix(namespace, fname, matrix)
Helper routine to output a 2D matrix.
Definition: isdf_utils.F90:134
subroutine, public weighted_kmeans(space, mesh, weight, centroids, n_iter, centroid_tol, discretize, inertia)
Weighted K-means clustering.
subroutine, public sample_initial_centroids(mesh, centroids, seed_value)
Sample initial centroids from the full mesh.
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
Definition: lcao.F90:811
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_new_line()
Definition: messages.F90:1134
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
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
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
Integration tests for ISDF.
Definition: test_isdf.F90:114
subroutine, public test_isdf(namespace, serial)
Set up an electron system, compute some optimal centroid positions, and use these to build a set of I...
Definition: test_isdf.F90:154
subroutine test_isdf_get_interpolation_points(namespace, isdf_options, sys, centroids)
Initialise a system and centroids for molecular ISDF application testing.
Definition: test_isdf.F90:256
type(type_t), public type_float
Definition: types.F90:133
This module defines the unit system, used for input and output.
type(unit_t), public unit_angstrom
For XYZ files.
type(unit_t), public unit_one
some special units required for particular quantities
Encapsulate centroid points and their indexing across a domain-decomposed mesh.
Definition: centroids.F90:155
Class describing the electron system.
Definition: electrons.F90:218
int true(void)