22 use,
intrinsic :: iso_fortran_env
61 type(namespace_t),
intent(in) :: namespace
62 logical,
optional,
intent(in) :: serial
64 type(electrons_t),
pointer :: sys
65 real(real64),
parameter :: LMM_R_SINGLE_ATOM = 100.0_real64
68 type(centroids_t) :: centroids
69 integer(int64),
allocatable :: indices(:)
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(:, :)
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)
101 call isdf_options%init(namespace, sys%st)
106 output_cubes = .false.
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.'
116 indices = centroids%global_mesh_indices()
117 n_int =
size(indices)
120 sys%gr, sys%st, indices, phi, isdf_vectors)
122 safe_allocate(gram_matrix(1:n_int, 1:n_int))
128 sys%ions, phi, indices, isdf_vectors, output_cubes)
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()
145 safe_allocate(gram_matrix(1:n_int, 1:n_int))
151 safe_deallocate_a(gram_matrix)
152 safe_deallocate_a(isdf_vectors)
154 safe_deallocate_p(sys)
169 integer(int64) :: rng_seed
170 character(len=1),
allocatable :: dummy_species(:)
171 real(real64),
allocatable :: centroid_positions(:, :)
177 safe_allocate(centroid_positions(1:sys%space%dim, 1:isdf_options%n_interp))
182 allocate(dummy_species(isdf_options%n_interp), source=
'A')
187 call weighted_kmeans(sys%space, sys%gr, sys%st%rho(1:sys%gr%np, 1), centroid_positions)
195 call centroids%init(sys%gr, centroid_positions, check_duplicates=.
true.)
196 safe_deallocate_a(centroid_positions)
199 call centroids%output_all_indices(
"centroid_indices.dat")
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
character(len= *), parameter, public static_dir
subroutine, public write_standard_xyz_file(namespace, fname, pos, species, header)
Write a standard xyz file with atom labels and positions (in Angstrom).
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Interoperable Separable Density Fitting (ISDF) molecular implementation.
subroutine, public interpolative_separable_density_fitting_vectors(namespace, mesh, st, centroids, isdf_vectors)
Top-level routine for computing ISDF vectors.
subroutine, public interpolation_vector_gram_matrix(mesh, isdf_vectors, gram_matrix)
Compute the Gram matrix for the ISDF interpolation vectors.
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.
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)
This module defines the meshes, which are used in Octopus.
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
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.
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...
subroutine test_isdf_get_interpolation_points(namespace, isdf_options, sys, centroids)
Initialise a system and centroids for molecular ISDF application testing.
type(type_t), public type_float
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.
Class describing the electron system.