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+1, +6
10point = getUnifRand(1, 10, ndim, npnt)
11point
12+4.00000000, +10.0000000, +10.0000000, +4.00000000, +3.00000000, +1.00000000
13call setResized(distance, [npnt, npnt])
14call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
15distance
16+0.00000000, +6.00000000, +6.00000000, +0.00000000, +1.00000000, +3.00000000
17+6.00000000, +0.00000000, +0.00000000, +6.00000000, +7.00000000, +9.00000000
18+6.00000000, +0.00000000, +0.00000000, +6.00000000, +7.00000000, +9.00000000
19+0.00000000, +6.00000000, +6.00000000, +0.00000000, +1.00000000, +3.00000000
20+1.00000000, +7.00000000, +7.00000000, +1.00000000, +0.00000000, +2.00000000
21+3.00000000, +9.00000000, +9.00000000, +3.00000000, +2.00000000, +0.00000000
22call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
23distance(1:npnt-1, 1:npnt)
24+6.00000000, +6.00000000, +6.00000000, +0.00000000, +1.00000000, +3.00000000
25+6.00000000, +0.00000000, +0.00000000, +6.00000000, +7.00000000, +9.00000000
26+0.00000000, +6.00000000, +6.00000000, +6.00000000, +7.00000000, +9.00000000
27+1.00000000, +7.00000000, +7.00000000, +1.00000000, +1.00000000, +3.00000000
28+3.00000000, +9.00000000, +9.00000000, +3.00000000, +2.00000000, +2.00000000
29call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
30distance
31+0.00000000, +6.00000000, +6.00000000, +0.00000000, +1.00000000, +3.00000000
32+6.00000000, +0.00000000, +0.00000000, +6.00000000, +7.00000000, +9.00000000
33+6.00000000, +0.00000000, +0.00000000, +6.00000000, +7.00000000, +9.00000000
34+0.00000000, +6.00000000, +6.00000000, +0.00000000, +1.00000000, +3.00000000
35+1.00000000, +7.00000000, +7.00000000, +1.00000000, +0.00000000, +2.00000000
36+3.00000000, +9.00000000, +9.00000000, +3.00000000, +2.00000000, +0.00000000
37call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
38distance(1:npnt-1, 1:npnt)
39+36.0000000, +36.0000000, +36.0000000, +0.00000000, +1.00000000, +9.00000000
40+36.0000000, +0.00000000, +0.00000000, +36.0000000, +49.0000000, +81.0000000
41+0.00000000, +36.0000000, +36.0000000, +36.0000000, +49.0000000, +81.0000000
42+1.00000000, +49.0000000, +49.0000000, +1.00000000, +1.00000000, +9.00000000
43+9.00000000, +81.0000000, +81.0000000, +9.00000000, +4.00000000, +4.00000000
44call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
45distance
46+0.00000000, +36.0000000, +36.0000000, +0.00000000, +1.00000000, +9.00000000
47+36.0000000, +0.00000000, +0.00000000, +36.0000000, +49.0000000, +81.0000000
48+36.0000000, +0.00000000, +0.00000000, +36.0000000, +49.0000000, +81.0000000
49+0.00000000, +36.0000000, +36.0000000, +0.00000000, +1.00000000, +9.00000000
50+1.00000000, +49.0000000, +49.0000000, +1.00000000, +0.00000000, +4.00000000
51+9.00000000, +81.0000000, +81.0000000, +9.00000000, +4.00000000, +0.00000000
52call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
53distance(1:npnt-1, 1:npnt)
54+36.0000000, +36.0000000, +36.0000000, +0.00000000, +1.00000000, +9.00000000
55+36.0000000, +0.00000000, +0.00000000, +36.0000000, +49.0000000, +81.0000000
56+0.00000000, +36.0000000, +36.0000000, +36.0000000, +49.0000000, +81.0000000
57+1.00000000, +49.0000000, +49.0000000, +1.00000000, +1.00000000, +9.00000000
58+9.00000000, +81.0000000, +81.0000000, +9.00000000, +4.00000000, +4.00000000
59
60
61ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
62[ndim, npnt]
63+3, +1
64point = getUnifRand(1, 10, ndim, npnt)
65point
66+9.00000000
67+9.00000000
68+7.00000000
69call setResized(distance, [npnt, npnt])
70call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
71distance
72+0.00000000
73call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
74distance(1:npnt-1, 1:npnt)
75call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
76distance
77+0.00000000
78call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
79distance(1:npnt-1, 1:npnt)
80call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
81distance
82+0.00000000
83call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
84distance(1:npnt-1, 1:npnt)
85
86
87ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
88[ndim, npnt]
89+2, +7
90point = getUnifRand(1, 10, ndim, npnt)
91point
92+2.00000000, +2.00000000, +4.00000000, +6.00000000, +8.00000000, +1.00000000, +2.00000000
93+3.00000000, +4.00000000, +4.00000000, +1.00000000, +8.00000000, +3.00000000, +9.00000000
94call setResized(distance, [npnt, npnt])
95call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
96distance
97+0.00000000, +1.00000000, +2.23606801, +4.47213602, +7.81025028, +1.00000000, +6.00000000
98+1.00000000, +0.00000000, +2.00000000, +5.00000000, +7.21110249, +1.41421354, +5.00000000
99+2.23606801, +2.00000000, +0.00000000, +3.60555124, +5.65685415, +3.16227770, +5.38516474
100+4.47213602, +5.00000000, +3.60555124, +0.00000000, +7.28010988, +5.38516474, +8.94427204
101+7.81025028, +7.21110249, +5.65685415, +7.28010988, +0.00000000, +8.60232544, +6.08276224
102+1.00000000, +1.41421354, +3.16227770, +5.38516474, +8.60232544, +0.00000000, +6.08276224
103+6.00000000, +5.00000000, +5.38516474, +8.94427204, +6.08276224, +6.08276224, +0.00000000
104call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
105distance(1:npnt-1, 1:npnt)
106+1.00000000, +1.00000000, +2.23606801, +4.47213602, +7.81025028, +1.00000000, +6.00000000
107+2.23606801, +2.00000000, +2.00000000, +5.00000000, +7.21110249, +1.41421354, +5.00000000
108+4.47213602, +5.00000000, +3.60555124, +3.60555124, +5.65685415, +3.16227770, +5.38516474
109+7.81025028, +7.21110249, +5.65685415, +7.28010988, +7.28010988, +5.38516474, +8.94427204
110+1.00000000, +1.41421354, +3.16227770, +5.38516474, +8.60232544, +8.60232544, +6.08276224
111+6.00000000, +5.00000000, +5.38516474, +8.94427204, +6.08276224, +6.08276224, +6.08276224
112call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
113distance
114+0.00000000, +1.00000000, +2.23606801, +4.47213602, +7.81024981, +1.00000000, +6.00000000
115+1.00000000, +0.00000000, +2.00000000, +5.00000000, +7.21110249, +1.41421354, +5.00000000
116+2.23606801, +2.00000000, +0.00000000, +3.60555124, +5.65685415, +3.16227770, +5.38516474
117+4.47213602, +5.00000000, +3.60555124, +0.00000000, +7.28010988, +5.38516474, +8.94427204
118+7.81024981, +7.21110249, +5.65685415, +7.28010988, +0.00000000, +8.60232544, +6.08276272
119+1.00000000, +1.41421354, +3.16227770, +5.38516474, +8.60232544, +0.00000000, +6.08276272
120+6.00000000, +5.00000000, +5.38516474, +8.94427204, +6.08276272, +6.08276272, +0.00000000
121call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
122distance(1:npnt-1, 1:npnt)
123+1.00000000, +1.00000000, +5.00000000, +20.0000000, +61.0000000, +1.00000000, +36.0000000
124+5.00000000, +4.00000000, +4.00000000, +25.0000000, +52.0000000, +2.00000000, +25.0000000
125+20.0000000, +25.0000000, +13.0000000, +13.0000000, +32.0000000, +10.0000000, +29.0000000
126+61.0000000, +52.0000000, +32.0000000, +53.0000000, +53.0000000, +29.0000000, +80.0000000
127+1.00000000, +2.00000000, +10.0000000, +29.0000000, +74.0000000, +74.0000000, +37.0000000
128+36.0000000, +25.0000000, +29.0000000, +80.0000000, +37.0000000, +37.0000000, +37.0000000
129call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
130distance
131+0.00000000, +1.00000000, +5.00000000, +20.0000000, +61.0000000, +1.00000000, +36.0000000
132+1.00000000, +0.00000000, +4.00000000, +25.0000000, +52.0000000, +2.00000000, +25.0000000
133+5.00000000, +4.00000000, +0.00000000, +13.0000000, +32.0000000, +10.0000000, +29.0000000
134+20.0000000, +25.0000000, +13.0000000, +0.00000000, +53.0000000, +29.0000000, +80.0000000
135+61.0000000, +52.0000000, +32.0000000, +53.0000000, +0.00000000, +74.0000000, +37.0000000
136+1.00000000, +2.00000000, +10.0000000, +29.0000000, +74.0000000, +0.00000000, +37.0000000
137+36.0000000, +25.0000000, +29.0000000, +80.0000000, +37.0000000, +37.0000000, +0.00000000
138call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
139distance(1:npnt-1, 1:npnt)
140+1.00000000, +1.00000000, +5.00000000, +20.0000000, +61.0000000, +1.00000000, +36.0000000
141+5.00000000, +4.00000000, +4.00000000, +25.0000000, +52.0000000, +2.00000000, +25.0000000
142+20.0000000, +25.0000000, +13.0000000, +13.0000000, +32.0000000, +10.0000000, +29.0000000
143+61.0000000, +52.0000000, +32.0000000, +53.0000000, +53.0000000, +29.0000000, +80.0000000
144+1.00000000, +2.00000000, +10.0000000, +29.0000000, +74.0000000, +74.0000000, +37.0000000
145+36.0000000, +25.0000000, +29.0000000, +80.0000000, +37.0000000, +37.0000000, +37.0000000
146
147
148ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
149[ndim, npnt]
150+2, +4
151point = getUnifRand(1, 10, ndim, npnt)
152point
153+10.0000000, +1.00000000, +5.00000000, +1.00000000
154+5.00000000, +5.00000000, +3.00000000, +1.00000000
155call setResized(distance, [npnt, npnt])
156call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
157distance
158+0.00000000, +9.00000000, +5.38516474, +9.84885788
159+9.00000000, +0.00000000, +4.47213602, +4.00000000
160+5.38516474, +4.47213602, +0.00000000, +4.47213602
161+9.84885788, +4.00000000, +4.47213602, +0.00000000
162call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
163distance(1:npnt-1, 1:npnt)
164+9.00000000, +9.00000000, +5.38516474, +9.84885788
165+5.38516474, +4.47213602, +4.47213602, +4.00000000
166+9.84885788, +4.00000000, +4.47213602, +4.47213602
167call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
168distance
169+0.00000000, +9.00000000, +5.38516474, +9.84885788
170+9.00000000, +0.00000000, +4.47213602, +4.00000000
171+5.38516474, +4.47213602, +0.00000000, +4.47213602
172+9.84885788, +4.00000000, +4.47213602, +0.00000000
173call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
174distance(1:npnt-1, 1:npnt)
175+81.0000000, +81.0000000, +29.0000000, +97.0000000
176+29.0000000, +20.0000000, +20.0000000, +16.0000000
177+97.0000000, +16.0000000, +20.0000000, +20.0000000
178call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
179distance
180+0.00000000, +81.0000000, +29.0000000, +97.0000000
181+81.0000000, +0.00000000, +20.0000000, +16.0000000
182+29.0000000, +20.0000000, +0.00000000, +20.0000000
183+97.0000000, +16.0000000, +20.0000000, +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+81.0000000, +81.0000000, +29.0000000, +97.0000000
187+29.0000000, +20.0000000, +20.0000000, +16.0000000
188+97.0000000, +16.0000000, +20.0000000, +20.0000000
189
190
191ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
192[ndim, npnt]
193+3, +3
194point = getUnifRand(1, 10, ndim, npnt)
195point
196+4.00000000, +9.00000000, +5.00000000
197+4.00000000, +7.00000000, +2.00000000
198+5.00000000, +8.00000000, +1.00000000
199call setResized(distance, [npnt, npnt])
200call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
201distance
202+0.00000000, +6.55743837, +4.58257580
203+6.55743837, +0.00000000, +9.48683357
204+4.58257580, +9.48683357, +0.00000000
205call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
206distance(1:npnt-1, 1:npnt)
207+6.55743837, +6.55743837, +4.58257580
208+4.58257580, +9.48683357, +9.48683357
209call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
210distance
211+0.00000000, +6.55743837, +4.58257580
212+6.55743837, +0.00000000, +9.48683262
213+4.58257580, +9.48683262, +0.00000000
214call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
215distance(1:npnt-1, 1:npnt)
216+43.0000000, +43.0000000, +21.0000000
217+21.0000000, +90.0000000, +90.0000000
218call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
219distance
220+0.00000000, +43.0000000, +21.0000000
221+43.0000000, +0.00000000, +90.0000000
222+21.0000000, +90.0000000, +0.00000000
223call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
224distance(1:npnt-1, 1:npnt)
225+43.0000000, +43.0000000, +21.0000000
226+21.0000000, +90.0000000, +90.0000000
227
228
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: