21 use,
intrinsic :: iso_fortran_env
67 class(mesh_t),
intent(in) :: mesh
68 real(real64),
intent(in) :: centroids(:, :)
69 integer,
intent(out) :: ip_to_ic(:)
71 integer :: ip, ic, icen
72 integer :: n_centroids
73 real(real64) :: min_dist, dist
74 real(real64),
allocatable :: point(:)
79 assert(
size(ip_to_ic) == mesh%np)
81 n_centroids =
size(centroids, 2)
82 safe_allocate(point(1:
size(centroids, 1)))
90 min_dist = sum((centroids(:, 1) - point(:))**2)
91 do ic = 2, n_centroids
92 dist = sum((centroids(:, ic) - point(:))**2)
93 if (dist < min_dist)
then
102 safe_deallocate_a(point)
120 class(mesh_t),
intent(in) :: mesh
121 real(real64),
intent(in) :: weight(:)
122 integer,
intent(in) :: ip_to_ic(:)
123 real(real64),
contiguous,
intent(inout) :: centroids(:, :)
126 integer :: n_centroids, ic, ip
127 real(real64) :: one_over_denom
128 real(real64),
allocatable :: denominator(:)
134 assert(mesh%np ==
size(weight))
136 n_centroids =
size(centroids, 2)
137 safe_allocate(denominator(1:n_centroids))
139 do ic = 1, n_centroids
140 centroids(:, ic) = 0._real64
141 denominator(ic) = 0._real64
148 centroids(:, ic) = centroids(:, ic) + (mesh%x(ip, :) * weight(ip))
149 denominator(ic) = denominator(ic) + weight(ip)
154 call mesh%allreduce(centroids)
155 call mesh%allreduce(denominator)
161 do ic = 1, n_centroids
162 one_over_denom = 1._real64 / denominator(ic)
163 centroids(:, ic) = centroids(:, ic) * one_over_denom
167 safe_deallocate_a(denominator)
178 real(real64),
intent(in) :: points(:, :)
179 real(real64),
intent(in) :: updated_points(:, :)
180 real(real64),
intent(in) :: tol
181 logical,
intent(out) :: points_differ(:)
184 real(real64),
allocatable :: diff(:)
188 n_dim =
size(points, 1)
189 allocate(diff(n_dim))
192 do ip = 1,
size(points, 2)
193 diff(:) = abs(updated_points(:, ip) - points(:, ip))
194 points_differ(ip) = any(diff > tol)
209 real(real64),
intent(in) :: points(:, :)
210 real(real64),
intent(in) :: updated_points(:, :)
211 real(real64),
intent(in) :: tol
212 logical,
intent(in) :: points_differ(:)
214 integer,
allocatable :: indices(:)
215 integer :: i, j, n_unconverged, ndim
216 character(len=50) :: f_string
217 real(real64),
allocatable :: diff(:)
221 indices = pack([(i, i=1,
size(points_differ))], points_differ)
222 n_unconverged =
size(indices)
223 ndim =
size(points, 1)
226 write(f_string,
'(A, I1, A, I1, A, I1, A)')
'(', &
227 & ndim,
'(F16.10, X), ', &
228 & ndim,
'(F16.10, X), ', &
229 & ndim,
'(F16.10, X), F16.10)'
231 write(
message(1),
'(a)')
"# Current Point , Prior Point , |ri - r_{i-1}| , tol"
233 do j = 1, n_unconverged
235 diff(:) = abs(updated_points(:, i) - points(:, i))
236 write(
message(1), f_string) updated_points(:, i), points(:, i), diff, tol
239 write(
message(1), *)
"Summary:", n_unconverged,
"of out",
size(points, 2),
"are not converged"
272 subroutine weighted_kmeans(space, mesh, weight, centroids, n_iter, centroid_tol, discretize)
273 class(
space_t),
intent(in) :: space
274 class(
mesh_t),
intent(in) :: mesh
275 real(real64),
intent(in) :: weight(:)
276 real(real64),
contiguous,
intent(inout) :: centroids(:, :)
278 integer,
optional,
intent(in) :: n_iter
279 real(real64),
optional,
intent(in) :: centroid_tol
280 logical,
optional,
intent(in) :: discretize
282 logical :: discretize_centroids
283 integer :: n_iterations, n_centroid, i
285 integer,
allocatable :: ip_to_ic(:)
286 real(real64),
allocatable :: prior_centroids(:, :)
287 logical,
allocatable :: points_differ(:)
296 assert(n_iterations >= 1)
298 assert(
size(weight) == mesh%np)
300 assert(
size(centroids, 1) == space%dim)
302 assert(.not. space%is_periodic())
305 n_centroid =
size(centroids, 2)
306 safe_allocate_source(prior_centroids(space%dim,
size(centroids, 2)), centroids)
307 safe_allocate(ip_to_ic(1:mesh%np))
308 safe_allocate(points_differ(1:n_centroid))
311 write(
message(1),
'(a)')
'Performing weighted Kmeans clustering '
315 do i = 1, n_iterations
317 write(
message(1),
'(a, I3)')
'Iteration ', i
326 if (any(points_differ))
then
327 prior_centroids = centroids
330 write(
message(1),
'(a)')
'All centroid points converged'
339 if (discretize_centroids)
then
343 safe_deallocate_a(prior_centroids)
344 safe_deallocate_a(ip_to_ic)
345 safe_deallocate_a(points_differ)
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 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 weighted_kmeans(space, mesh, weight, centroids, n_iter, centroid_tol, discretize)
Weighted K-means clustering.
subroutine, public update_centroids(mesh, weight, ip_to_ic, centroids)
Compute a new set of 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.
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
This module contains some common usage patterns of MPI routines.
Describes mesh distribution to nodes.