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

Return the full or a subset of the Euclidean (squared) distance matrix of the input set of npnt points in ndim dimensions. More...

Detailed Description

Return the full or a subset of the Euclidean (squared) distance matrix of the input set of npnt points in ndim dimensions.

Parameters
[in]pack: The input scalar that can be:
  1. the constant rdpack or an object of type rdpack_type, implying the use of Rectangular Default Packing format for the output matrix.
[in]subset: The input scalar that can be:
  1. the constant uppLowDia or an object of type uppLowDia_type, indicating that the output distance must contain the full distance matrix of shape (1:npnt, 1:npnt) including the zero-valued diagonals.
  2. the constant uppLow or an object of type uppLow_type, indicating that the output distance must exclude the zero-valued diagonals from the distance matrix yielding a distance matrix of shape (1:npnt - 1, 1:npnt).
    Motivation: The zero-valued diagonal elements of the distance matrix are are frequently troubling for subsequent vector operations on the output distance matrix.
    Such vector operations include but are not limited to finding the extrema of distances, for example, the nearest and farthest neighbors.
    This subset value offers a fast convenient method of excluding self-distance values from the output distance matrix such that each column (1:npnt-1 , i) of the distance matrix contains only the distances of point(1:ndim, i) with all other npnt - 1 points in point.
    For example, finding the nearest neighbor of the points using the output distance matrix would be as simple as minval(distance, 1).
    Finding the actual index of the point that is the nearest neighbor to each point would be slightly more involved as a two-step process:
    nn1loc(1 : npnt) = minloc(distance(1 : npnt - 1, 1 : npnt), 1)
    nn1loc = merge(nn1loc, nn1loc + 1, getRange(1, npnt) <= nn1loc)
    where nn1loc is the vector of indices of the first nearest neighbors such that point(:,nn1loc(i)) is the nearest neighbor to point(:,i).
[in]point: The input contiguous matrix of shape (1:ndim, 1:npnt) of,
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing npnt points in the ndim-dimensional Euclidean space whose distances with respect to each other must be computed and returned.
[in]method: The input scalar that can be,
  1. the constant euclid or an object of type euclid_type, implying that all distance calculations must be done without undue numerical overflow.
    This option is computationally the most expensive method.
  2. the constant euclidu or an object of type euclidu_type, implying that all distance calculations must be without runtime checks for numerical overflow.
    This option is computationally faster than the euclid method.
  3. the constant euclidsq or an object of type euclidsq_type implying that the squared values of all distance calculations must be returned without runtime checks for numerical overflow.
    This option is computationally the fastest approach to constructing the distance matrix because it avoid costly sqrt() operations and runtime overflow checks.
(optional, default = euclid)
Returns
distance : The output contiguous array of rank 2 of the same type and kind as the input argument point.
On output, it contains the requested subset of the (squared) distance matrix in the specified packing format pack.
Any element of distance that is not included in the specified subset will remain intact, if any such element exists.


Possible calling interfaces

