ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_clusKmeans::setKmeans Interface Reference

Compute and return an iteratively-refined set of cluster centers given the input sample using the k-means approach.
More...

Detailed Description

Compute and return an iteratively-refined set of cluster centers given the input sample using the k-means approach.

See the documentation of pm_clusKmeans for more information on the Kmeans clustering algorithm.
The metric used within this generic interface is the Euclidean distance.

Parameters
[in,out]rng: The input/output scalar that can be an object of,
  1. type rngf_type, implying the use of intrinsic Fortran uniform RNG.
  2. type xoshiro256ssw_type, implying the use of xoshiro256** uniform RNG.
(optional. If this argument is present, then all intent(inout) arguments below have intent(out) argument and will be initialized using the k-means++ algorithm.
If this argument is missing, then the user must initialize all of the following arguments with intent(inout) via k-means++ or any other method before passing them to this generic interface.)
[in,out]membership: The input/output vector of shape (1:nsam) of type integer of default kind IK, containing the membership of each input sample in sample from its nearest cluster center, such that cluster(membership(i)) is the nearest cluster center to the ith sample sample(:, i) at a squared-distance of disq(i).
If the argument rng is missing, then membership has intent(inout) and must be properly initialized before calling this routine.
If the argument rng is present, then membership has intent(out) and will contain the final cluster memberships on output.
[in,out]disq: The input/output vector of shape (1:nsam) of the same type and kind as the input argument sample, containing the Euclidean squared distance of each input sample in sample from its nearest cluster center.
If the argument rng is missing, then disq has intent(inout) and must be properly initialized before calling this routine.
If the argument rng is present, then disq has intent(out) and will contain the final squared-distances from cluster centers on output.
[in]sample: The input vector, or matrix of,
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the sample of nsam points in a ndim-dimensional space whose corresponding cluster centers must be computed.
  1. If sample is a vector of shape (1 : nsam) and center is a vector of shape (1 : ncls), then the input sample must be a collection of nsam points (in univariate space).
  2. If sample is a matrix of shape (1 : ndim, 1 : nsam) and center is a matrix of shape (1 : ndim, 1 : ncls), then the input sample must be a collection of nsam points (in ndim-dimensional space).
[in,out]center: The input/output vector of shape (1:ncls) or matrix of shape (1 : ndim, 1 : ncls) of the same type and kind as the input argument sample, containing the set of ncls unique cluster centers (centroids) computed based on the input sample memberships and minimum distances.
If the argument rng is missing, then center has intent(inout) and must be properly initialized before calling this routine.
If the argument rng is present, then center has intent(out) and will contain the final cluster centers on output.
[in,out]size: The input/output vector of shape (1:ncls) type integer of default kind IK, containing the sizes (number of members) of the clusters with the corresponding centers output in the argument center.
If the argument rng is missing, then size has intent(inout) and must be properly initialized before calling this routine.
If the argument rng is present, then size has intent(out) and will contain the final cluster sizes (member counts) on output.
[in,out]potential: The input/output vector of shape (1:ncls) of the same type and kind as the input argument sample, the ith element of which contains the sum of squared distances of all members of the ith cluster from the cluster center as output in the ith element of center.
If the argument rng is missing, then potential has intent(inout) and must be properly initialized before calling this routine (although its values are not explicitly referenced).
If the argument rng is present, then potential has intent(out) and will contain the final cluster potentials (sums of squared distances from cluster centers) on output.
[out]failed: The output scalar of type logical of default kind LK that is .true. if and only if the algorithm fails to converge within the user-specified or default criteria for convergence.
Failure occurs only if any(size < minsize) .or. maxniter < niter.
[out]niter: The output scalar of type integer of default kind IK, containing the number of refinement iterations performed within the algorithm to achieve convergence.
An output niter value larger than the input maxniter implies lack of convergence before return.
(optional. If missing, the number of refinement iterations will not be output.)
[in]maxniter: The input non-negative scalar of type integer of default kind IK, containing the maximum number of refinement iterations allowed within the algorithm to achieve convergence.
If convergence does not occur within the maximum specified value, the output arguments can be passed again as is to the generic interface (without the optional rng argument) to continue the refinement iterations until convergence.
A reasonable choice can be 300 or comparable values.
(optional, default = 300.)
[in]minsize: The input non-negative scalar of type integer of default kind IK, containing the minimum allowed size of each cluster.
If any cluster has any number of members below the specified minsize, the algorithm will return without achieving convergence.
The situation can be detected of any element of the output size is smaller than the specified minsize.
A reasonable choice can be ndim = ubound(sample, 1) or comparable values although any non-negative value including zero is possible.
(optional, default = 1.)
[in]nfail: The input non-negative scalar of type integer of default kind IK, containing the number of times the k-means algorithm is allowed to fail before returning without convergence.
(optional. It can be present only if the argument rng is also present, allowing random initializations in case of failures.)


