21 use,
intrinsic :: iso_fortran_env
73 class(mesh_t),
intent(in) :: mesh
74 real(real64),
intent(in) :: centroids(:, :)
75 integer,
intent(out) :: ip_to_ic(:)
77 integer :: ip, ic, icen
78 integer :: n_centroids
79 real(real64) :: min_dist, dist
80 real(real64),
allocatable :: point(:)
85 real(real64),
parameter :: tol = 1.0e-13_real64
90 assert(
size(ip_to_ic) == mesh%np)
92 n_centroids =
size(centroids, 2)
93 safe_allocate(point(1:
size(centroids, 1)))
101 min_dist = sum((centroids(:, 1) - point(:))**2)
102 do ic = 2, n_centroids
103 dist = sum((centroids(:, ic) - point(:))**2)
104 if (dist < min_dist - tol)
then
113 safe_deallocate_a(point)
131 class(mesh_t),
intent(in) :: mesh
132 real(real64),
intent(in) :: weight(:)
133 integer,
intent(in) :: ip_to_ic(:)
134 real(real64),
contiguous,
intent(inout) :: centroids(:, :)
137 integer :: n_centroids, ic, ip
138 real(real64) :: one_over_denom
139 real(real64),
allocatable :: denominator(:)
145 assert(mesh%np ==
size(weight))
147 n_centroids =
size(centroids, 2)
148 safe_allocate(denominator(1:n_centroids))
150 do ic = 1, n_centroids
151 centroids(:, ic) = 0._real64
152 denominator(ic) = 0._real64
159 centroids(:, ic) = centroids(:, ic) + (mesh%x(ip, :) * weight(ip))
160 denominator(ic) = denominator(ic) + weight(ip)
165 call mesh%allreduce(centroids)
166 call mesh%allreduce(denominator)
172 do ic = 1, n_centroids
173 one_over_denom = 1._real64 / denominator(ic)
174 centroids(:, ic) = centroids(:, ic) * one_over_denom
178 safe_deallocate_a(denominator)
189 real(real64),
intent(in) :: points(:, :)
190 real(real64),
intent(in) :: updated_points(:, :)
191 real(real64),
intent(in) :: tol
192 logical,
intent(out) :: points_differ(:)
195 real(real64),
allocatable :: diff(:)
199 n_dim =
size(points, 1)
200 allocate(diff(n_dim))
203 do ip = 1,
size(points, 2)
204 diff(:) = abs(updated_points(:, ip) - points(:, ip))
205 points_differ(ip) = any(diff > tol)
220 real(real64),
intent(in) :: points(:, :)
221 real(real64),
intent(in) :: updated_points(:, :)
222 real(real64),
intent(in) :: tol
223 logical,
intent(in) :: points_differ(:)
225 integer,
allocatable :: indices(:)
226 integer :: i, j, n_unconverged, ndim
227 character(len=50) :: f_string
228 real(real64),
allocatable :: diff(:)
232 indices = pack([(i, i=1,
size(points_differ))], points_differ)
233 n_unconverged =
size(indices)
234 ndim =
size(points, 1)
237 write(f_string,
'(A, I1, A, I1, A, I1, A)')
'(', &
238 & ndim,
'(F16.10, X), ', &
239 & ndim,
'(F16.10, X), ', &
240 & ndim,
'(F16.10, X), F16.10)'
242 write(
message(1),
'(a)')
"# Current Point , Prior Point , |ri - r_{i-1}| , tol"
244 do j = 1, n_unconverged
246 diff(:) = abs(updated_points(:, i) - points(:, i))
247 write(
message(1), f_string) updated_points(:, i), points(:, i), diff, tol
250 write(
message(1), *)
"Summary:", n_unconverged,
"of out",
size(points, 2),
"are not converged"
283 subroutine weighted_kmeans(space, mesh, weight, centroids, n_iter, centroid_tol, discretize, inertia)
284 class(
space_t),
intent(in) :: space
285 class(
mesh_t),
intent(in) :: mesh
286 real(real64),
intent(in) :: weight(:)
287 real(real64),
contiguous,
intent(inout) :: centroids(:, :)
289 integer,
optional,
intent(in ) :: n_iter
290 real(real64),
optional,
intent(in ) :: centroid_tol
291 logical,
optional,
intent(in ) :: discretize
292 real(real64),
optional,
intent(out) :: inertia
294 logical :: discretize_centroids
295 integer :: n_iterations, n_centroid, i
297 integer,
allocatable :: ip_to_ic(:)
298 real(real64),
allocatable :: prior_centroids(:, :)
299 logical,
allocatable :: points_differ(:)
308 assert(n_iterations >= 1)
310 assert(
size(weight) == mesh%np)
312 assert(
size(centroids, 1) == space%dim)
314 assert(.not. space%is_periodic())
317 n_centroid =
size(centroids, 2)
318 safe_allocate_source(prior_centroids(space%dim,
size(centroids, 2)), centroids)
319 safe_allocate(ip_to_ic(1:mesh%np))
320 safe_allocate(points_differ(1:n_centroid))
322 write(
message(1),
'(a)')
'Debug: Performing weighted Kmeans clustering '
325 do i = 1, n_iterations
326 write(
message(1),
'(a, I3)')
'Debug: Iteration ', i
334 if (any(points_differ))
then
335 prior_centroids = centroids
337 write(
message(1),
'(a)')
'Debug: All centroid points converged'
345 if (discretize_centroids)
then
349 if (
present(inertia))
then
353 safe_deallocate_a(prior_centroids)
354 safe_deallocate_a(ip_to_ic)
355 safe_deallocate_a(points_differ)
367 class(
mesh_t),
intent(in ) :: mesh
368 real(real64),
contiguous,
intent(out) :: centroids(:, :)
369 integer(int64),
intent(inout),
optional :: seed_value
372 integer(int32) :: n_centroids
373 integer(int64),
allocatable :: centroid_idx(:)
374 integer(int64) :: ipg
379 n_centroids =
size(centroids, 2)
380 safe_allocate(centroid_idx(1:n_centroids))
386 do ic = 1, n_centroids
387 ipg = centroid_idx(ic)
391 safe_deallocate_a(centroid_idx)
408 class(
mesh_t),
intent(in) :: mesh
409 real(real64),
intent(in) :: centroids(:, :)
410 real(real64),
intent(in) :: weight(:)
411 integer,
intent(in) :: ip_to_ic(:)
412 real(real64),
intent(out) :: inertia
420 inertia = inertia + weight(ip) * sum((centroids(:, ic) - mesh%x(ip, :))**2)
424 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.