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, +1
10point = getUnifRand(1, 10, ndim, npnt)
11point
12+4.00000000
13+3.00000000
14call setResized(distance, [npnt, npnt])
15call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
16distance
17+0.00000000
18call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
19distance(1:npnt-1, 1:npnt)
20call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
21distance
22+0.00000000
23call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
24distance(1:npnt-1, 1:npnt)
25call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
26distance
27+0.00000000
28call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
29distance(1:npnt-1, 1:npnt)
30
31
32ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
33[ndim, npnt]
34+1, +6
35point = getUnifRand(1, 10, ndim, npnt)
36point
37+3.00000000, +9.00000000, +8.00000000, +2.00000000, +7.00000000, +1.00000000
38call setResized(distance, [npnt, npnt])
39call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
40distance
41+0.00000000, +6.00000000, +5.00000000, +1.00000000, +4.00000000, +2.00000000
42+6.00000000, +0.00000000, +1.00000000, +7.00000000, +2.00000000, +8.00000000
43+5.00000000, +1.00000000, +0.00000000, +6.00000000, +1.00000000, +7.00000000
44+1.00000000, +7.00000000, +6.00000000, +0.00000000, +5.00000000, +1.00000000
45+4.00000000, +2.00000000, +1.00000000, +5.00000000, +0.00000000, +6.00000000
46+2.00000000, +8.00000000, +7.00000000, +1.00000000, +6.00000000, +0.00000000
47call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
48distance(1:npnt-1, 1:npnt)
49+6.00000000, +6.00000000, +5.00000000, +1.00000000, +4.00000000, +2.00000000
50+5.00000000, +1.00000000, +1.00000000, +7.00000000, +2.00000000, +8.00000000
51+1.00000000, +7.00000000, +6.00000000, +6.00000000, +1.00000000, +7.00000000
52+4.00000000, +2.00000000, +1.00000000, +5.00000000, +5.00000000, +1.00000000
53+2.00000000, +8.00000000, +7.00000000, +1.00000000, +6.00000000, +6.00000000
54call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
55distance
56+0.00000000, +6.00000000, +5.00000000, +1.00000000, +4.00000000, +2.00000000
57+6.00000000, +0.00000000, +1.00000000, +7.00000000, +2.00000000, +8.00000000
58+5.00000000, +1.00000000, +0.00000000, +6.00000000, +1.00000000, +7.00000000
59+1.00000000, +7.00000000, +6.00000000, +0.00000000, +5.00000000, +1.00000000
60+4.00000000, +2.00000000, +1.00000000, +5.00000000, +0.00000000, +6.00000000
61+2.00000000, +8.00000000, +7.00000000, +1.00000000, +6.00000000, +0.00000000
62call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
63distance(1:npnt-1, 1:npnt)
64+36.0000000, +36.0000000, +25.0000000, +1.00000000, +16.0000000, +4.00000000
65+25.0000000, +1.00000000, +1.00000000, +49.0000000, +4.00000000, +64.0000000
66+1.00000000, +49.0000000, +36.0000000, +36.0000000, +1.00000000, +49.0000000
67+16.0000000, +4.00000000, +1.00000000, +25.0000000, +25.0000000, +1.00000000
68+4.00000000, +64.0000000, +49.0000000, +1.00000000, +36.0000000, +36.0000000
69call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
70distance
71+0.00000000, +36.0000000, +25.0000000, +1.00000000, +16.0000000, +4.00000000
72+36.0000000, +0.00000000, +1.00000000, +49.0000000, +4.00000000, +64.0000000
73+25.0000000, +1.00000000, +0.00000000, +36.0000000, +1.00000000, +49.0000000
74+1.00000000, +49.0000000, +36.0000000, +0.00000000, +25.0000000, +1.00000000
75+16.0000000, +4.00000000, +1.00000000, +25.0000000, +0.00000000, +36.0000000
76+4.00000000, +64.0000000, +49.0000000, +1.00000000, +36.0000000, +0.00000000
77call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
78distance(1:npnt-1, 1:npnt)
79+36.0000000, +36.0000000, +25.0000000, +1.00000000, +16.0000000, +4.00000000
80+25.0000000, +1.00000000, +1.00000000, +49.0000000, +4.00000000, +64.0000000
81+1.00000000, +49.0000000, +36.0000000, +36.0000000, +1.00000000, +49.0000000
82+16.0000000, +4.00000000, +1.00000000, +25.0000000, +25.0000000, +1.00000000
83+4.00000000, +64.0000000, +49.0000000, +1.00000000, +36.0000000, +36.0000000
84
85
86ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
87[ndim, npnt]
88+2, +5
89point = getUnifRand(1, 10, ndim, npnt)
90point
91+2.00000000, +7.00000000, +9.00000000, +5.00000000, +4.00000000
92+5.00000000, +8.00000000, +9.00000000, +8.00000000, +10.0000000
93call setResized(distance, [npnt, npnt])
94call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
95distance
96+0.00000000, +5.83095169, +8.06225777, +4.24264050, +5.38516474
97+5.83095169, +0.00000000, +2.23606801, +2.00000000, +3.60555124
98+8.06225777, +2.23606801, +0.00000000, +4.12310553, +5.09901953
99+4.24264050, +2.00000000, +4.12310553, +0.00000000, +2.23606801
100+5.38516474, +3.60555124, +5.09901953, +2.23606801, +0.00000000
101call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
102distance(1:npnt-1, 1:npnt)
103+5.83095169, +5.83095169, +8.06225777, +4.24264050, +5.38516474
104+8.06225777, +2.23606801, +2.23606801, +2.00000000, +3.60555124
105+4.24264050, +2.00000000, +4.12310553, +4.12310553, +5.09901953
106+5.38516474, +3.60555124, +5.09901953, +2.23606801, +2.23606801
107call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
108distance
109+0.00000000, +5.83095169, +8.06225777, +4.24264050, +5.38516474
110+5.83095169, +0.00000000, +2.23606801, +2.00000000, +3.60555124
111+8.06225777, +2.23606801, +0.00000000, +4.12310553, +5.09901953
112+4.24264050, +2.00000000, +4.12310553, +0.00000000, +2.23606801
113+5.38516474, +3.60555124, +5.09901953, +2.23606801, +0.00000000
114call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
115distance(1:npnt-1, 1:npnt)
116+34.0000000, +34.0000000, +65.0000000, +18.0000000, +29.0000000
117+65.0000000, +5.00000000, +5.00000000, +4.00000000, +13.0000000
118+18.0000000, +4.00000000, +17.0000000, +17.0000000, +26.0000000
119+29.0000000, +13.0000000, +26.0000000, +5.00000000, +5.00000000
120call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
121distance
122+0.00000000, +34.0000000, +65.0000000, +18.0000000, +29.0000000
123+34.0000000, +0.00000000, +5.00000000, +4.00000000, +13.0000000
124+65.0000000, +5.00000000, +0.00000000, +17.0000000, +26.0000000
125+18.0000000, +4.00000000, +17.0000000, +0.00000000, +5.00000000
126+29.0000000, +13.0000000, +26.0000000, +5.00000000, +0.00000000
127call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
128distance(1:npnt-1, 1:npnt)
129+34.0000000, +34.0000000, +65.0000000, +18.0000000, +29.0000000
130+65.0000000, +5.00000000, +5.00000000, +4.00000000, +13.0000000
131+18.0000000, +4.00000000, +17.0000000, +17.0000000, +26.0000000
132+29.0000000, +13.0000000, +26.0000000, +5.00000000, +5.00000000
133
134
135ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
136[ndim, npnt]
137+1, +1
138point = getUnifRand(1, 10, ndim, npnt)
139point
140+3.00000000
141call setResized(distance, [npnt, npnt])
142call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
143distance
144+0.00000000
145call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
146distance(1:npnt-1, 1:npnt)
147call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
148distance
149+0.00000000
150call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
151distance(1:npnt-1, 1:npnt)
152call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
153distance
154+0.00000000
155call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
156distance(1:npnt-1, 1:npnt)
157
158
159ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 7)
160[ndim, npnt]
161+3, +6
162point = getUnifRand(1, 10, ndim, npnt)
163point
164+1.00000000, +6.00000000, +4.00000000, +2.00000000, +7.00000000, +5.00000000
165+1.00000000, +9.00000000, +5.00000000, +5.00000000, +7.00000000, +8.00000000
166+3.00000000, +8.00000000, +2.00000000, +2.00000000, +5.00000000, +10.0000000
167call setResized(distance, [npnt, npnt])
168call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclid)
169distance
170+0.00000000, +10.6770782, +5.09901953, +4.24264050, +8.71779823, +10.6770782
171+10.6770782, +0.00000000, +7.48331451, +8.24621105, +3.74165726, +2.44948983
172+5.09901953, +7.48331451, +0.00000000, +2.00000000, +4.69041586, +8.60232544
173+4.24264050, +8.24621105, +2.00000000, +0.00000000, +6.16441345, +9.05538559
174+8.71779823, +3.74165726, +4.69041586, +6.16441345, +0.00000000, +5.47722578
175+10.6770782, +2.44948983, +8.60232544, +9.05538559, +5.47722578, +0.00000000
176call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclid) ! drop the zero-valued diagonal elements of the distance matrix.
177distance(1:npnt-1, 1:npnt)
178+10.6770782, +10.6770782, +5.09901953, +4.24264050, +8.71779823, +10.6770782
179+5.09901953, +7.48331451, +7.48331451, +8.24621105, +3.74165726, +2.44948983
180+4.24264050, +8.24621105, +2.00000000, +2.00000000, +4.69041586, +8.60232544
181+8.71779823, +3.74165726, +4.69041586, +6.16441345, +6.16441345, +9.05538559
182+10.6770782, +2.44948983, +8.60232544, +9.05538559, +5.47722578, +5.47722578
183call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidu)
184distance
185+0.00000000, +10.6770782, +5.09901953, +4.24264050, +8.71779823, +10.6770782
186+10.6770782, +0.00000000, +7.48331499, +8.24621105, +3.74165750, +2.44948983
187+5.09901953, +7.48331499, +0.00000000, +2.00000000, +4.69041586, +8.60232544
188+4.24264050, +8.24621105, +2.00000000, +0.00000000, +6.16441393, +9.05538559
189+8.71779823, +3.74165750, +4.69041586, +6.16441393, +0.00000000, +5.47722578
190+10.6770782, +2.44948983, +8.60232544, +9.05538559, +5.47722578, +0.00000000
191call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
192distance(1:npnt-1, 1:npnt)
193+114.000000, +114.000000, +26.0000000, +18.0000000, +76.0000000, +114.000000
194+26.0000000, +56.0000000, +56.0000000, +68.0000000, +14.0000000, +6.00000000
195+18.0000000, +68.0000000, +4.00000000, +4.00000000, +22.0000000, +74.0000000
196+76.0000000, +14.0000000, +22.0000000, +38.0000000, +38.0000000, +82.0000000
197+114.000000, +6.00000000, +74.0000000, +82.0000000, +30.0000000, +30.0000000
198call setDisMatEuclid(distance, rdpack, uppLowDia, point, euclidsq)
199distance
200+0.00000000, +114.000000, +26.0000000, +18.0000000, +76.0000000, +114.000000
201+114.000000, +0.00000000, +56.0000000, +68.0000000, +14.0000000, +6.00000000
202+26.0000000, +56.0000000, +0.00000000, +4.00000000, +22.0000000, +74.0000000
203+18.0000000, +68.0000000, +4.00000000, +0.00000000, +38.0000000, +82.0000000
204+76.0000000, +14.0000000, +22.0000000, +38.0000000, +0.00000000, +30.0000000
205+114.000000, +6.00000000, +74.0000000, +82.0000000, +30.0000000, +0.00000000
206call setDisMatEuclid(distance(1:npnt-1, 1:npnt), rdpack, uppLow, point, euclidsq) ! drop the zero-valued diagonal elements of the distance matrix.
207distance(1:npnt-1, 1:npnt)
208+114.000000, +114.000000, +26.0000000, +18.0000000, +76.0000000, +114.000000
209+26.0000000, +56.0000000, +56.0000000, +68.0000000, +14.0000000, +6.00000000
210+18.0000000, +68.0000000, +4.00000000, +4.00000000, +22.0000000, +74.0000000
211+76.0000000, +14.0000000, +22.0000000, +38.0000000, +38.0000000, +82.0000000
212+114.000000, +6.00000000, +74.0000000, +82.0000000, +30.0000000, +30.0000000
213
214
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 3449 of file pm_distanceEuclid.F90.


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