Possible calling interfaces

call setKmeans(membership(1 : nsam), disq(1 : nsam), sample(1 : nsam) , center(1 : ncls) , size(1 : ncls), potential(1 : ncls), failed, niter = niter, maxniter = maxniter, minsize = minsize)
call setKmeans(membership(1 : nsam), disq(1 : nsam), sample(1 : ndim, 1 : nsam), center(1 : ndim, 1 : ncls), size(1 : ncls), potential(1 : ncls), failed, niter = niter, maxniter = maxniter, minsize = minsize)
call setKmeans(rng, membership(1 : nsam), disq(1 : nsam), csdisq(0 : nsam), sample(1 : nsam) , center(1 : ncls) , size(1 : ncls), potential(1 : ncls), failed, niter = niter, maxniter = maxniter, minsize = minsize, nfail = nfail)
call setKmeans(rng, membership(1 : nsam), disq(1 : nsam), csdisq(0 : nsam), sample(1 : ndim, 1 : nsam), center(1 : ndim, 1 : ncls), size(1 : ncls), potential(1 : ncls), failed, niter = niter, maxniter = maxniter, minsize = minsize, nfail = nfail)
Compute and return an iteratively-refined set of cluster centers given the input sample using the k-m...
This module contains procedures and routines for the computing the Kmeans clustering of a given set o...
Warning
If the argument rng is present, then all arguments associated with setKmeansPP equally apply to this generic interface.
The condition ubound(center, rank(center)) > 0 must hold for the corresponding input arguments.
The condition ubound(sample, rank(sample)) == ubound(disq, 1) must hold for the corresponding input arguments.
The condition ubound(center, rank(center)) == ubound(size, 1) must hold for the corresponding input arguments.
The condition ubound(center, rank(center)) == ubound(potential, 1) must hold for the corresponding input arguments.
The condition ubound(sample, rank(sample)) == ubound(membership, 1) must hold for the corresponding input arguments.
The condition ubound(center, rank(center)) <= ubound(sample, rank(sample)) must hold for the corresponding input arguments (the number of clusters must be less than or equal to the sample size).
The condition ubound(sample, 1) == ubound(center, 1) must hold for the corresponding input arguments.
The condition 0 <= maxniter must hold for the corresponding input arguments.
The condition 0 <= minsize must hold for the corresponding input arguments.
The condition 0 <= nfail must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds. The procedures of this generic interface are always impure when the input rng argument is set to an object of type rngf_type.
Remarks
The functionality of this generic interface is highly similar to setCenter, with the major difference being that setKmeans simultaneously computes the new cluster centers and sample memberships, whereas setCenter computes the new cluster centers based on given sample membership.
See also
setKmeans
setCenter
setMember
setKmeansPP
k-means clustering


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: RKG => RKS ! all other real kinds are also supported.
5 use pm_io, only: display_type
6 use pm_distUnif, only: getUnifRand
10
11 implicit none
12
13 logical(LK) :: failed
14 integer(IK) :: ndim, nsam, ncls, itry, niter
15 real(RKG) , allocatable :: sample(:,:), center(:,:), disq(:), csdisq(:), potential(:)
16 integer(IK) , allocatable :: membership(:), size(:)
17 type(display_type) :: disp
18
19 disp = display_type(file = "main.out.F90")
20
21 call disp%skip
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%show("! Compute cluster centers based on an input sample and cluster memberships and member-center distances.")
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%skip
26
27 do itry = 1, 10
28 call disp%skip()
29 call disp%show("ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);")
30 ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
31 call disp%show("[ndim, nsam, ncls]")
32 call disp%show( [ndim, nsam, ncls] )
33 call disp%show("sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.")
34 sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
35 call disp%show("sample")
36 call disp%show( sample )
37 call disp%show("call setResized(disq, nsam)")
38 call setResized(disq, nsam)
39 call disp%show("call setResized(csdisq, nsam + 1_IK)")
40 call setResized(csdisq, nsam + 1_IK)
41 call disp%show("call setResized(membership, nsam)")
42 call setResized(membership, nsam)
43 call disp%show("call setResized(center, [ndim, ncls])")
44 call setResized(center, [ndim, ncls])
45 call disp%show("call setResized(potential, ncls)")
46 call setResized(potential, ncls)
47 call disp%show("call setResized(size, ncls)")
48 call setResized(size, ncls)
49 call disp%skip()
50
51 call disp%show("call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.")
52 call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
53 call disp%show("failed")
54 call disp%show( failed )
55 call disp%show("niter")
56 call disp%show( niter )
57 call disp%show("disq")
58 call disp%show( disq )
59 call disp%show("csdisq")
60 call disp%show( csdisq )
61 call disp%show("membership")
62 call disp%show( membership )
63 call disp%show("potential")
64 call disp%show( potential )
65 call disp%show("center")
66 call disp%show( center )
67 call disp%show("size")
68 call disp%show( size )
69 call disp%skip()
70 end do
71
72 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
73 ! Output an example for visualization.
74 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
75
76 block
77 integer(IK) :: funit, i
78 ndim = 2
79 ncls = 5
80 nsam = 5000
81 center = getUnifRand(0., 1., ndim, ncls)
82 sample = getUnifRand(0., 1., ndim, nsam)
83 call setResized(disq, nsam)
84 call setResized(csdisq, nsam + 1_IK)
85 call setResized(membership, nsam)
86 call setResized(potential, ncls)
87 call setResized(size, ncls)
88 ! output the sample and the initial centers.
89 call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential)
90 open(newunit = funit, file = "setKmeansPP.center.txt")
91 do i = 1, ncls
92 write(funit, "(*(g0,:,','))") i, center(:,i)
93 end do
94 close(funit)
95 open(newunit = funit, file = "setKmeansPP.sample.txt")
96 do i = 1, nsam
97 write(funit, "(*(g0,:,','))") membership(i), sample(:,i)
98 end do
99 close(funit)
100 call setKmeans(membership, disq, sample, center, size, potential, failed)
101 open(newunit = funit, file = "setKmeans.center.txt")
102 do i = 1, ncls
103 write(funit, "(*(g0,:,','))") i, center(:,i)
104 end do
105 close(funit)
106 open(newunit = funit, file = "setKmeans.sample.txt")
107 do i = 1, nsam
108 write(funit, "(*(g0,:,','))") membership(i), sample(:,i)
109 end do
110 close(funit)
111 end block
112
113end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Compute and return an asymptotically optimal set of cluster centers for the input sample,...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains classes and procedures for computing various statistical quantities related to t...
type(rngf_type) rngf
The scalar constant object of type rngf_type whose presence signified the use of the Fortran intrinsi...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Compute cluster centers based on an input sample and cluster memberships and member-center distances.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
8[ndim, nsam, ncls]
9+3, +1, +1
10sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
11sample
12+3.83674002
13+2.53701353
14+3.95865202
15call setResized(disq, nsam)
16call setResized(csdisq, nsam + 1_IK)
17call setResized(membership, nsam)
18call setResized(center, [ndim, ncls])
19call setResized(potential, ncls)
20call setResized(size, ncls)
21
22call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
23failed
24F
25niter
26+1
27disq
28+0.00000000
29csdisq
30+0.00000000, +0.00000000
31membership
32+1
33potential
34+0.00000000
35center
36+3.83674002
37+2.53701353
38+3.95865202
39size
40+1
41
42
43ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
44[ndim, nsam, ncls]
45+1, +8, +5
46sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
47sample
48+3.45510840, +2.03410673, +2.91745758, +0.876081884, +4.99067736, +0.193027854, +1.06150985, +4.52673864
49call setResized(disq, nsam)
50call setResized(csdisq, nsam + 1_IK)
51call setResized(membership, nsam)
52call setResized(center, [ndim, ncls])
53call setResized(potential, ncls)
54call setResized(size, ncls)
55
56call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
57failed
58F
59niter
60+1
61disq
62+0.00000000, +0.00000000, +0.00000000, +0.275146402E-1, +0.538097806E-1, +0.267473757, +0.123414040, +0.538097806E-1
63csdisq
64+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.343835279E-1, +2.39235544, +3.14661646, +3.14661646, +4.29500771
65membership
66+4, +3, +1, +2, +5, +2, +2, +5
67potential
68+0.00000000, +0.788644493, +0.00000000, +0.00000000, +0.215239123
69center
70+2.91745758, +0.710206509, +2.03410673, +3.45510840, +4.75870800
71size
72+1, +3, +1, +1, +2
73
74
75ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
76[ndim, nsam, ncls]
77+5, +7, +4
78sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
79sample
80+0.469803810E-1, +2.73589039, +4.41585493, +0.727153420, +2.62369585, +1.11789441, +1.87587357
81+3.46518493, +4.19699478, +2.18554473, +4.56806946, +3.66272402, +3.48705959, +3.84169102
82+1.99363613, +3.11102891, +0.148693025, +2.53067231, +1.17997670, +0.909391344, +4.08227491
83+2.98503137, +4.04010487, +2.36521268, +0.914018452, +2.98129821, +4.13123798, +3.59145927
84+2.41045833, +4.36546135, +0.135650039, +1.33044147, +4.11749935, +0.595386624, +2.83048511
85call setResized(disq, nsam)
86call setResized(csdisq, nsam + 1_IK)
87call setResized(membership, nsam)
88call setResized(center, [ndim, ncls])
89call setResized(potential, ncls)
90call setResized(size, ncls)
91
92call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
93failed
94F
95niter
96+1
97disq
98+1.62158608, +1.30238795, +0.00000000, +0.00000000, +1.30238819, +4.17535830, +4.66269016
99csdisq
100+0.00000000, +0.00000000, +5.20955276, +27.9016762, +27.9016762, +27.9016762, +34.8328743, +43.2260704
101membership
102+1, +3, +4, +2, +3, +1, +1
103potential
104+15.3243933, +0.00000000, +5.20955276, +0.00000000
105center
106+1.01358283, +0.727153420, +2.67979312, +4.41585493
107+3.59797859, +4.56806946, +3.92985940, +2.18554473
108+2.32843423, +2.53067231, +2.14550281, +0.148693025
109+3.56924295, +0.914018452, +3.51070166, +2.36521268
110+1.94544351, +1.33044147, +4.24148035, +0.135650039
111size
112+3, +1, +2, +1
113
114
115ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
116[ndim, nsam, ncls]
117+5, +3, +2
118sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
119sample
120+1.09774446, +4.89227819, +4.58643341
121+1.28245616, +4.16129065, +3.24558687
122+3.62755585, +4.83825779, +2.18615103
123+2.41229582, +0.319649577, +2.79636765
124+3.34439254, +0.498896837E-1, +3.59288716
125call setResized(disq, nsam)
126call setResized(csdisq, nsam + 1_IK)
127call setResized(membership, nsam)
128call setResized(center, [ndim, ncls])
129call setResized(potential, ncls)
130call setResized(size, ncls)
131
132call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
133failed
134F
135niter
136+1
137disq
138+4.57793522, +0.00000000, +4.57793570
139csdisq
140+0.00000000, +18.3117428, +44.9644318, +44.9644318
141membership
142+1, +2, +1
143potential
144+18.3117428, +0.00000000
145center
146+2.84208894, +4.89227819
147+2.26402140, +4.16129065
148+2.90685344, +4.83825779
149+2.60433173, +0.319649577
150+3.46863985, +0.498896837E-1
151size
152+2, +1
153
154
155ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
156[ndim, nsam, ncls]
157+4, +1, +1
158sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
159sample
160+4.11432219
161+3.73472953
162+0.461795926
163+3.41511512
164call setResized(disq, nsam)
165call setResized(csdisq, nsam + 1_IK)
166call setResized(membership, nsam)
167call setResized(center, [ndim, ncls])
168call setResized(potential, ncls)
169call setResized(size, ncls)
170
171call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
172failed
173F
174niter
175+1
176disq
177+0.00000000
178csdisq
179+0.00000000, +0.00000000
180membership
181+1
182potential
183+0.00000000
184center
185+4.11432219
186+3.73472953
187+0.461795926
188+3.41511512
189size
190+1
191
192
193ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
194[ndim, nsam, ncls]
195+1, +3, +2
196sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
197sample
198+4.75599241, +2.13216352, +3.73829889
199call setResized(disq, nsam)
200call setResized(csdisq, nsam + 1_IK)
201call setResized(membership, nsam)
202call setResized(center, [ndim, ncls])
203call setResized(potential, ncls)
204call setResized(size, ncls)
205
206call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
207failed
208F
209niter
210+1
211disq
212+0.258925021, +0.00000000, +0.258925021
213csdisq
214+0.00000000, +0.00000000, +6.88447809, +7.92017841
215membership
216+1, +2, +1
217potential
218+1.03570008, +0.00000000
219center
220+4.24714565, +2.13216352
221size
222+2, +1
223
224
225ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
226[ndim, nsam, ncls]
227+4, +8, +5
228sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
229sample
230+0.902241170, +1.21660054, +0.526398122, +0.639582574, +2.76477456, +4.75375938, +0.833216310E-1, +1.46838212
231+4.17369270, +0.762449801, +4.25817156, +2.11539268, +2.62207699, +4.54568768, +0.680652559, +4.71295786
232+0.565861464, +3.02655268, +4.04550982, +0.945316553, +2.94427657, +1.17018247, +4.12222910, +4.81162930
233+1.57048488, +0.627495050, +0.624884665, +3.18371582, +0.962584615, +1.04376769, +0.693204105, +4.24750137
234call setResized(disq, nsam)
235call setResized(csdisq, nsam + 1_IK)
236call setResized(membership, nsam)
237call setResized(center, [ndim, ncls])
238call setResized(potential, ncls)
239call setResized(size, ncls)
240
241call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
242failed
243F
244niter
245+1
246disq
247+1.76302242, +0.623959064, +2.25347304, +1.76302207, +2.25347233, +0.00000000, +0.623959124, +0.00000000
248csdisq
249+0.00000000, +7.05208874, +9.54792500, +22.5534363, +22.5534363, +33.3637772, +33.3637772, +33.3637772, +33.3637772
250membership
251+4, +1, +5, +4, +5, +2, +1, +3
252potential
253+2.49583626, +0.00000000, +0.00000000, +7.05208874, +9.01389027
254center
255+0.649961114, +4.75375938, +1.46838212, +0.770911872, +1.64558637
256+0.721551180, +4.54568768, +4.71295786, +3.14454269, +3.44012427
257+3.57439089, +1.17018247, +4.81162930, +0.755589008, +3.49489307
258+0.660349607, +1.04376769, +4.24750137, +2.37710047, +0.793734670
259size
260+2, +1, +1, +2, +2
261
262
263ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
264[ndim, nsam, ncls]
265+4, +2, +1
266sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
267sample
268+1.38009012, +2.38166952
269+3.22566700, +4.63457537
270+2.31684113, +4.32816744
271+4.27526045, +1.37325311
272call setResized(disq, nsam)
273call setResized(csdisq, nsam + 1_IK)
274call setResized(membership, nsam)
275call setResized(center, [ndim, ncls])
276call setResized(potential, ncls)
277call setResized(size, ncls)
278
279call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
280failed
281F
282niter
283+1
284disq
285+3.86381578, +3.86381626
286csdisq
287+0.00000000, +0.00000000, +15.4552650
288membership
289+1, +1
290potential
291+15.4552650
292center
293+1.88087988
294+3.93012118
295+3.32250428
296+2.82425690
297size
298+2
299
300
301ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
302[ndim, nsam, ncls]
303+3, +5, +5
304sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
305sample
306+4.17415190, +1.61431932, +0.246993899, +3.67080092, +0.259948969
307+1.99562073, +2.48159909, +1.11275339, +2.78574252, +0.413348973
308+4.22451878, +0.876610279, +4.62803841, +2.45714283, +3.37406683
309call setResized(disq, nsam)
310call setResized(csdisq, nsam + 1_IK)
311call setResized(membership, nsam)
312call setResized(center, [ndim, ncls])
313call setResized(potential, ncls)
314call setResized(size, ncls)
315
316call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
317failed
318F
319niter
320+1
321disq
322+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
323csdisq
324+0.00000000, +0.00000000, +0.00000000, +2.06177902, +2.06177902, +2.06177902
325membership
326+2, +3, +5, +4, +1
327potential
328+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
329center
330+0.259948969, +4.17415190, +1.61431932, +3.67080092, +0.246993899
331+0.413348973, +1.99562073, +2.48159909, +2.78574252, +1.11275339
332+3.37406683, +4.22451878, +0.876610279, +2.45714283, +4.62803841
333size
334+1, +1, +1, +1, +1
335
336
337ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
338[ndim, nsam, ncls]
339+2, +7, +5
340sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
341sample
342+2.92859769, +2.81345630, +1.68383336, +4.10582447, +4.04435062, +2.28779602, +3.89948130
343+1.53934121, +0.998155773, +1.97992420, +2.51024294, +4.32186508, +0.838081241, +1.66061461
344call setResized(disq, nsam)
345call setResized(csdisq, nsam + 1_IK)
346call setResized(membership, nsam)
347call setResized(center, [ndim, ncls])
348call setResized(potential, ncls)
349call setResized(size, ncls)
350
351call setKmeans(rngf, membership, disq, csdisq, sample, center, size, potential, failed, niter) ! compute the new clusters and memberships.
352failed
353F
354niter
355+1
356disq
357+0.235013276, +0.348635092E-1, +0.00000000, +0.00000000, +0.00000000, +0.233614594, +0.00000000
358csdisq
359+0.00000000, +0.306139201, +0.306139201, +0.306139201, +1.07058501, +1.07058501, +1.37252760, +1.37252760
360membership
361+3, +3, +4, +5, +2, +3, +1
362potential
363+0.00000000, +0.00000000, +0.608081818, +0.00000000, +0.00000000
364center
365+3.89948130, +4.04435062, +2.67661667, +1.68383336, +4.10582447
366+1.66061461, +4.32186508, +1.12519288, +1.97992420, +2.51024294
367size
368+1, +1, +3, +1, +1
369
370

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6import glob
7import sys
8
9
10fontsize = 17
11
12prefixes = { "setKmeansPP" : "Memberships based on k-means++."
13 , "setKmeans" : "k-means refinement of k-means++ memberships."
14 }
15for prefix in list(prefixes.keys()):
16
17 legends = []
18 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
19 ax = plt.subplot()
20
21 fileList = glob.glob(prefix + ".*.txt")
22 if len(fileList) == 2:
23
24 for file in fileList:
25
26 kind = file.split(".")[1]
27 df = pd.read_csv(file, delimiter = ",", header = None)
28
29 if kind == "center":
30 ax.scatter ( df.values[:, 1]
31 , df.values[:,2]
32 , zorder = 100
33 , marker = "*"
34 , c = "red"
35 , s = 50
36 )
37 legends.append("center")
38 elif kind == "sample":
39 ax.scatter ( df.values[:, 1]
40 , df.values[:,2]
41 , c = df.values[:, 0]
42 , s = 10
43 )
44 legends.append("sample")
45 else:
46 sys.exit("Ambiguous file exists: {}".format(file))
47
48 ax.legend(legends, fontsize = fontsize)
49 plt.xticks(fontsize = fontsize - 2)
50 plt.yticks(fontsize = fontsize - 2)
51 ax.set_xlabel("X", fontsize = 17)
52 ax.set_ylabel("Y", fontsize = 17)
53 ax.set_title(prefixes[prefix], fontsize = fontsize)
54
55 plt.axis('equal')
56 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
57 ax.tick_params(axis = "y", which = "minor")
58 ax.tick_params(axis = "x", which = "minor")
59 ax.set_axisbelow(True)
60 plt.tight_layout()
61
62 plt.savefig(prefix + ".png")
63
64 else:
65
66 sys.exit("Ambiguous file list exists.")

Visualization of the example output


Test:
test_pm_clusKmeans


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1661 of file pm_clusKmeans.F90.


The documentation for this interface was generated from the following file: