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

Compute and return an asymptotically optimal set of cluster centers for the input sample, cluster membership IDs, and sample distances-squared from their corresponding cluster centers.
More...

Detailed Description

Compute and return an asymptotically optimal set of cluster centers for the input sample, cluster membership IDs, and sample distances-squared from their corresponding cluster centers.

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.
[out]membership: The 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).
[out]disq: The 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.
[out]csdisq: The output vector of shape (1:nsam+1) of the same type and kind as the input argument sample, containing the cumulative sum of the Euclidean squared distance of each input sample in sample from its nearest cluster center.
While the output contents are mostly useless, this argument can aid the algorithm efficiency to resolving the need for internal space allocation.
This potential speed-up is particularly relevant when the procedure is called repeatedly many times on samples of the same size.
[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]ncls: The input scalar of type integer of default kind IK, containing the number of the desired clusters to be identified in the sample.
(optional, default = ubound(center, 2). It must be present if and only if the output arguments center, size, and potential are all missing.)
[out]center: The 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 random cluster centers (centroids) selected from the input sample based on the computed memberships and minimum sample-cluster distances disq.
(optional. If missing, no cluster center information will be output.)
[out]size: The 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.
(optional. If missing, no cluster size information will be output.)
[out]potential: The 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.
(optional. If missing, no cluster potential information will be output.)


Possible calling interfaces

