Octopus
isdf_options.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
22
24 use debug_oct_m
25 use global_oct_m
26 use io_oct_m
28 use mesh_oct_m
31 use parser_oct_m
33 use space_oct_m
35 implicit none
36 private
37
43 type, public :: isdf_options_t
44 logical :: check_n_interp
45 integer :: n_interp
46 integer :: n_ks_states = 0
47 ! This should be consistent with the subspace of ACE
48 integer(int64) :: seed
49 logical :: use_serial
50 type(centroids_t) :: centroids
51 real(real64), allocatable :: interpolation_coords(:, :)
52
53 contains
54 procedure :: init => isdf_options_init
55 procedure :: write_info => isdf_options_write_info
56 procedure :: get_interpolation_points => isdf_get_interpolation_points
57 procedure :: end => isdf_options_end
59 end type isdf_options_t
60
62 character(len=0), parameter, private :: NULL_FILE = ""
63
64contains
65
71 subroutine isdf_options_init(this, namespace, space, gr, n_ks_states, rng_seed)
72 class(isdf_options_t), intent(inout) :: this
73 type(namespace_t), intent(in ) :: namespace
74 class(space_t), intent(in ) :: space
75 class(mesh_t), intent(in ) :: gr
76 integer, intent(in ) :: n_ks_states
77 ! If ISDF is used with ACE, the number of states should be greater or equal to the subspace in which ACE
78 ! is exact. Else the density fitting will not capture the features of these states.
79 integer(int64), intent(inout), optional :: rng_seed
80
81 integer, parameter :: NULL_N_INTERP = 0
82 integer(int64), parameter :: default_rng_seed = 101_int64
83
84 push_sub(isdf_options_init)
85
86 this%seed = optional_default(rng_seed, default_rng_seed)
87 this%n_ks_states = n_ks_states
88
89 !%Variable ISDFNpoints
90 !%Type integer
91 !%Default 0
92 !%Section ISDF
93 !%Description
94 !% Specify the number of interpolation points, and therefore the number density fitting vectors.
95 !% The number of interpolation points should be considerably smaller than the number
96 !% of grid points. ISDFNpoints MUST be set. A sensible guess is 10 times the number of occupied states.
97 !%
98 !%Variable CheckISDFNpoints
99 !%Type logical
100 !%Default no
101 !%Section ISDF
102 !%Description
103 !% We can get a very good indication of the optimal number of interpolation points to use in ISDF
104 !% by choosing an excessive number,and computing the rank of the fitting coefficient Gram matrix, $CC^T$.
105 !% If the rank is equal to ISDFNpoints, this suggests one can choose more, or has the optimal number.
106 !% If the rank is less than ISDFNpoints, then it defines the upper bound beyond which the exact exchange energy
107 !% will not increase in precision. The same argument is not true of the two-electron integrals, which can be
108 !% converged further at the expense of an increase in the condition number of $CC^T$.
109 !%
110 !%Variable ISDFRunSerial
111 !%Type logical
112 !%Default no
113 !%Section ISDF
114 !%Description
115 !% Use the reference serial implementation of ISDF.
116 !%
117 !%End
118 call parse_variable(namespace, 'ISDFRunSerial', .false., this%use_serial)
119 call parse_variable(namespace, 'CheckISDFNpoints', .false., this%check_n_interp)
120
121 ! Set number of interpolation points
122 call parse_variable(namespace, 'ISDFNpoints', null_n_interp, this%n_interp)
123 if (this%n_interp == null_n_interp) then
124 message(1) = "ISDFNpoints must be specified"
125 call messages_fatal(1, namespace=namespace)
126 endif
127
128 safe_allocate(this%interpolation_coords(1:space%dim, 1:this%n_interp))
129 call sample_initial_centroids(gr, this%interpolation_coords, seed_value=this%seed)
130
131 pop_sub(isdf_options_init)
132
133 end subroutine isdf_options_init
134
135
137 subroutine isdf_options_end(this)
138 class(isdf_options_t) :: this
142 safe_deallocate_a(this%interpolation_coords)
143 call this%centroids%end()
147 end subroutine isdf_options_end
148
151 subroutine isdf_options_finalise(this)
152 type(isdf_options_t) :: this
154 push_sub(isdf_options_finalise)
155 call this%end()
156 pop_sub(isdf_options_finalise)
158 end subroutine isdf_options_finalise
159
160
162 subroutine isdf_options_write_info(this, namespace)
163 class(isdf_options_t), intent(in) :: this
164 type(namespace_t), intent(in) :: namespace
165
167
168 call messages_print_var_value("Interpolation Point Algorithm", "Weighted Kmeans Clustering", &
169 namespace=namespace)
170
171 call messages_print_var_value("Number of Interpolation Points (ISDFNpoints)", this%n_interp, &
172 namespace=namespace)
173
174 call messages_print_var_value("Number of states used in the ISDF vector construction", &
175 this%n_ks_states, namespace=namespace)
176
178
179 end subroutine isdf_options_write_info
180
181
183 subroutine isdf_get_interpolation_points(this, namespace, space, mesh, weight)
184 class(isdf_options_t), intent(inout) :: this
185 type(namespace_t), intent(in ) :: namespace
186 class(space_t), intent(in ) :: space
187 class(mesh_t), intent(in ) :: mesh
188 real(real64), intent(in ) :: weight(:)
189
190 logical :: check_duplicates
191
192 push_sub_with_profile(isdf_get_interpolation_points)
193
194 check_duplicates = debug%info .eqv. .true.
195
196 write(message(1),'(a, I3)') 'Info: Performing weighted KMeans'
197 call messages_info(1, namespace=namespace)
198 call weighted_kmeans(space, mesh, weight, this%interpolation_coords)
199 call this%centroids%end()
200 call this%centroids%init(mesh, this%interpolation_coords, check_duplicates = check_duplicates)
201
202 pop_sub_with_profile(isdf_get_interpolation_points)
203
204 end subroutine isdf_get_interpolation_points
205
206end module isdf_options_oct_m
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
type(debug_t), save, public debug
Definition: debug.F90:158
Definition: io.F90:116
subroutine isdf_get_interpolation_points(this, namespace, space, mesh, weight)
Populate the attribute centroids, with interpolation points.
subroutine isdf_options_finalise(this)
Automatically invoked finalise.
subroutine isdf_options_write_info(this, namespace)
Print ISDF settings via octopus messages.
subroutine isdf_options_init(this, namespace, space, gr, n_ks_states, rng_seed)
Initialise isdf_inp_options_t instance.
subroutine isdf_options_end(this)
Type-bound routine to manually free the attribute data.
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.
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: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
Encapsulate centroid points and their indexing across a domain-decomposed mesh.
Definition: centroids.F90:159
Describes mesh distribution to nodes.
Definition: mesh.F90:187
int true(void)