21 use,
intrinsic :: iso_fortran_env
44 logical :: check_n_interp
46 integer :: n_ks_states = 0
48 integer(int64) :: seed
50 type(centroids_t) :: centroids
51 real(real64),
allocatable :: interpolation_coords(:, :)
62 character(len=0),
parameter,
private :: NULL_FILE =
""
71 subroutine isdf_options_init(this, namespace, space, st, 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 type(states_elec_t),
intent(in ) :: st
76 class(mesh_t),
intent(in ) :: gr
77 integer,
intent(in ) :: n_ks_states
80 integer(int64),
intent(inout),
optional :: rng_seed
82 integer,
parameter :: NULL_N_INTERP = 0
83 integer(int64),
parameter :: default_rng_seed = 101_int64
88 this%n_ks_states = n_ks_states
119 call parse_variable(namespace,
'ISDFRunSerial', .false., this%use_serial)
120 call parse_variable(namespace,
'CheckISDFNpoints', .false., this%check_n_interp)
123 call parse_variable(namespace,
'ISDFNpoints', null_n_interp, this%n_interp)
124 if (this%n_interp == null_n_interp)
then
125 message(1) =
"ISDFNpoints must be specified"
129 safe_allocate(this%interpolation_coords(1:space%dim, 1:this%n_interp))
143 safe_deallocate_a(this%interpolation_coords)
144 call this%centroids%end()
176 this%n_ks_states, namespace=namespace)
187 class(
space_t),
intent(in ) :: space
188 class(
mesh_t),
intent(in ) :: mesh
189 real(real64),
intent(in ) :: weight(:)
191 logical :: check_duplicates
195 check_duplicates =
debug%info .eqv. .
true.
197 write(
message(1),
'(a, I3)')
'Info: Performing weighted KMeans'
200 call this%centroids%end()
201 call this%centroids%init(mesh, this%interpolation_coords, check_duplicates = check_duplicates)
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
type(debug_t), save, public debug
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_init(this, namespace, space, st, gr, n_ks_states, rng_seed)
Initialise isdf_inp_options_t instance.
subroutine isdf_options_write_info(this, namespace)
Print ISDF settings via octopus messages.
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.
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)
Encapsulate centroid points and their indexing across a domain-decomposed mesh.
Describes mesh distribution to nodes.