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

Generate and return the (squared) Euclidean distance of a (set of) point(s) in ndim-dimensions from a reference point (possibly origin), optionally robustly without underflow or overflow.
More...

Detailed Description

Generate and return the (squared) Euclidean distance of a (set of) point(s) in ndim-dimensions from a reference point (possibly origin), optionally robustly without underflow or overflow.

Parameters
[out]distance: The output object of,
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the requested Euclidean (squared) distance(s).
The rank and shape of the output distance follows that of the interfaces illustrated below.
[in]point: The input contiguous vector of shape (1:ndim) or matrix of shape (1:ndim, 1:npnt) of the same type and kind as the output distance, containing a (set of npnt) point(s) in the ndim-dimensional Euclidean space whose distances with respect to the input reference point ref must be returned.
[in]ref: The input contiguous vector of shape (1:ndim) or matrix of shape (1:ndim, 1:nref) of the same type and kind as point, containing the (set of nref) reference point(s) from which the distance(s) of point must be computed.
(optional, default = [(0., i = 1, size(point, 1))].)
[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 euclidv or an object of type euclidv_type, implying that the volumes corresponding to all Euclidean distances must be returned without runtime checks for numerical overflow.
    This option is computationally the fastest approach to computing the distances because it avoid costly sqrt() operations and runtime overflow
  4. 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 computing the distances because it avoid costly sqrt() operations and runtime overflow checks.
(optional, default = euclid)


Possible calling interfaces

! distance with respect to origin.
call setDisEuclid(distance, point(1:ndim), method)
call setDisEuclid(distance(1:npnt), point(1:ndim, 1:npnt), method)
! distance with respect to custom reference.
call setDisEuclid(distance, point(1:ndim), ref(1:ndim), method)
call setDisEuclid(distance(1:nref), point(1:ndim), ref(1:ndim, 1:nref), method)
call setDisEuclid(distance(1:npnt), point(1:ndim, 1:npnt), ref(1:ndim), method)
call setDisEuclid(distance(1:npnt, 1:nref), point(1:npnt), ref(1:nref), method)
call setDisEuclid(distance(1:npnt, 1:nref), point(1:ndim, 1:npnt), ref(1:ndim, 1:nref), method)
!
Generate and return the (squared) Euclidean distance of a (set of) point(s) in ndim-dimensions from a...
This module contains procedures and generic interfaces for computing the Euclidean norm of a single p...
Warning
The shapes of distance, point, and ref must be consistent as in the above interface at all times.
The condition size(point, 1) == size(ref, 1) .or. rank(distance) == 2 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.
BLAS/LAPACK equivalent:
The procedures under discussion combine, modernize, and extend the interface and functionalities of Version 3.11 of BLAS/LAPACK routine(s): dlapy3
See also
getDisEuclid
setDisEuclid
getDisMatEuclid
setDisMatEuclid
Intrinsic Fortran procedure hypot(x, y) (robust)
Intrinsic Fortran procedure norm2(x(:)) (unsafe)


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKH
4 use pm_kind, only: TKG => RKS ! all processor kinds are supported.
5 use pm_io, only: display_type
9 use pm_distUnif, only: getUnifRand
10
11 implicit none
12
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 of a point with respect to a reference in arbitrary dimensions without undue overflow/underflow.")
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%skip()
21
22 block
23 real(TKG) :: distance
24 real(TKG), allocatable :: point(:), ref(:)
25 call disp%skip()
26 call disp%show("point = getUnifRand(1, 10, 3_IK)")
27 point = getUnifRand(1, 10, 3_IK)
28 call disp%show("point")
29 call disp%show( point )
30 call disp%show("ref = getUnifRand(1, 10, 3_IK)")
31 ref = getUnifRand(1, 10, 3_IK)
32 call disp%show("ref")
33 call disp%show( ref )
34 call disp%show("call setDisEuclid(distance, point, ref, euclidsq)")
35 call setDisEuclid(distance, point, ref, euclidsq)
36 call disp%show("[norm2(point - ref), sqrt(distance), distance]")
37 call disp%show( [norm2(point - ref), sqrt(distance), distance] )
38 call disp%skip()
39 end block
40
41 call disp%skip()
42 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
43 call disp%show("! Compute the asymmetric matrix of (squared) distances of a set of points from a set of reference points with or without undue overflow/underflow.")
44 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
45 call disp%skip()
46
47 block
48 real(TKG), allocatable :: point(:,:), ref(:,:), distance(:,:)
49 integer(IK) :: ndim, npnt, nref
50 call disp%skip()
51 call disp%show("ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 4); nref = getUnifRand(1, 3)")
52 ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 4); nref = getUnifRand(1, 3)
53 call disp%show("[ndim, npnt, nref]")
54 call disp%show( [ndim, npnt, nref] )
55 call disp%show("point = getUnifRand(1, 10, ndim, npnt)")
56 point = getUnifRand(1, 10, ndim, npnt)
57 call disp%show("point")
58 call disp%show( point )
59 call disp%show("ref = getUnifRand(1, 10, ndim, nref)")
60 ref = getUnifRand(1, 10, ndim, nref)
61 call disp%show("ref")
62 call disp%show( ref )
63 call disp%show("call setResized(distance, [npnt, nref])")
64 call setResized(distance, [npnt, nref])
65 call disp%show("call setDisEuclid(distance, point, ref, euclidsq)")
66 call setDisEuclid(distance, point, ref, euclidsq)
67 call disp%show("distance")
68 call disp%show( distance )
69 call disp%show("call setDisEuclid(distance(1,:), point(:,1), ref, euclidsq)")
70 call setDisEuclid(distance(1,:), point(:,1), ref, euclidsq)
71 call disp%show("distance(1,:)")
72 call disp%show( distance(1,:) )
73 call disp%show("call setDisEuclid(distance(:,1), point, ref(:,1), euclidsq)")
74 call setDisEuclid(distance(:,1), point, ref(:,1), euclidsq)
75 call disp%show("distance(:,1)")
76 call disp%show( distance(:,1) )
77 call disp%show("call setDisEuclid(distance(1,1), point(:,1), ref(:,1), euclidsq)")
78 call setDisEuclid(distance(1,1), point(:,1), ref(:,1), euclidsq)
79 call disp%show("distance(1,1)")
80 call disp%show( distance(1,1) )
81 call disp%show("call setDisEuclid(distance(1,1), point(:,1), euclidsq) ! `ref` is the origin.")
82 call setDisEuclid(distance(1,1), point(:,1), euclidsq)
83 call disp%show("distance(1,1)")
84 call disp%show( distance(1,1) )
85 call disp%show("call setDisEuclid(distance(:,1), point, euclidsq) ! `ref` is the origin.")
86 call setDisEuclid(distance(:,1), point, euclidsq)
87 call disp%show("distance(:,1)")
88 call disp%show( distance(:,1) )
89 call disp%show("call setDisEuclid(distance, point(1,:), ref(1,:), euclid) ! 1D point and ref (faster than below).")
90 call setDisEuclid(distance, point(1,:), ref(1,:), euclid)
91 call disp%show("distance")
92 call disp%show( distance )
93 call disp%show("call setDisEuclid(distance, point(1:1,:), ref(1:1,:), euclid) ! For comparison with the above.")
94 call setDisEuclid(distance, point(1:1,:), ref(1:1,:), euclid)
95 call disp%show("distance")
96 call disp%show( distance )
97 call disp%show("call setDisEuclid(distance, point(1,:), ref(1,:), euclidsq) ! 1D point and ref (faster than below).")
98 call setDisEuclid(distance, point(1,:), ref(1,:), euclidsq)
99 call disp%show("distance")
100 call disp%show( distance )
101 call disp%show("call setDisEuclid(distance, point(1:1,:), ref(1:1,:), euclidsq) ! For comparison with the above.")
102 call setDisEuclid(distance, point(1:1,:), ref(1:1,:), euclidsq)
103 call disp%show("distance")
104 call disp%show( distance )
105 call disp%show("call setDisEuclid(distance, point(1,:), ref(1,:), euclidv) ! volume.")
106 call setDisEuclid(distance, point(1,:), ref(1,:), euclidv)
107 call disp%show("distance")
108 call disp%show( distance )
109 call disp%show("call setDisEuclid(distance, point(1,:), ref(1,:), euclid) ! reference volume.")
110 call setDisEuclid(distance, point(1,:), ref(1,:), euclid)
111 call disp%show("getVolUnitBall(1._TKG) * distance")
112 call disp%show( getVolUnitBall(1._TKG) * distance )
113 call disp%skip()
114 call disp%show("call setDisEuclid(distance, point, ref, euclidv) ! volume")
115 call setDisEuclid(distance, point, ref, euclidv) ! volume
116 call disp%show("distance")
117 call disp%show( distance )
118 call disp%show("call setDisEuclid(distance, point, ref, euclid) ! reference volume")
119 call setDisEuclid(distance, point, ref, euclid) ! reference volume
120 call disp%show("getVolUnitBall(real(ndim, TKG)) * distance**ndim")
121 call disp%show( getVolUnitBall(real(ndim, TKG)) * distance**ndim )
122 end block
123
124end 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...
Generate and return the natural logarithm of the volume of an -dimensional ball of unit-radius.
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(euclidv_type), parameter euclidv
This is a scalar parameter object of type euclidv_typethat is exclusively used to request computing E...
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 setting up and computing the properties of the hyper-...
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 of a point with respect to a reference in arbitrary dimensions without undue overflow/underflow.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7point = getUnifRand(1, 10, 3_IK)
8point
9+5.00000000, +10.0000000, +9.00000000
10ref = getUnifRand(1, 10, 3_IK)
11ref
12+5.00000000, +9.00000000, +3.00000000
13call setDisEuclid(distance, point, ref, euclidsq)
14[norm2(point - ref), sqrt(distance), distance]
15+6.08276224, +6.08276272, +37.0000000
16
17
18!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
19! Compute the asymmetric matrix of (squared) distances of a set of points from a set of reference points with or without undue overflow/underflow.
20!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
21
22
23ndim = getUnifRand(1, 3); npnt = getUnifRand(1, 4); nref = getUnifRand(1, 3)
24[ndim, npnt, nref]
25+3, +1, +2
26point = getUnifRand(1, 10, ndim, npnt)
27point
28+7.00000000
29+9.00000000
30+7.00000000
31ref = getUnifRand(1, 10, ndim, nref)
32ref
33+5.00000000, +8.00000000
34+6.00000000, +2.00000000
35+6.00000000, +7.00000000
36call setResized(distance, [npnt, nref])
37call setDisEuclid(distance, point, ref, euclidsq)
38distance
39+14.0000000, +50.0000000
40call setDisEuclid(distance(1,:), point(:,1), ref, euclidsq)
41distance(1,:)
42+14.0000000, +50.0000000
43call setDisEuclid(distance(:,1), point, ref(:,1), euclidsq)
44distance(:,1)
45+14.0000000
46call setDisEuclid(distance(1,1), point(:,1), ref(:,1), euclidsq)
47distance(1,1)
48+14.0000000
49call setDisEuclid(distance(1,1), point(:,1), euclidsq) ! `ref` is the origin.
50distance(1,1)
51+179.000000
52call setDisEuclid(distance(:,1), point, euclidsq) ! `ref` is the origin.
53distance(:,1)
54+179.000000
55call setDisEuclid(distance, point(1,:), ref(1,:), euclid) ! 1D point and ref (faster than below).
56distance
57+2.00000000, +1.00000000
58call setDisEuclid(distance, point(1:1,:), ref(1:1,:), euclid) ! For comparison with the above.
59distance
60+2.00000000, +1.00000000
61call setDisEuclid(distance, point(1,:), ref(1,:), euclidsq) ! 1D point and ref (faster than below).
62distance
63+4.00000000, +1.00000000
64call setDisEuclid(distance, point(1:1,:), ref(1:1,:), euclidsq) ! For comparison with the above.
65distance
66+4.00000000, +1.00000000
67call setDisEuclid(distance, point(1,:), ref(1,:), euclidv) ! volume.
68distance
69+4.00000000, +2.00000000
70call setDisEuclid(distance, point(1,:), ref(1,:), euclid) ! reference volume.
71getVolUnitBall(1._TKG) * distance
72+4.00000000, +2.00000000
73
74call setDisEuclid(distance, point, ref, euclidv) ! volume
75distance
76+219.422241, +1480.96106
77call setDisEuclid(distance, point, ref, euclid) ! reference volume
78getVolUnitBall(real(ndim, TKG)) * distance**ndim
79+219.422226, +1480.96106
80
Test:
test_pm_distanceEuclid
Internal naming convention:
The following illustrates the internal naming convention used for the procedures within this generic interface.
setDE_MEQ_D1_D1_D2_RK5()
|| ||| || || || |||
|| ||| || || || The type and kind parameters.
|| ||| || || The dimension of `ref` array: D1 => vector, D2 => matrix, XX => `ref` missing.
|| ||| || The dimension of `point` array: D1 => vector, D2 => matrix.
|| ||| The dimension of the output `distance`: D0 => scalar, D1 => vector, D2 => matrix.
|| The Method of Euclidean distance computation: MED => euclid_type(default/safe), MEU => eulidu_type(unsafe), MEV => eulidv_type(volume), MEQ => eulidsq_type(squared)
The abbreviation for DisEuclid to shorten procedure names.
This is a concrete derived type whose instances are exclusively used to request safe method of comput...
Todo:
Normal Priority: A benchmark comparison with the equivalent compiler implementations would be informative.


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 1065 of file pm_distanceEuclid.F90.


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