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, +5
10point = getUnifRand(1, 10, ndim, npnt)
11point
12+9.00000000, +8.00000000, +1.00000000, +5.00000000, +2.00000000
13distance = getDisMatEuclid(point)
14distance
15+0.00000000, +1.00000000, +8.00000000, +4.00000000, +7.00000000
16+1.00000000, +0.00000000, +7.00000000, +3.00000000, +6.00000000
17+8.00000000, +7.00000000, +0.00000000, +4.00000000, +1.00000000
18+4.00000000, +3.00000000, +4.00000000, +0.00000000, +3.00000000
19+7.00000000, +6.00000000, +1.00000000, +3.00000000, +0.00000000
20distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
21distance
22+1.00000000, +1.00000000, +8.00000000, +4.00000000, +7.00000000
23+8.00000000, +7.00000000, +7.00000000, +3.00000000, +6.00000000
24+4.00000000, +3.00000000, +4.00000000, +4.00000000, +1.00000000
25+7.00000000, +6.00000000, +1.00000000, +3.00000000, +3.00000000
26distance = getDisMatEuclid(point, euclid)
27distance
28+0.00000000, +1.00000000, +8.00000000, +4.00000000, +7.00000000
29+1.00000000, +0.00000000, +7.00000000, +3.00000000, +6.00000000
30+8.00000000, +7.00000000, +0.00000000, +4.00000000, +1.00000000
31+4.00000000, +3.00000000, +4.00000000, +0.00000000, +3.00000000
32+7.00000000, +6.00000000, +1.00000000, +3.00000000, +0.00000000
33distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
34distance
35+1.00000000, +1.00000000, +8.00000000, +4.00000000, +7.00000000
36+8.00000000, +7.00000000, +7.00000000, +3.00000000, +6.00000000
37+4.00000000, +3.00000000, +4.00000000, +4.00000000, +1.00000000
38+7.00000000, +6.00000000, +1.00000000, +3.00000000, +3.00000000
39distance = getDisMatEuclid(point, euclidu)
40distance
41+0.00000000, +1.00000000, +8.00000000, +4.00000000, +7.00000000
42+1.00000000, +0.00000000, +7.00000000, +3.00000000, +6.00000000
43+8.00000000, +7.00000000, +0.00000000, +4.00000000, +1.00000000
44+4.00000000, +3.00000000, +4.00000000, +0.00000000, +3.00000000
45+7.00000000, +6.00000000, +1.00000000, +3.00000000, +0.00000000
46distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
47distance
48+1.00000000, +1.00000000, +64.0000000, +16.0000000, +49.0000000
49+64.0000000, +49.0000000, +49.0000000, +9.00000000, +36.0000000
50+16.0000000, +9.00000000, +16.0000000, +16.0000000, +1.00000000
51+49.0000000, +36.0000000, +1.00000000, +9.00000000, +9.00000000
52distance = getDisMatEuclid(point, euclidsq)
53distance
54+0.00000000, +1.00000000, +64.0000000, +16.0000000, +49.0000000
55+1.00000000, +0.00000000, +49.0000000, +9.00000000, +36.0000000
56+64.0000000, +49.0000000, +0.00000000, +16.0000000, +1.00000000
57+16.0000000, +9.00000000, +16.0000000, +0.00000000, +9.00000000
58+49.0000000, +36.0000000, +1.00000000, +9.00000000, +0.00000000
59distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
60distance
61+1.00000000, +1.00000000, +64.0000000, +16.0000000, +49.0000000
62+64.0000000, +49.0000000, +49.0000000, +9.00000000, +36.0000000
63+16.0000000, +9.00000000, +16.0000000, +16.0000000, +1.00000000
64+49.0000000, +36.0000000, +1.00000000, +9.00000000, +9.00000000
65
66
67ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
68[ndim, npnt]
69+3, +5
70point = getUnifRand(1, 10, ndim, npnt)
71point
72+4.00000000, +3.00000000, +9.00000000, +9.00000000, +9.00000000
73+7.00000000, +2.00000000, +8.00000000, +7.00000000, +1.00000000
74+10.0000000, +8.00000000, +5.00000000, +5.00000000, +4.00000000
75distance = getDisMatEuclid(point)
76distance
77+0.00000000, +5.47722530, +7.14142847, +7.07106781, +9.84885788
78+5.47722530, +0.00000000, +9.00000000, +8.36660004, +7.28011036
79+7.14142847, +9.00000000, +0.00000000, +1.00000000, +7.07106781
80+7.07106781, +8.36660004, +1.00000000, +0.00000000, +6.08276224
81+9.84885788, +7.28011036, +7.07106781, +6.08276224, +0.00000000
82distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
83distance
84+5.47722530, +5.47722530, +7.14142847, +7.07106781, +9.84885788
85+7.14142847, +9.00000000, +9.00000000, +8.36660004, +7.28011036
86+7.07106781, +8.36660004, +1.00000000, +1.00000000, +7.07106781
87+9.84885788, +7.28011036, +7.07106781, +6.08276224, +6.08276224
88distance = getDisMatEuclid(point, euclid)
89distance
90+0.00000000, +5.47722530, +7.14142847, +7.07106781, +9.84885788
91+5.47722530, +0.00000000, +9.00000000, +8.36660004, +7.28011036
92+7.14142847, +9.00000000, +0.00000000, +1.00000000, +7.07106781
93+7.07106781, +8.36660004, +1.00000000, +0.00000000, +6.08276224
94+9.84885788, +7.28011036, +7.07106781, +6.08276224, +0.00000000
95distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
96distance
97+5.47722530, +5.47722530, +7.14142847, +7.07106781, +9.84885788
98+7.14142847, +9.00000000, +9.00000000, +8.36660004, +7.28011036
99+7.07106781, +8.36660004, +1.00000000, +1.00000000, +7.07106781
100+9.84885788, +7.28011036, +7.07106781, +6.08276224, +6.08276224
101distance = getDisMatEuclid(point, euclidu)
102distance
103+0.00000000, +5.47722578, +7.14142847, +7.07106781, +9.84885788
104+5.47722578, +0.00000000, +9.00000000, +8.36660004, +7.28010988
105+7.14142847, +9.00000000, +0.00000000, +1.00000000, +7.07106781
106+7.07106781, +8.36660004, +1.00000000, +0.00000000, +6.08276272
107+9.84885788, +7.28010988, +7.07106781, +6.08276272, +0.00000000
108distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
109distance
110+30.0000000, +30.0000000, +51.0000000, +50.0000000, +97.0000000
111+51.0000000, +81.0000000, +81.0000000, +70.0000000, +53.0000000
112+50.0000000, +70.0000000, +1.00000000, +1.00000000, +50.0000000
113+97.0000000, +53.0000000, +50.0000000, +37.0000000, +37.0000000
114distance = getDisMatEuclid(point, euclidsq)
115distance
116+0.00000000, +30.0000000, +51.0000000, +50.0000000, +97.0000000
117+30.0000000, +0.00000000, +81.0000000, +70.0000000, +53.0000000
118+51.0000000, +81.0000000, +0.00000000, +1.00000000, +50.0000000
119+50.0000000, +70.0000000, +1.00000000, +0.00000000, +37.0000000
120+97.0000000, +53.0000000, +50.0000000, +37.0000000, +0.00000000
121distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
122distance
123+30.0000000, +30.0000000, +51.0000000, +50.0000000, +97.0000000
124+51.0000000, +81.0000000, +81.0000000, +70.0000000, +53.0000000
125+50.0000000, +70.0000000, +1.00000000, +1.00000000, +50.0000000
126+97.0000000, +53.0000000, +50.0000000, +37.0000000, +37.0000000
127
128
129ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
130[ndim, npnt]
131+1, +6
132point = getUnifRand(1, 10, ndim, npnt)
133point
134+8.00000000, +8.00000000, +6.00000000, +8.00000000, +5.00000000, +5.00000000
135distance = getDisMatEuclid(point)
136distance
137+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
138+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
139+2.00000000, +2.00000000, +0.00000000, +2.00000000, +1.00000000, +1.00000000
140+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
141+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
142+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
143distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
144distance
145+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
146+2.00000000, +2.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
147+0.00000000, +0.00000000, +2.00000000, +2.00000000, +1.00000000, +1.00000000
148+3.00000000, +3.00000000, +1.00000000, +3.00000000, +3.00000000, +3.00000000
149+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
150distance = getDisMatEuclid(point, euclid)
151distance
152+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
153+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
154+2.00000000, +2.00000000, +0.00000000, +2.00000000, +1.00000000, +1.00000000
155+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
156+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
157+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
158distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
159distance
160+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
161+2.00000000, +2.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
162+0.00000000, +0.00000000, +2.00000000, +2.00000000, +1.00000000, +1.00000000
163+3.00000000, +3.00000000, +1.00000000, +3.00000000, +3.00000000, +3.00000000
164+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
165distance = getDisMatEuclid(point, euclidu)
166distance
167+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
168+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
169+2.00000000, +2.00000000, +0.00000000, +2.00000000, +1.00000000, +1.00000000
170+0.00000000, +0.00000000, +2.00000000, +0.00000000, +3.00000000, +3.00000000
171+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
172+3.00000000, +3.00000000, +1.00000000, +3.00000000, +0.00000000, +0.00000000
173distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
174distance
175+0.00000000, +0.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
176+4.00000000, +4.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
177+0.00000000, +0.00000000, +4.00000000, +4.00000000, +1.00000000, +1.00000000
178+9.00000000, +9.00000000, +1.00000000, +9.00000000, +9.00000000, +9.00000000
179+9.00000000, +9.00000000, +1.00000000, +9.00000000, +0.00000000, +0.00000000
180distance = getDisMatEuclid(point, euclidsq)
181distance
182+0.00000000, +0.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
183+0.00000000, +0.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
184+4.00000000, +4.00000000, +0.00000000, +4.00000000, +1.00000000, +1.00000000
185+0.00000000, +0.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
186+9.00000000, +9.00000000, +1.00000000, +9.00000000, +0.00000000, +0.00000000
187+9.00000000, +9.00000000, +1.00000000, +9.00000000, +0.00000000, +0.00000000
188distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
189distance
190+0.00000000, +0.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
191+4.00000000, +4.00000000, +4.00000000, +0.00000000, +9.00000000, +9.00000000
192+0.00000000, +0.00000000, +4.00000000, +4.00000000, +1.00000000, +1.00000000
193+9.00000000, +9.00000000, +1.00000000, +9.00000000, +9.00000000, +9.00000000
194+9.00000000, +9.00000000, +1.00000000, +9.00000000, +0.00000000, +0.00000000
195
196
197ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
198[ndim, npnt]
199+2, +5
200point = getUnifRand(1, 10, ndim, npnt)
201point
202+7.00000000, +5.00000000, +3.00000000, +8.00000000, +1.00000000
203+2.00000000, +1.00000000, +2.00000000, +7.00000000, +6.00000000
204distance = getDisMatEuclid(point)
205distance
206+0.00000000, +2.23606801, +4.00000000, +5.09901953, +7.21110249
207+2.23606801, +0.00000000, +2.23606801, +6.70820427, +6.40312433
208+4.00000000, +2.23606801, +0.00000000, +7.07106781, +4.47213602
209+5.09901953, +6.70820427, +7.07106781, +0.00000000, +7.07106781
210+7.21110249, +6.40312433, +4.47213602, +7.07106781, +0.00000000
211distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
212distance
213+2.23606801, +2.23606801, +4.00000000, +5.09901953, +7.21110249
214+4.00000000, +2.23606801, +2.23606801, +6.70820427, +6.40312433
215+5.09901953, +6.70820427, +7.07106781, +7.07106781, +4.47213602
216+7.21110249, +6.40312433, +4.47213602, +7.07106781, +7.07106781
217distance = getDisMatEuclid(point, euclid)
218distance
219+0.00000000, +2.23606801, +4.00000000, +5.09901953, +7.21110249
220+2.23606801, +0.00000000, +2.23606801, +6.70820427, +6.40312433
221+4.00000000, +2.23606801, +0.00000000, +7.07106781, +4.47213602
222+5.09901953, +6.70820427, +7.07106781, +0.00000000, +7.07106781
223+7.21110249, +6.40312433, +4.47213602, +7.07106781, +0.00000000
224distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
225distance
226+2.23606801, +2.23606801, +4.00000000, +5.09901953, +7.21110249
227+4.00000000, +2.23606801, +2.23606801, +6.70820427, +6.40312433
228+5.09901953, +6.70820427, +7.07106781, +7.07106781, +4.47213602
229+7.21110249, +6.40312433, +4.47213602, +7.07106781, +7.07106781
230distance = getDisMatEuclid(point, euclidu)
231distance
232+0.00000000, +2.23606801, +4.00000000, +5.09901953, +7.21110249
233+2.23606801, +0.00000000, +2.23606801, +6.70820379, +6.40312433
234+4.00000000, +2.23606801, +0.00000000, +7.07106781, +4.47213602
235+5.09901953, +6.70820379, +7.07106781, +0.00000000, +7.07106781
236+7.21110249, +6.40312433, +4.47213602, +7.07106781, +0.00000000
237distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
238distance
239+5.00000000, +5.00000000, +16.0000000, +26.0000000, +52.0000000
240+16.0000000, +5.00000000, +5.00000000, +45.0000000, +41.0000000
241+26.0000000, +45.0000000, +50.0000000, +50.0000000, +20.0000000
242+52.0000000, +41.0000000, +20.0000000, +50.0000000, +50.0000000
243distance = getDisMatEuclid(point, euclidsq)
244distance
245+0.00000000, +5.00000000, +16.0000000, +26.0000000, +52.0000000
246+5.00000000, +0.00000000, +5.00000000, +45.0000000, +41.0000000
247+16.0000000, +5.00000000, +0.00000000, +50.0000000, +20.0000000
248+26.0000000, +45.0000000, +50.0000000, +0.00000000, +50.0000000
249+52.0000000, +41.0000000, +20.0000000, +50.0000000, +0.00000000
250distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
251distance
252+5.00000000, +5.00000000, +16.0000000, +26.0000000, +52.0000000
253+16.0000000, +5.00000000, +5.00000000, +45.0000000, +41.0000000
254+26.0000000, +45.0000000, +50.0000000, +50.0000000, +20.0000000
255+52.0000000, +41.0000000, +20.0000000, +50.0000000, +50.0000000
256
257
258ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
259[ndim, npnt]
260+3, +1
261point = getUnifRand(1, 10, ndim, npnt)
262point
263+10.0000000
264+2.00000000
265+1.00000000
266distance = getDisMatEuclid(point)
267distance
268+0.00000000
269distance = getDisMatEuclid(rdpack, uppLow, point) ! drop the zero-valued diagonal elements of the distance matrix.
270distance
271distance = getDisMatEuclid(point, euclid)
272distance
273+0.00000000
274distance = getDisMatEuclid(rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
275distance
276distance = getDisMatEuclid(point, euclidu)
277distance
278+0.00000000
279distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
280distance
281distance = getDisMatEuclid(point, euclidsq)
282distance
283+0.00000000
284distance = getDisMatEuclid(rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
285distance
286
287
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 Austin

Definition at line 3110 of file pm_distanceEuclid.F90.


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