21 use,
intrinsic :: iso_fortran_env
77 class(mesh_t),
intent(in) :: mesh
78 real(real64),
intent(in) :: centroids(:, :)
79 integer,
intent(out) :: ip_to_ic(:)
81 integer :: ip, ic, icen
82 integer :: n_centroids
83 real(real64) :: min_dist, dist
84 real(real64),
allocatable :: point(:)
89 real(real64),
parameter :: tol = 1.0e-13_real64
94 assert(
size(ip_to_ic) == mesh%np)
96 n_centroids =
size(centroids, 2)
97 safe_allocate(point(1:
size(centroids, 1)))
103 point = mesh%x(ip, :)
105 min_dist = sum((centroids(:, 1) - point(:))**2)
106 do ic = 2, n_centroids
107 dist = sum((centroids(:, ic) - point(:))**2)
108 if (dist < min_dist - tol)
then
117 safe_deallocate_a(point)
137 class(mesh_t),
intent(in) :: mesh
138 real(real64),
intent(in) :: weight(:)
139 integer,
intent(in) :: ip_to_ic(:)
140 real(real64),
contiguous,
intent(inout) :: centroids(:, :)
143 integer :: n_centroids, ic, ip
144 real(real64) :: one_over_denom
145 real(real64),
allocatable :: denominator(:)
151 assert(mesh%np ==
size(weight))
153 n_centroids =
size(centroids, 2)
154 safe_allocate(denominator(1:n_centroids))
156 do ic = 1, n_centroids
157 centroids(:, ic) = 0._real64
158 denominator(ic) = 0._real64
165 centroids(:, ic) = centroids(:, ic) + (mesh%x(ip, :) * weight(ip))
166 denominator(ic) = denominator(ic) + weight(ip)
171 call mesh%allreduce(centroids)
172 call mesh%allreduce(denominator)
178 do ic = 1, n_centroids
179 one_over_denom = 1._real64 / denominator(ic)
180 centroids(:, ic) = centroids(:, ic) * one_over_denom
184 safe_deallocate_a(denominator)
195 real(real64),
intent(in) :: points(:, :)
196 real(real64),
intent(in) :: updated_points(:, :)
197 real(real64),
intent(in) :: tol
198 logical,
intent(out) :: points_differ(:)
201 real(real64),
allocatable :: diff(:)
205 n_dim =
size(points, 1)
206 allocate(diff(n_dim))
209 do ip = 1,
size(points, 2)
210 diff(:) = abs(updated_points(:, ip) - points(:, ip))
211 points_differ(ip) = any(diff > tol)
226 real(real64),
intent(in) :: points(:, :)
227 real(real64),
intent(in) :: updated_points(:, :)
228 real(real64),
intent(in) :: tol
229 logical,
intent(in) :: points_differ(:)
231 integer,
allocatable :: indices(:)
232 integer :: i, j, n_unconverged, ndim
233 character(len=50) :: f_string
234 real(real64),
allocatable :: diff(:)
238 indices = pack([(i, i=1,
size(points_differ))], points_differ)
239 n_unconverged =
size(indices)
240 ndim =
size(points, 1)
243 write(f_string,
'(A, I1, A, I1, A, I1, A)')
'(', &
244 & ndim,
'(F16.10, X), ', &
245 & ndim,
'(F16.10, X), ', &
246 & ndim,
'(F16.10, X), F16.10)'
248 write(
message(1),
'(a)')
"# Current Point , Prior Point , |ri - r_{i-1}| , tol"
250 do j = 1, n_unconverged
252 diff(:) = abs(updated_points(:, i) - points(:, i))
253 write(
message(1), f_string) updated_points(:, i), points(:, i), diff, tol
256 write(
message(1), *)
"Summary:", n_unconverged,
"of out",
size(points, 2),
"are not converged"
293 subroutine weighted_kmeans(space, mesh, weight, centroids, n_iter, centroid_tol, discretize, inertia)
294 class(
space_t),
intent(in) :: space
295 class(
mesh_t),
intent(in) :: mesh
296 real(real64),
intent(in) :: weight(:)
297 real(real64),
contiguous,
intent(inout) :: centroids(:, :)
299 integer,
optional,
intent(in ) :: n_iter
300 real(real64),
optional,
intent(in ) :: centroid_tol
301 logical,
optional,
intent(in ) :: discretize
302 real(real64),
optional,
intent(out) :: inertia
304 logical :: discretize_centroids
305 integer :: n_iterations, n_centroid, i
307 integer,
allocatable :: ip_to_ic(:)
308 real(real64),
allocatable :: prior_centroids(:, :)
309 logical,
allocatable :: points_differ(:)
318 assert(n_iterations >= 1)
320 assert(
size(weight) == mesh%np)
322 assert(
size(centroids, 1) == space%dim)
324 assert(.not. space%is_periodic())
327 n_centroid =
size(centroids, 2)
328 safe_allocate_source(prior_centroids(space%dim,
size(centroids, 2)), centroids)
329 safe_allocate(ip_to_ic(1:mesh%np))
330 safe_allocate(points_differ(1:n_centroid))
332 write(
message(1),
'(a)')
'Debug: Performing weighted Kmeans clustering '
335 do i = 1, n_iterations
336 write(
message(1),
'(a, I3)')
'Debug: Iteration ', i
344 if (any(points_differ))
then
345 prior_centroids = centroids
347 write(
message(1),
'(a)')
'Debug: All centroid points converged'
355 if (discretize_centroids)
then
359 if (
present(inertia))
then
363 safe_deallocate_a(prior_centroids)
364 safe_deallocate_a(ip_to_ic)
365 safe_deallocate_a(points_differ)
377 class(
mesh_t),
intent(in ) :: mesh
378 real(real64),
contiguous,
intent(out) :: centroids(:, :)
379 integer(int64),
intent(inout),
optional :: seed_value
382 integer(int32) :: n_centroids
383 integer(int64),
allocatable :: centroid_idx(:)
384 integer(int64) :: ipg
389 n_centroids =
size(centroids, 2)
390 safe_allocate(centroid_idx(1:n_centroids))
396 do ic = 1, n_centroids
397 ipg = centroid_idx(ic)
401 safe_deallocate_a(centroid_idx)
418 class(
mesh_t),
intent(in) :: mesh
419 real(real64),
intent(in) :: centroids(:, :)
420 real(real64),
intent(in) :: weight(:)
421 integer,
intent(in) :: ip_to_ic(:)
422 real(real64),
intent(out) :: inertia
430 inertia = inertia + weight(ip) * sum((centroids(:, ic) - mesh%x(ip, :))**2)
434 call mesh%allreduce(inertia)
type(debug_t), save, public debug
subroutine report_differences_in_grids(points, updated_points, tol, points_differ)
Report differences returned from compute_grid_difference.
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 assign_points_to_centroids_finite_bc(mesh, centroids, ip_to_ic)
Assign each grid point to the closest centroid. A centroid and its set of nearest grid points defines...
subroutine, public update_centroids(mesh, weight, ip_to_ic, centroids)
Compute a new set of centroids.
subroutine compute_centroid_inertia(mesh, centroids, weight, ip_to_ic, inertia)
Compute the inertia of all centroids.
subroutine compute_grid_difference(points, updated_points, tol, points_differ)
Compute the difference in two grids as .
This module defines the meshes, which are used in Octopus.
subroutine, public mesh_discretize_values_to_mesh(mesh, values)
Assign a set of values to their nearest discrete points on the mesh.
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
This module contains some common usage patterns of MPI routines.
subroutine, public fisher_yates_shuffle(m, n, seed, values)
Return m random numbers from a range of with no replacement.
This module is intended to contain "only mathematical" functions and procedures.
Describes mesh distribution to nodes.