call setKmeansPP(rng, membership(1 : nsam), disq(1 : nsam), csdisq(0 : nsam), sample(1 : ndim, 1 : nsam), ncls)
call setKmeansPP(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))
Compute and return an asymptotically optimal set of cluster centers for the input sample,...
This module contains procedures and routines for the computing the Kmeans clustering of a given set o...
Warning
The condition ubound(center, rank(center)) > 0 must hold for the corresponding input arguments.
The condition ubound(sample, rank(sample)) == size(disq, 1) must hold for the corresponding input arguments.
The condition ubound(sample, rank(sample)) == size(csdisq, 1) - 1 must hold for the corresponding input arguments.
The condition ubound(center, rank(center)) == size(size, 1) must hold for the corresponding input arguments.
The condition ubound(center, rank(center)) == size(potential, 1) must hold for the corresponding input arguments.
The condition ubound(sample, rank(sample)) == size(membership, 1) must hold for the corresponding input arguments.
The condition ubound(center, rank(center)) <= size(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.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
By definition, the number of points in the input sample must be larger than the specified number of clusters.
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 similar to setMember, with the major difference being that setKmeansPP simultaneously computes the new cluster centers and sample memberships, whereas setMember computes the new sample memberships based on a given set of cluster centers.
The output of setKmeansPP can be directly passed to setCenter to the learn the new updated cluster centers and their sizes.
Note
Dropping the optional arguments can aid runtime performance.
This is particularly relevant when the output of this generic interface is directly passed to the k-means algorithm.
See also
setKmeans
setCenter
setMember
setKmeansPP
Arthur, D.; Vassilvitskii, S. (2007). k-means++: the advantages of careful seeding


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
9
10 implicit none
11
12 integer(IK) :: ndim, nsam, ncls, itry
13 real(RKG) , allocatable :: sample(:,:), center(:,:), disq(:), csdisq(:), potential(:)
14 integer(IK) , allocatable :: membership(:), size(:)
15 type(display_type) :: disp
16
17 disp = display_type(file = "main.out.F90")
18
19 call disp%skip
20 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%show("! Compute cluster centers based on an input sample and cluster memberships and member-center distances.")
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%skip
24
25 do itry = 1, 10
26 call disp%skip()
27 call disp%show("ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);")
28 ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
29 call disp%show("[ndim, nsam, ncls]")
30 call disp%show( [ndim, nsam, ncls] )
31 call disp%show("sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.")
32 sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
33 call disp%show("sample")
34 call disp%show( sample )
35 call disp%show("call setResized(disq, nsam)")
36 call setResized(disq, nsam)
37 call disp%show("call setResized(csdisq, nsam + 1_IK)")
38 call setResized(csdisq, nsam + 1_IK)
39 call disp%show("call setResized(membership, nsam)")
40 call setResized(membership, nsam)
41 call disp%show("call setResized(center, [ndim, ncls])")
42 call setResized(center, [ndim, ncls])
43 call disp%show("call setResized(potential, ncls)")
44 call setResized(potential, ncls)
45 call disp%show("call setResized(size, ncls)")
46 call setResized(size, ncls)
47 call disp%skip()
48
49 call disp%show("call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.")
50 call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
51 call disp%show("disq")
52 call disp%show( disq )
53 call disp%show("csdisq")
54 call disp%show( csdisq )
55 call disp%show("membership")
56 call disp%show( membership )
57 call disp%show("potential")
58 call disp%show( potential )
59 call disp%show("center")
60 call disp%show( center )
61 call disp%show("size")
62 call disp%show( size )
63 call disp%skip()
64 end do
65
66 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
67 ! Output an example for visualization.
68 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
69
70 block
71 integer(IK) :: funit, i
72 ndim = 2
73 ncls = 5
74 nsam = 5000
75 center = getUnifRand(0., 1., ndim, ncls)
76 sample = getUnifRand(0., 1., ndim, nsam)
77 call setResized(csdisq, nsam + 1_IK)
78 call setResized(disq, nsam)
79 call setResized(membership, nsam)
80 call setResized(center, [ndim, ncls])
81 call setResized(potential, ncls)
82 call setResized(size, ncls)
83 call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential)
84 open(newunit = funit, file = "setKmeansPP.center.txt")
85 do i = 1, ncls
86 write(funit, "(*(g0,:,','))") i, center(:,i)
87 end do
88 close(funit)
89 open(newunit = funit, file = "setKmeansPP.sample.txt")
90 do i = 1, nsam
91 write(funit, "(*(g0,:,','))") membership(i), sample(:,i)
92 end do
93 close(funit)
94 end block
95
96end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
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+5, +5, +3
10sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
11sample
12+2.61626530, +2.49449944, +2.88740635, +0.938760340, +3.10102725
13+1.56099772, +2.95452976, +1.49571240, +3.10879445, +2.55161071
14+3.58589792, +0.811363161, +1.74135542, +0.506854951, +2.02165461
15+1.31248772, +0.703998208, +1.62052417, +2.34893680, +0.116751492
16+4.57497692, +3.05059314, +3.31685781, +3.02566719, +0.682346523
17call setResized(disq, nsam)
18call setResized(csdisq, nsam + 1_IK)
19call setResized(membership, nsam)
20call setResized(center, [ndim, ncls])
21call setResized(potential, ncls)
22call setResized(size, ncls)
23
24call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
25disq
26+5.15786695, +4.05832624, +0.00000000, +0.00000000, +0.00000000
27csdisq
28+0.00000000, +5.15786695, +9.21619320, +9.21619320, +17.7548161, +17.7548161
29membership
30+1, +1, +1, +3, +2
31potential
32+9.21619320, +0.00000000, +0.00000000
33center
34+2.88740635, +3.10102725, +0.938760340
35+1.49571240, +2.55161071, +3.10879445
36+1.74135542, +2.02165461, +0.506854951
37+1.62052417, +0.116751492, +2.34893680
38+3.31685781, +0.682346523, +3.02566719
39size
40+3, +1, +1
41
42
43ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
44[ndim, nsam, ncls]
45+2, +3, +3
46sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
47sample
48+1.68761706, +0.359880924, +1.75377822
49+4.90363646, +1.16889060, +2.63341904
50call setResized(disq, nsam)
51call setResized(csdisq, nsam + 1_IK)
52call setResized(membership, nsam)
53call setResized(center, [ndim, ncls])
54call setResized(potential, ncls)
55call setResized(size, ncls)
56
57call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
58disq
59+0.00000000, +0.00000000, +0.00000000
60csdisq
61+0.00000000, +0.00000000, +0.00000000, +4.08779335
62membership
63+2, +1, +3
64potential
65+0.00000000, +0.00000000, +0.00000000
66center
67+0.359880924, +1.68761706, +1.75377822
68+1.16889060, +4.90363646, +2.63341904
69size
70+1, +1, +1
71
72
73ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
74[ndim, nsam, ncls]
75+2, +4, +3
76sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
77sample
78+1.10077357, +0.920268893, +2.86390972, +1.84875309
79+2.30324936, +3.61482835, +4.15435028, +1.69605255
80call setResized(disq, nsam)
81call setResized(csdisq, nsam + 1_IK)
82call setResized(membership, nsam)
83call setResized(center, [ndim, ncls])
84call setResized(potential, ncls)
85call setResized(size, ncls)
86
87call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
88disq
89+0.928161263, +0.00000000, +0.00000000, +0.00000000
90csdisq
91+0.00000000, +0.928161263, +4.99698496, +4.99698496, +4.99698496
92membership
93+1, +3, +2, +1
94potential
95+0.928161263, +0.00000000, +0.00000000
96center
97+1.84875309, +2.86390972, +0.920268893
98+1.69605255, +4.15435028, +3.61482835
99size
100+2, +1, +1
101
102
103ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
104[ndim, nsam, ncls]
105+1, +4, +4
106sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
107sample
108+2.21469545, +2.30991077, +3.10806012, +3.73275876
109call setResized(disq, nsam)
110call setResized(csdisq, nsam + 1_IK)
111call setResized(membership, nsam)
112call setResized(center, [ndim, ncls])
113call setResized(potential, ncls)
114call setResized(size, ncls)
115
116call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
117disq
118+0.00000000, +0.00000000, +0.00000000, +0.00000000
119csdisq
120+0.00000000, +0.906595774E-2, +0.906595774E-2, +0.906595774E-2, +0.906595774E-2
121membership
122+4, +1, +3, +2
123potential
124+0.00000000, +0.00000000, +0.00000000, +0.00000000
125center
126+2.30991077, +3.73275876, +3.10806012, +2.21469545
127size
128+1, +1, +1, +1
129
130
131ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
132[ndim, nsam, ncls]
133+5, +4, +4
134sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
135sample
136+0.615631938, +1.07353866, +1.92217255, +3.95547628
137+2.04954028, +3.28535891, +4.07291746, +1.32997990
138+1.57302237, +1.27022123, +2.39697719, +0.370814502
139+1.63032830, +4.09974623, +1.69385195, +1.72365546
140+4.81616783, +4.82765388, +3.71876240, +4.30589342
141call setResized(disq, nsam)
142call setResized(csdisq, nsam + 1_IK)
143call setResized(membership, nsam)
144call setResized(center, [ndim, ncls])
145call setResized(potential, ncls)
146call setResized(size, ncls)
147
148call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
149disq
150+0.00000000, +0.00000000, +0.00000000, +0.00000000
151csdisq
152+0.00000000, +0.00000000, +0.00000000, +7.68833923, +7.68833923
153membership
154+1, +3, +4, +2
155potential
156+0.00000000, +0.00000000, +0.00000000, +0.00000000
157center
158+0.615631938, +3.95547628, +1.07353866, +1.92217255
159+2.04954028, +1.32997990, +3.28535891, +4.07291746
160+1.57302237, +0.370814502, +1.27022123, +2.39697719
161+1.63032830, +1.72365546, +4.09974623, +1.69385195
162+4.81616783, +4.30589342, +4.82765388, +3.71876240
163size
164+1, +1, +1, +1
165
166
167ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
168[ndim, nsam, ncls]
169+3, +6, +3
170sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
171sample
172+3.47035861, +1.23233342, +0.206613839, +4.36533022, +2.02480912, +1.81483841
173+4.39511585, +4.52171373, +0.604673028, +4.33753872, +1.73817015, +0.811356902
174+4.07394934, +0.347383320, +1.02019429, +4.81816816, +0.164374709E-1, +3.34907794
175call setResized(disq, nsam)
176call setResized(csdisq, nsam + 1_IK)
177call setResized(membership, nsam)
178call setResized(center, [ndim, ncls])
179call setResized(potential, ncls)
180call setResized(size, ncls)
181
182call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
183disq
184+1.35815096, +0.00000000, +0.00000000, +0.00000000, +5.59817791, +8.05280304
185csdisq
186+0.00000000, +1.35815096, +18.2061348, +18.2061348, +18.2061348, +23.8043137, +31.8571167
187membership
188+2, +3, +1, +2, +1, +1
189potential
190+13.6509809, +1.35815096, +0.00000000
191center
192+0.206613839, +4.36533022, +1.23233342
193+0.604673028, +4.33753872, +4.52171373
194+1.02019429, +4.81816816, +0.347383320
195size
196+3, +2, +1
197
198
199ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
200[ndim, nsam, ncls]
201+5, +7, +5
202sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
203sample
204+3.01253295, +4.42124319, +0.192546248, +1.75626183, +3.62300348, +2.25239635, +4.98670197
205+2.40991235, +0.520440340, +3.22777724, +0.323043168, +4.63835859, +0.551990569, +3.64964962
206+4.85423803, +0.503692031, +4.59394169, +1.41908431, +0.408355296, +3.99546576, +3.17770386
207+2.04802322, +4.97627258, +0.434995592, +3.86973143, +3.43769431, +0.598253608, +3.80662227
208+3.37010598, +4.15741968, +3.55607939, +3.75267959, +2.18137693, +3.98265505, +0.906212628
209call setResized(disq, nsam)
210call setResized(csdisq, nsam + 1_IK)
211call setResized(membership, nsam)
212call setResized(center, [ndim, ncls])
213call setResized(potential, ncls)
214call setResized(size, ncls)
215
216call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
217disq
218+0.00000000, +0.00000000, +11.3254271, +0.00000000, +0.00000000, +7.24421930, +0.00000000
219csdisq
220+0.00000000, +0.00000000, +9.36728191, +20.6927090, +20.6927090, +20.6927090, +27.9369278, +27.9369278
221membership
222+1, +5, +1, +2, +3, +1, +4
223potential
224+18.5696468, +0.00000000, +0.00000000, +0.00000000, +0.00000000
225center
226+3.01253295, +1.75626183, +3.62300348, +4.98670197, +4.42124319
227+2.40991235, +0.323043168, +4.63835859, +3.64964962, +0.520440340
228+4.85423803, +1.41908431, +0.408355296, +3.17770386, +0.503692031
229+2.04802322, +3.86973143, +3.43769431, +3.80662227, +4.97627258
230+3.37010598, +3.75267959, +2.18137693, +0.906212628, +4.15741968
231size
232+3, +1, +1, +1, +1
233
234
235ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
236[ndim, nsam, ncls]
237+4, +4, +2
238sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
239sample
240+4.25245190, +3.03357100, +2.39169240, +2.59643316
241+3.91888857, +2.00881863, +0.877840519, +0.764200091E-1
242+2.38328004, +3.87669468, +1.32884407, +4.41034889
243+3.53272057, +3.48084307, +4.12918758, +3.49399662
244call setResized(disq, nsam)
245call setResized(csdisq, nsam + 1_IK)
246call setResized(membership, nsam)
247call setResized(center, [ndim, ncls])
248call setResized(potential, ncls)
249call setResized(size, ncls)
250
251call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
252disq
253+7.36701679, +0.00000000, +0.00000000, +4.21021366
254csdisq
255+0.00000000, +7.36701679, +7.36701679, +15.9700308, +20.1802444
256membership
257+1, +1, +2, +1
258potential
259+11.5772305, +0.00000000
260center
261+3.03357100, +2.39169240
262+2.00881863, +0.877840519
263+3.87669468, +1.32884407
264+3.48084307, +4.12918758
265size
266+3, +1
267
268
269ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
270[ndim, nsam, ncls]
271+3, +4, +4
272sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
273sample
274+2.72424674, +4.91062355, +2.37947321, +2.49682641
275+1.08171606, +3.72852230, +2.85086894, +4.93677711
276+2.00969887, +2.22752786, +0.800605416, +2.70743179
277call setResized(disq, nsam)
278call setResized(csdisq, nsam + 1_IK)
279call setResized(membership, nsam)
280call setResized(center, [ndim, ncls])
281call setResized(potential, ncls)
282call setResized(size, ncls)
283
284call setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
285disq
286+0.00000000, +0.00000000, +0.00000000, +0.00000000
287csdisq
288+0.00000000, +0.00000000, +0.00000000, +4.71067762, +4.71067762
289membership
290+1, +2, +4, +3
291potential
292+0.00000000, +0.00000000, +0.00000000, +0.00000000
293center
294+2.72424674, +4.91062355, +2.49682641, +2.37947321
295+1.08171606, +3.72852230, +4.93677711, +2.85086894
296+2.00969887, +2.22752786, +2.70743179, +0.800605416
297size
298+1, +1, +1, +1
299
300
301ndim = getUnifRand(1, 5); ncls = getUnifRand(1, 5); nsam = getUnifRand(ncls, 2 * ncls);
302[ndim, nsam, ncls]
303+3, +10, +5
304sample = getUnifRand(0., 5., ndim, nsam) ! Create a random sample.
305sample
306+2.80958200, +2.73564911, +4.47637177, +1.69878721, +1.40944481, +1.50349712, +4.23320818, +4.12399673, +2.96782517, +2.52528906
307+4.34577370, +1.60551965, +1.52460754, +2.11627674, +2.19805050, +4.96329308, +1.74042082, +2.37020826, +3.42810845, +2.34280229
308+2.90490389, +4.74109697, +0.929734111, +0.979961157, +1.89935803, +3.81704617, +1.64498448, +3.47271585, +3.30787802, +0.431272686
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 setKmeansPP(rngf, membership, disq, csdisq, sample, center, size, potential) ! compute the new clusters and memberships.
317disq
318+0.00000000, +0.00000000, +0.617286980, +0.00000000, +0.935696602, +2.91919136, +0.00000000, +0.00000000, +1.02953863, +1.03547812
319csdisq
320+0.00000000, +0.00000000, +0.00000000, +7.30596399, +7.30596399, +8.24166107, +11.1608524, +14.9100132, +14.9100132, +15.9395523, +16.9750309
321membership
322+4, +1, +5, +3, +3, +4, +5, +2, +4, +3
323potential
324+0.00000000, +0.00000000, +1.97117472, +3.94872999, +0.617286980
325center
326+2.73564911, +4.12399673, +1.69878721, +2.80958200, +4.23320818
327+1.60551965, +2.37020826, +2.11627674, +4.34577370, +1.74042082
328+4.74109697, +3.47271585, +0.979961157, +2.90490389, +1.64498448
329size
330+1, +1, +3, +3, +2
331
332

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
8import os
9
10fontsize = 17
11fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
12ax = plt.subplot()
13
14parent = os.path.basename(os.path.dirname(__file__))
15pattern = parent + "*.txt"
16
17fileList = glob.glob(pattern)
18legends = []
19if len(fileList) == 2:
20 for file in fileList:
21
22 kind = file.split(".")[1]
23 prefix = file.split(".")[0]
24 df = pd.read_csv(file, delimiter = ",", header = None)
25
26 if kind == "center":
27 ax.scatter ( df.values[:, 1]
28 , df.values[:,2]
29 , zorder = 100
30 , marker = "*"
31 , c = "red"
32 , s = 50
33 )
34 legends.append("center")
35 elif kind == "sample":
36 ax.scatter ( df.values[:, 1]
37 , df.values[:,2]
38 , c = df.values[:, 0]
39 , s = 10
40 )
41 legends.append("sample")
42 else:
43 sys.exit("Ambiguous file exists: {}".format(file))
44
45 ax.legend(legends, fontsize = fontsize)
46 plt.xticks(fontsize = fontsize - 2)
47 plt.yticks(fontsize = fontsize - 2)
48 ax.set_xlabel("X", fontsize = 17)
49 ax.set_ylabel("Y", fontsize = 17)
50 ax.set_title("Membership Scatter Plot", fontsize = fontsize)
51
52 plt.axis('equal')
53 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
54 ax.tick_params(axis = "y", which = "minor")
55 ax.tick_params(axis = "x", which = "minor")
56 ax.set_axisbelow(True)
57 plt.tight_layout()
58
59 plt.savefig(prefix + ".png")
60else:
61 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 at Austin

Definition at line 1181 of file pm_clusKmeans.F90.


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