distance(1:npnt, 1:npnt) = getDisMatEuclid(pack, subset, point(1:ndim, 1:npnt), method) ! subset = uppLowDia, pack = rdpack
distance(1:npnt-1, 1:npnt) = getDisMatEuclid(pack, subset, point(1:ndim, 1:npnt), method) ! subset = uppLow, pack = rdpack
!
Return the full or a subset of the Euclidean (squared) distance matrix of the input set of npnt point...
This module contains procedures and generic interfaces for computing the Euclidean norm of a single p...
Warning
The condition size(point, 1) == size(point, 2) must hold for the corresponding input arguments.
The condition shape(distance) == [size(point, 1), size(point, 1)] .or. .not. same_type_as(subset, uppLowDia) must hold for the corresponding input arguments.
The condition shape(distance) == [size(point, 1) - 1, size(point, 1)] .or. .not. same_type_as(subset, uppLow) 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.
Developer Remark:
The input arguments pack, subset appear first for a good reason: To allow the possibility of adding of similarly-named arguments for the input point matrix.
See also
euclid
euclidu
euclidsq
euclid_type
euclidu_type
euclidsq_type
getDisEuclid
setDisEuclid
getDisMatEuclid
setDisMatEuclid


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKH
4 use pm_kind, only: RKG => RKS ! all processor kinds are supported.
5 use pm_io, only: display_type
6 use pm_distanceEuclid, only: getDisMatEuclid, rdpack, uppLow, uppLowDia, euclid, euclidu, euclidsq
8 use pm_distUnif, only: getUnifRand
9
10 implicit none
11
12 integer(IK) :: ndim, npnt, itry, ntry = 5
13 type(display_type) :: disp
14 disp = display_type(file = "main.out.F90")
15
16 call disp%skip()
17 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
18 call disp%show("! Compute the distance matrix of a set of points.")
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%skip()
21
22 block
23 real(RKG), allocatable :: distance(:,:), point(:,:)
24 do itry = 1, ntry
25 call disp%skip()
26 call disp%show("ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)")
27 ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
28 call disp%show("[ndim, npnt]")
29 call disp%show( [ndim, npnt] )
30 call disp%show("point = getUnifRand(1, 10, ndim, npnt)")
31 point = getUnifRand(1, 10, ndim, npnt)
32 call disp%show("point")
33 call disp%show( point )
34 call disp%show("distance = getDisMatEuclid(point)")
35 distance = getDisMatEuclid(point)
36 call disp%show("distance")
37 call disp%show( distance )
38 call disp%show("distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.")
39 distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
40 call disp%show("distance")
41 call disp%show( distance )
42 call disp%show("distance = getDisMatEuclid(point, euclid)")
43 distance = getDisMatEuclid(point, euclid)
44 call disp%show("distance")
45 call disp%show( distance )
46 call disp%show("distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.")
47 distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
48 call disp%show("distance")
49 call disp%show( distance )
50 call disp%show("distance = getDisMatEuclid(point, euclidu)")
51 distance = getDisMatEuclid(point, euclidu)
52 call disp%show("distance")
53 call disp%show( distance )
54 call disp%show("distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.")
55 distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
56 call disp%show("distance")
57 call disp%show( distance )
58 call disp%show("distance = getDisMatEuclid(point, euclidsq)")
59 distance = getDisMatEuclid(point, euclidsq)
60 call disp%show("distance")
61 call disp%show( distance )
62 call disp%show("distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.")
63 distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
64 call disp%show("distance")
65 call disp%show( distance )
66 call disp%skip()
67 end do
68 end block
69
70end 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(euclidu_type), parameter euclidu
This is a scalar parameter object of type euclidu_typethat is exclusively used to request unsafe meth...
type(euclid_type), parameter euclid
This is a scalar parameter object of type euclid_type that is exclusively used to request safe method...
type(euclidsq_type), parameter euclidsq
This is a scalar parameter object of type euclidsq_typethat is exclusively used to request computing ...
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 RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
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 the distance matrix of a set of points.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
8[ndim, npnt]
9+1, +4
10point = getUnifRand(1, 10, ndim, npnt)
11point
12+6.00000000, +5.00000000, +2.00000000, +6.00000000
13distance = getDisMatEuclid(point)
14distance
15+0.00000000, +1.00000000, +4.00000000, +0.00000000
16+1.00000000, +0.00000000, +3.00000000, +1.00000000
17+4.00000000, +3.00000000, +0.00000000, +4.00000000
18+0.00000000, +1.00000000, +4.00000000, +0.00000000
19distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
20distance
21+1.00000000, +1.00000000, +4.00000000, +0.00000000
22+4.00000000, +3.00000000, +3.00000000, +1.00000000
23+0.00000000, +1.00000000, +4.00000000, +4.00000000
24distance = getDisMatEuclid(point, euclid)
25distance
26+0.00000000, +1.00000000, +4.00000000, +0.00000000
27+1.00000000, +0.00000000, +3.00000000, +1.00000000
28+4.00000000, +3.00000000, +0.00000000, +4.00000000
29+0.00000000, +1.00000000, +4.00000000, +0.00000000
30distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
31distance
32+1.00000000, +1.00000000, +4.00000000, +0.00000000
33+4.00000000, +3.00000000, +3.00000000, +1.00000000
34+0.00000000, +1.00000000, +4.00000000, +4.00000000
35distance = getDisMatEuclid(point, euclidu)
36distance
37+0.00000000, +1.00000000, +4.00000000, +0.00000000
38+1.00000000, +0.00000000, +3.00000000, +1.00000000
39+4.00000000, +3.00000000, +0.00000000, +4.00000000
40+0.00000000, +1.00000000, +4.00000000, +0.00000000
41distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
42distance
43+1.00000000, +1.00000000, +16.0000000, +0.00000000
44+16.0000000, +9.00000000, +9.00000000, +1.00000000
45+0.00000000, +1.00000000, +16.0000000, +16.0000000
46distance = getDisMatEuclid(point, euclidsq)
47distance
48+0.00000000, +1.00000000, +16.0000000, +0.00000000
49+1.00000000, +0.00000000, +9.00000000, +1.00000000
50+16.0000000, +9.00000000, +0.00000000, +16.0000000
51+0.00000000, +1.00000000, +16.0000000, +0.00000000
52distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
53distance
54+1.00000000, +1.00000000, +16.0000000, +0.00000000
55+16.0000000, +9.00000000, +9.00000000, +1.00000000
56+0.00000000, +1.00000000, +16.0000000, +16.0000000
57
58
59ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
60[ndim, npnt]
61+3, +5
62point = getUnifRand(1, 10, ndim, npnt)
63point
64+1.00000000, +6.00000000, +1.00000000, +9.00000000, +9.00000000
65+3.00000000, +1.00000000, +9.00000000, +7.00000000, +9.00000000
66+9.00000000, +8.00000000, +2.00000000, +8.00000000, +9.00000000
67distance = getDisMatEuclid(point)
68distance
69+0.00000000, +5.47722530, +9.21954536, +9.00000000, +10.0000000
70+5.47722530, +0.00000000, +11.1803398, +6.70820427, +8.60232544
71+9.21954536, +11.1803398, +0.00000000, +10.1980391, +10.6301460
72+9.00000000, +6.70820427, +10.1980391, +0.00000000, +2.23606801
73+10.0000000, +8.60232544, +10.6301460, +2.23606801, +0.00000000
74distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
75distance
76+5.47722530, +5.47722530, +9.21954536, +9.00000000, +10.0000000
77+9.21954536, +11.1803398, +11.1803398, +6.70820427, +8.60232544
78+9.00000000, +6.70820427, +10.1980391, +10.1980391, +10.6301460
79+10.0000000, +8.60232544, +10.6301460, +2.23606801, +2.23606801
80distance = getDisMatEuclid(point, euclid)
81distance
82+0.00000000, +5.47722530, +9.21954536, +9.00000000, +10.0000000
83+5.47722530, +0.00000000, +11.1803398, +6.70820427, +8.60232544
84+9.21954536, +11.1803398, +0.00000000, +10.1980391, +10.6301460
85+9.00000000, +6.70820427, +10.1980391, +0.00000000, +2.23606801
86+10.0000000, +8.60232544, +10.6301460, +2.23606801, +0.00000000
87distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
88distance
89+5.47722530, +5.47722530, +9.21954536, +9.00000000, +10.0000000
90+9.21954536, +11.1803398, +11.1803398, +6.70820427, +8.60232544
91+9.00000000, +6.70820427, +10.1980391, +10.1980391, +10.6301460
92+10.0000000, +8.60232544, +10.6301460, +2.23606801, +2.23606801
93distance = getDisMatEuclid(point, euclidu)
94distance
95+0.00000000, +5.47722578, +9.21954441, +9.00000000, +10.0000000
96+5.47722578, +0.00000000, +11.1803398, +6.70820379, +8.60232544
97+9.21954441, +11.1803398, +0.00000000, +10.1980391, +10.6301460
98+9.00000000, +6.70820379, +10.1980391, +0.00000000, +2.23606801
99+10.0000000, +8.60232544, +10.6301460, +2.23606801, +0.00000000
100distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
101distance
102+30.0000000, +30.0000000, +85.0000000, +81.0000000, +100.000000
103+85.0000000, +125.000000, +125.000000, +45.0000000, +74.0000000
104+81.0000000, +45.0000000, +104.000000, +104.000000, +113.000000
105+100.000000, +74.0000000, +113.000000, +5.00000000, +5.00000000
106distance = getDisMatEuclid(point, euclidsq)
107distance
108+0.00000000, +30.0000000, +85.0000000, +81.0000000, +100.000000
109+30.0000000, +0.00000000, +125.000000, +45.0000000, +74.0000000
110+85.0000000, +125.000000, +0.00000000, +104.000000, +113.000000
111+81.0000000, +45.0000000, +104.000000, +0.00000000, +5.00000000
112+100.000000, +74.0000000, +113.000000, +5.00000000, +0.00000000
113distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
114distance
115+30.0000000, +30.0000000, +85.0000000, +81.0000000, +100.000000
116+85.0000000, +125.000000, +125.000000, +45.0000000, +74.0000000
117+81.0000000, +45.0000000, +104.000000, +104.000000, +113.000000
118+100.000000, +74.0000000, +113.000000, +5.00000000, +5.00000000
119
120
121ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
122[ndim, npnt]
123+2, +2
124point = getUnifRand(1, 10, ndim, npnt)
125point
126+3.00000000, +7.00000000
127+4.00000000, +4.00000000
128distance = getDisMatEuclid(point)
129distance
130+0.00000000, +4.00000000
131+4.00000000, +0.00000000
132distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
133distance
134+4.00000000, +4.00000000
135distance = getDisMatEuclid(point, euclid)
136distance
137+0.00000000, +4.00000000
138+4.00000000, +0.00000000
139distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
140distance
141+4.00000000, +4.00000000
142distance = getDisMatEuclid(point, euclidu)
143distance
144+0.00000000, +4.00000000
145+4.00000000, +0.00000000
146distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
147distance
148+16.0000000, +16.0000000
149distance = getDisMatEuclid(point, euclidsq)
150distance
151+0.00000000, +16.0000000
152+16.0000000, +0.00000000
153distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
154distance
155+16.0000000, +16.0000000
156
157
158ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
159[ndim, npnt]
160+3, +5
161point = getUnifRand(1, 10, ndim, npnt)
162point
163+8.00000000, +2.00000000, +5.00000000, +6.00000000, +4.00000000
164+9.00000000, +5.00000000, +10.0000000, +1.00000000, +3.00000000
165+3.00000000, +3.00000000, +5.00000000, +5.00000000, +9.00000000
166distance = getDisMatEuclid(point)
167distance
168+0.00000000, +7.21110249, +3.74165726, +8.48528099, +9.38083172
169+7.21110249, +0.00000000, +6.16441345, +6.00000000, +6.63324928
170+3.74165726, +6.16441345, +0.00000000, +9.05538559, +8.12403870
171+8.48528099, +6.00000000, +9.05538559, +0.00000000, +4.89897966
172+9.38083172, +6.63324928, +8.12403870, +4.89897966, +0.00000000
173distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
174distance
175+7.21110249, +7.21110249, +3.74165726, +8.48528099, +9.38083172
176+3.74165726, +6.16441345, +6.16441345, +6.00000000, +6.63324928
177+8.48528099, +6.00000000, +9.05538559, +9.05538559, +8.12403870
178+9.38083172, +6.63324928, +8.12403870, +4.89897966, +4.89897966
179distance = getDisMatEuclid(point, euclid)
180distance
181+0.00000000, +7.21110249, +3.74165726, +8.48528099, +9.38083172
182+7.21110249, +0.00000000, +6.16441345, +6.00000000, +6.63324928
183+3.74165726, +6.16441345, +0.00000000, +9.05538559, +8.12403870
184+8.48528099, +6.00000000, +9.05538559, +0.00000000, +4.89897966
185+9.38083172, +6.63324928, +8.12403870, +4.89897966, +0.00000000
186distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
187distance
188+7.21110249, +7.21110249, +3.74165726, +8.48528099, +9.38083172
189+3.74165726, +6.16441345, +6.16441345, +6.00000000, +6.63324928
190+8.48528099, +6.00000000, +9.05538559, +9.05538559, +8.12403870
191+9.38083172, +6.63324928, +8.12403870, +4.89897966, +4.89897966
192distance = getDisMatEuclid(point, euclidu)
193distance
194+0.00000000, +7.21110249, +3.74165750, +8.48528099, +9.38083172
195+7.21110249, +0.00000000, +6.16441393, +6.00000000, +6.63324976
196+3.74165750, +6.16441393, +0.00000000, +9.05538559, +8.12403870
197+8.48528099, +6.00000000, +9.05538559, +0.00000000, +4.89897966
198+9.38083172, +6.63324976, +8.12403870, +4.89897966, +0.00000000
199distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
200distance
201+52.0000000, +52.0000000, +14.0000000, +72.0000000, +88.0000000
202+14.0000000, +38.0000000, +38.0000000, +36.0000000, +44.0000000
203+72.0000000, +36.0000000, +82.0000000, +82.0000000, +66.0000000
204+88.0000000, +44.0000000, +66.0000000, +24.0000000, +24.0000000
205distance = getDisMatEuclid(point, euclidsq)
206distance
207+0.00000000, +52.0000000, +14.0000000, +72.0000000, +88.0000000
208+52.0000000, +0.00000000, +38.0000000, +36.0000000, +44.0000000
209+14.0000000, +38.0000000, +0.00000000, +82.0000000, +66.0000000
210+72.0000000, +36.0000000, +82.0000000, +0.00000000, +24.0000000
211+88.0000000, +44.0000000, +66.0000000, +24.0000000, +0.00000000
212distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
213distance
214+52.0000000, +52.0000000, +14.0000000, +72.0000000, +88.0000000
215+14.0000000, +38.0000000, +38.0000000, +36.0000000, +44.0000000
216+72.0000000, +36.0000000, +82.0000000, +82.0000000, +66.0000000
217+88.0000000, +44.0000000, +66.0000000, +24.0000000, +24.0000000
218
219
220ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
221[ndim, npnt]
222+1, +5
223point = getUnifRand(1, 10, ndim, npnt)
224point
225+10.0000000, +3.00000000, +1.00000000, +9.00000000, +10.0000000
226distance = getDisMatEuclid(point)
227distance
228+0.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
229+7.00000000, +0.00000000, +2.00000000, +6.00000000, +7.00000000
230+9.00000000, +2.00000000, +0.00000000, +8.00000000, +9.00000000
231+1.00000000, +6.00000000, +8.00000000, +0.00000000, +1.00000000
232+0.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
233distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
234distance
235+7.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
236+9.00000000, +2.00000000, +2.00000000, +6.00000000, +7.00000000
237+1.00000000, +6.00000000, +8.00000000, +8.00000000, +9.00000000
238+0.00000000, +7.00000000, +9.00000000, +1.00000000, +1.00000000
239distance = getDisMatEuclid(point, euclid)
240distance
241+0.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
242+7.00000000, +0.00000000, +2.00000000, +6.00000000, +7.00000000
243+9.00000000, +2.00000000, +0.00000000, +8.00000000, +9.00000000
244+1.00000000, +6.00000000, +8.00000000, +0.00000000, +1.00000000
245+0.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
246distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
247distance
248+7.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
249+9.00000000, +2.00000000, +2.00000000, +6.00000000, +7.00000000
250+1.00000000, +6.00000000, +8.00000000, +8.00000000, +9.00000000
251+0.00000000, +7.00000000, +9.00000000, +1.00000000, +1.00000000
252distance = getDisMatEuclid(point, euclidu)
253distance
254+0.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
255+7.00000000, +0.00000000, +2.00000000, +6.00000000, +7.00000000
256+9.00000000, +2.00000000, +0.00000000, +8.00000000, +9.00000000
257+1.00000000, +6.00000000, +8.00000000, +0.00000000, +1.00000000
258+0.00000000, +7.00000000, +9.00000000, +1.00000000, +0.00000000
259distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
260distance
261+49.0000000, +49.0000000, +81.0000000, +1.00000000, +0.00000000
262+81.0000000, +4.00000000, +4.00000000, +36.0000000, +49.0000000
263+1.00000000, +36.0000000, +64.0000000, +64.0000000, +81.0000000
264+0.00000000, +49.0000000, +81.0000000, +1.00000000, +1.00000000
265distance = getDisMatEuclid(point, euclidsq)
266distance
267+0.00000000, +49.0000000, +81.0000000, +1.00000000, +0.00000000
268+49.0000000, +0.00000000, +4.00000000, +36.0000000, +49.0000000
269+81.0000000, +4.00000000, +0.00000000, +64.0000000, +81.0000000
270+1.00000000, +36.0000000, +64.0000000, +0.00000000, +1.00000000
271+0.00000000, +49.0000000, +81.0000000, +1.00000000, +0.00000000
272distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
273distance
274+49.0000000, +49.0000000, +81.0000000, +1.00000000, +0.00000000
275+81.0000000, +4.00000000, +4.00000000, +36.0000000, +49.0000000
276+1.00000000, +36.0000000, +64.0000000, +64.0000000, +81.0000000
277+0.00000000, +49.0000000, +81.0000000, +1.00000000, +1.00000000
278
279
Test:
test_pm_distanceEuclid
Todo:
High Priority: This generic interface must be extended to allow other packing and subsets of the output distance matrix.


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, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin

Definition at line 3110 of file pm_distanceEuclid.F90.


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