ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distanceEuclid::setDisMatEuclid 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,out]distance: The input/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.
[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.


Possible calling interfaces

call setDisMatEuclid(distance(1:npnt, 1:npnt), pack, subset, point(1:ndim, 1:npnt), method) ! subset /= uppLow
call setDisMatEuclid(distance(1:npnt-1, 1:npnt), pack, subset, point(1:ndim, 1:npnt), method) ! subset = uppLow
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.
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: setDisMatEuclid, 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("call setResized(distance, [npnt, npnt])")
35 call setResized(distance, [npnt, npnt])
36 call disp%show("call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)")
37 call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
38 call disp%show("distance")
39 call disp%show( distance )
40 call disp%show("call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.")
41 call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
42 call disp%show("distance(1:npnt-1, 1:npnt)")
43 call disp%show( distance(1:npnt-1, 1:npnt) )
44 call disp%show("call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)")
45 call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
46 call disp%show("distance")
47 call disp%show( distance )
48 call disp%show("call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.")
49 call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
50 call disp%show("distance(1:npnt-1, 1:npnt)")
51 call disp%show( distance(1:npnt-1, 1:npnt) )
52 call disp%show("call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)")
53 call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
54 call disp%show("distance")
55 call disp%show( distance )
56 call disp%show("call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.")
57 call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
58 call disp%show("distance(1:npnt-1, 1:npnt)")
59 call disp%show( distance(1:npnt-1, 1:npnt) )
60 call disp%skip()
61 end do
62 end block
63
64end 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+2, +7
10point = getUnifRand(1, 10, ndim, npnt)
11point
12+4.00000000, +8.00000000, +9.00000000, +10.0000000, +8.00000000, +9.00000000, +10.0000000
13+6.00000000, +10.0000000, +5.00000000, +2.00000000, +3.00000000, +2.00000000, +3.00000000
14call setResized(distance, [npnt, npnt])
15call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
16distance
17+0.00000000, +5.65685415, +5.09901953, +7.21110249, +5.00000000, +6.40312433, +6.70820427
18+5.65685415, +0.00000000, +5.09901953, +8.24621105, +7.00000000, +8.06225777, +7.28010988
19+5.09901953, +5.09901953, +0.00000000, +3.16227770, +2.23606801, +3.00000000, +2.23606801
20+7.21110249, +8.24621105, +3.16227770, +0.00000000, +2.23606801, +1.00000000, +1.00000000
21+5.00000000, +7.00000000, +2.23606801, +2.23606801, +0.00000000, +1.41421354, +2.00000000
22+6.40312433, +8.06225777, +3.00000000, +1.00000000, +1.41421354, +0.00000000, +1.41421354
23+6.70820427, +7.28010988, +2.23606801, +1.00000000, +2.00000000, +1.41421354, +0.00000000
24call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
25distance(1:npnt-1, 1:npnt)
26+5.65685415, +5.65685415, +5.09901953, +7.21110249, +5.00000000, +6.40312433, +6.70820427
27+5.09901953, +5.09901953, +5.09901953, +8.24621105, +7.00000000, +8.06225777, +7.28010988
28+7.21110249, +8.24621105, +3.16227770, +3.16227770, +2.23606801, +3.00000000, +2.23606801
29+5.00000000, +7.00000000, +2.23606801, +2.23606801, +2.23606801, +1.00000000, +1.00000000
30+6.40312433, +8.06225777, +3.00000000, +1.00000000, +1.41421354, +1.41421354, +2.00000000
31+6.70820427, +7.28010988, +2.23606801, +1.00000000, +2.00000000, +1.41421354, +1.41421354
32call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
33distance
34+0.00000000, +5.65685415, +5.09901953, +7.21110249, +5.00000000, +6.40312433, +6.70820379
35+5.65685415, +0.00000000, +5.09901953, +8.24621105, +7.00000000, +8.06225777, +7.28010988
36+5.09901953, +5.09901953, +0.00000000, +3.16227770, +2.23606801, +3.00000000, +2.23606801
37+7.21110249, +8.24621105, +3.16227770, +0.00000000, +2.23606801, +1.00000000, +1.00000000
38+5.00000000, +7.00000000, +2.23606801, +2.23606801, +0.00000000, +1.41421354, +2.00000000
39+6.40312433, +8.06225777, +3.00000000, +1.00000000, +1.41421354, +0.00000000, +1.41421354
40+6.70820379, +7.28010988, +2.23606801, +1.00000000, +2.00000000, +1.41421354, +0.00000000
41call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
42distance(1:npnt-1, 1:npnt)
43+32.0000000, +32.0000000, +26.0000000, +52.0000000, +25.0000000, +41.0000000, +45.0000000
44+26.0000000, +26.0000000, +26.0000000, +68.0000000, +49.0000000, +65.0000000, +53.0000000
45+52.0000000, +68.0000000, +10.0000000, +10.0000000, +5.00000000, +9.00000000, +5.00000000
46+25.0000000, +49.0000000, +5.00000000, +5.00000000, +5.00000000, +1.00000000, +1.00000000
47+41.0000000, +65.0000000, +9.00000000, +1.00000000, +2.00000000, +2.00000000, +4.00000000
48+45.0000000, +53.0000000, +5.00000000, +1.00000000, +4.00000000, +2.00000000, +2.00000000
49call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
50distance
51+0.00000000, +32.0000000, +26.0000000, +52.0000000, +25.0000000, +41.0000000, +45.0000000
52+32.0000000, +0.00000000, +26.0000000, +68.0000000, +49.0000000, +65.0000000, +53.0000000
53+26.0000000, +26.0000000, +0.00000000, +10.0000000, +5.00000000, +9.00000000, +5.00000000
54+52.0000000, +68.0000000, +10.0000000, +0.00000000, +5.00000000, +1.00000000, +1.00000000
55+25.0000000, +49.0000000, +5.00000000, +5.00000000, +0.00000000, +2.00000000, +4.00000000
56+41.0000000, +65.0000000, +9.00000000, +1.00000000, +2.00000000, +0.00000000, +2.00000000
57+45.0000000, +53.0000000, +5.00000000, +1.00000000, +4.00000000, +2.00000000, +0.00000000
58call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
59distance(1:npnt-1, 1:npnt)
60+32.0000000, +32.0000000, +26.0000000, +52.0000000, +25.0000000, +41.0000000, +45.0000000
61+26.0000000, +26.0000000, +26.0000000, +68.0000000, +49.0000000, +65.0000000, +53.0000000
62+52.0000000, +68.0000000, +10.0000000, +10.0000000, +5.00000000, +9.00000000, +5.00000000
63+25.0000000, +49.0000000, +5.00000000, +5.00000000, +5.00000000, +1.00000000, +1.00000000
64+41.0000000, +65.0000000, +9.00000000, +1.00000000, +2.00000000, +2.00000000, +4.00000000
65+45.0000000, +53.0000000, +5.00000000, +1.00000000, +4.00000000, +2.00000000, +2.00000000
66
67
68ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
69[ndim, npnt]
70+3, +2
71point = getUnifRand(1, 10, ndim, npnt)
72point
73+9.00000000, +3.00000000
74+8.00000000, +7.00000000
75+7.00000000, +10.0000000
76call setResized(distance, [npnt, npnt])
77call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
78distance
79+0.00000000, +6.78233051
80+6.78233051, +0.00000000
81call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
82distance(1:npnt-1, 1:npnt)
83+6.78233051, +6.78233051
84call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
85distance
86+0.00000000, +6.78233004
87+6.78233004, +0.00000000
88call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
89distance(1:npnt-1, 1:npnt)
90+46.0000000, +46.0000000
91call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
92distance
93+0.00000000, +46.0000000
94+46.0000000, +0.00000000
95call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
96distance(1:npnt-1, 1:npnt)
97+46.0000000, +46.0000000
98
99
100ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
101[ndim, npnt]
102+2, +7
103point = getUnifRand(1, 10, ndim, npnt)
104point
105+10.0000000, +10.0000000, +4.00000000, +7.00000000, +9.00000000, +5.00000000, +6.00000000
106+7.00000000, +9.00000000, +7.00000000, +2.00000000, +8.00000000, +7.00000000, +3.00000000
107call setResized(distance, [npnt, npnt])
108call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
109distance
110+0.00000000, +2.00000000, +6.00000000, +5.83095169, +1.41421354, +5.00000000, +5.65685415
111+2.00000000, +0.00000000, +6.32455540, +7.61577320, +1.41421354, +5.38516474, +7.21110249
112+6.00000000, +6.32455540, +0.00000000, +5.83095169, +5.09901953, +1.00000000, +4.47213602
113+5.83095169, +7.61577320, +5.83095169, +0.00000000, +6.32455540, +5.38516474, +1.41421354
114+1.41421354, +1.41421354, +5.09901953, +6.32455540, +0.00000000, +4.12310553, +5.83095169
115+5.00000000, +5.38516474, +1.00000000, +5.38516474, +4.12310553, +0.00000000, +4.12310553
116+5.65685415, +7.21110249, +4.47213602, +1.41421354, +5.83095169, +4.12310553, +0.00000000
117call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
118distance(1:npnt-1, 1:npnt)
119+2.00000000, +2.00000000, +6.00000000, +5.83095169, +1.41421354, +5.00000000, +5.65685415
120+6.00000000, +6.32455540, +6.32455540, +7.61577320, +1.41421354, +5.38516474, +7.21110249
121+5.83095169, +7.61577320, +5.83095169, +5.83095169, +5.09901953, +1.00000000, +4.47213602
122+1.41421354, +1.41421354, +5.09901953, +6.32455540, +6.32455540, +5.38516474, +1.41421354
123+5.00000000, +5.38516474, +1.00000000, +5.38516474, +4.12310553, +4.12310553, +5.83095169
124+5.65685415, +7.21110249, +4.47213602, +1.41421354, +5.83095169, +4.12310553, +4.12310553
125call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
126distance
127+0.00000000, +2.00000000, +6.00000000, +5.83095169, +1.41421354, +5.00000000, +5.65685415
128+2.00000000, +0.00000000, +6.32455540, +7.61577320, +1.41421354, +5.38516474, +7.21110249
129+6.00000000, +6.32455540, +0.00000000, +5.83095169, +5.09901953, +1.00000000, +4.47213602
130+5.83095169, +7.61577320, +5.83095169, +0.00000000, +6.32455540, +5.38516474, +1.41421354
131+1.41421354, +1.41421354, +5.09901953, +6.32455540, +0.00000000, +4.12310553, +5.83095169
132+5.00000000, +5.38516474, +1.00000000, +5.38516474, +4.12310553, +0.00000000, +4.12310553
133+5.65685415, +7.21110249, +4.47213602, +1.41421354, +5.83095169, +4.12310553, +0.00000000
134call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
135distance(1:npnt-1, 1:npnt)
136+4.00000000, +4.00000000, +36.0000000, +34.0000000, +2.00000000, +25.0000000, +32.0000000
137+36.0000000, +40.0000000, +40.0000000, +58.0000000, +2.00000000, +29.0000000, +52.0000000
138+34.0000000, +58.0000000, +34.0000000, +34.0000000, +26.0000000, +1.00000000, +20.0000000
139+2.00000000, +2.00000000, +26.0000000, +40.0000000, +40.0000000, +29.0000000, +2.00000000
140+25.0000000, +29.0000000, +1.00000000, +29.0000000, +17.0000000, +17.0000000, +34.0000000
141+32.0000000, +52.0000000, +20.0000000, +2.00000000, +34.0000000, +17.0000000, +17.0000000
142call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
143distance
144+0.00000000, +4.00000000, +36.0000000, +34.0000000, +2.00000000, +25.0000000, +32.0000000
145+4.00000000, +0.00000000, +40.0000000, +58.0000000, +2.00000000, +29.0000000, +52.0000000
146+36.0000000, +40.0000000, +0.00000000, +34.0000000, +26.0000000, +1.00000000, +20.0000000
147+34.0000000, +58.0000000, +34.0000000, +0.00000000, +40.0000000, +29.0000000, +2.00000000
148+2.00000000, +2.00000000, +26.0000000, +40.0000000, +0.00000000, +17.0000000, +34.0000000
149+25.0000000, +29.0000000, +1.00000000, +29.0000000, +17.0000000, +0.00000000, +17.0000000
150+32.0000000, +52.0000000, +20.0000000, +2.00000000, +34.0000000, +17.0000000, +0.00000000
151call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
152distance(1:npnt-1, 1:npnt)
153+4.00000000, +4.00000000, +36.0000000, +34.0000000, +2.00000000, +25.0000000, +32.0000000
154+36.0000000, +40.0000000, +40.0000000, +58.0000000, +2.00000000, +29.0000000, +52.0000000
155+34.0000000, +58.0000000, +34.0000000, +34.0000000, +26.0000000, +1.00000000, +20.0000000
156+2.00000000, +2.00000000, +26.0000000, +40.0000000, +40.0000000, +29.0000000, +2.00000000
157+25.0000000, +29.0000000, +1.00000000, +29.0000000, +17.0000000, +17.0000000, +34.0000000
158+32.0000000, +52.0000000, +20.0000000, +2.00000000, +34.0000000, +17.0000000, +17.0000000
159
160
161ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
162[ndim, npnt]
163+3, +3
164point = getUnifRand(1, 10, ndim, npnt)
165point
166+2.00000000, +6.00000000, +4.00000000
167+6.00000000, +1.00000000, +1.00000000
168+8.00000000, +10.0000000, +6.00000000
169call setResized(distance, [npnt, npnt])
170call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
171distance
172+0.00000000, +6.70820427, +5.74456263
173+6.70820427, +0.00000000, +4.47213602
174+5.74456263, +4.47213602, +0.00000000
175call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
176distance(1:npnt-1, 1:npnt)
177+6.70820427, +6.70820427, +5.74456263
178+5.74456263, +4.47213602, +4.47213602
179call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
180distance
181+0.00000000, +6.70820379, +5.74456263
182+6.70820379, +0.00000000, +4.47213602
183+5.74456263, +4.47213602, +0.00000000
184call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
185distance(1:npnt-1, 1:npnt)
186+45.0000000, +45.0000000, +33.0000000
187+33.0000000, +20.0000000, +20.0000000
188call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
189distance
190+0.00000000, +45.0000000, +33.0000000
191+45.0000000, +0.00000000, +20.0000000
192+33.0000000, +20.0000000, +0.00000000
193call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
194distance(1:npnt-1, 1:npnt)
195+45.0000000, +45.0000000, +33.0000000
196+33.0000000, +20.0000000, +20.0000000
197
198
199ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
200[ndim, npnt]
201+3, +2
202point = getUnifRand(1, 10, ndim, npnt)
203point
204+10.0000000, +8.00000000
205+10.0000000, +4.00000000
206+10.0000000, +8.00000000
207call setResized(distance, [npnt, npnt])
208call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
209distance
210+0.00000000, +6.63324928
211+6.63324928, +0.00000000
212call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
213distance(1:npnt-1, 1:npnt)
214+6.63324928, +6.63324928
215call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
216distance
217+0.00000000, +6.63324976
218+6.63324976, +0.00000000
219call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
220distance(1:npnt-1, 1:npnt)
221+44.0000000, +44.0000000
222call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
223distance
224+0.00000000, +44.0000000
225+44.0000000, +0.00000000
226call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
227distance(1:npnt-1, 1:npnt)
228+44.0000000, +44.0000000
229
230
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 3449 of file pm_distanceEuclid.F90.


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