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

Return a collection of random vectors of size ndim from the ndim-dimensional Multiple MultiVariate Uniform Ellipsoidal (MMVUE) distribution, with the specified input mean(1:ndim, 1:nell) and optionally the specified subset of the Cholesky Factorization of the Gramian matrices of the MMVUE distribution.
More...

Detailed Description

Return a collection of random vectors of size ndim from the ndim-dimensional Multiple MultiVariate Uniform Ellipsoidal (MMVUE) distribution, with the specified input mean(1:ndim, 1:nell) and optionally the specified subset of the Cholesky Factorization of the Gramian matrices of the MMVUE distribution.

Here,

  1. the variable ndim represents the number of dimensions of the domain of the distribution (that is, the number of dimensions of the ellipsoids),
  2. the variable nsam represents the number of randomly sampled points from the set of ellipsoids specified by the input arguments,
  3. the variable nell represents the number of ellipsoids in the distirbution.
Parameters
[in,out]rng: The input/output scalar that can be an object of,
  1. type rngf_type, implying the use of intrinsic Fortran uniform RNG.
  2. type xoshiro256ssw_type, implying the use of xoshiro256** uniform RNG.
[out]rand: The output contiguous vector of shape (1:ndim) (or matrix of shape (1:ndim, 1:nsam)) of
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the (nsam) random output vector(s).
[out]mahalSq: The output vector of shape (1:nell) (or matrix of shape (1:nell, 1:nsam)) of the same type and kind as the output rand, containing the squared Mahalanobis distance(s) of individual randomly sampled vector(s rand(:, isam)) from the centers of the specified ellipsoids in the distribution via the input argument mean(1:ndim, 1:nell).
  1. A value larger than 1 for the vector element mahalSq(iell) (or matrix element mahalSq(iell, isam)) implies that the random vector rand(1:ndim, isam) is outside ellipsoid iell specified by the input arguments.
  2. A value less than 1 implies that the random vector is within the ellipsoid.
Thus,
  1. for vector mahalSq(1:nell), the expression count(mahalSq(1:nell) <= 1) yields the membership count of rand(1:ndim) in all input ellipsoids.
  2. for matrix mahalSq(1:nell, 1:nsam), the expression count(mahalSq(1:nell, isam) <= 1) yields the membership count of rand(1:ndim, isam) in all input ellipsoids.
[out]invmul: The output scalar (or vector of shape (1:nsam)) of the same type and kind as the output rand, (each element of) which represents the inverse of the number of ellipsoids within which the corresponding random vector rand(1:ndim) (or rand(1:ndim, isam)) falls.
[out]membership: The output scalar (or vector of shape (1:nsam)) of type integer of default kind IK, (the isam element of) which contains the ID of the ellipsoid from which (the isam) sampled point has been generated.
[in]mean: The input contiguous matrix of shape (1:ndim, 1:nell), of the same type and kind as the output rand, representing the centers of the ellipsoids in the distribution.
[in]chol: The input contiguous array of shape (1:ndim, 1:ndim, 1:nell) whose specified triangular subset contains the Cholesky Factorization of the Gramian matrix of the MMVUE distribution.
(optional, the default is the Identity matrix of rank ndim. It must be present if and only if the input argument subset and invGram are also present.)
[in]subset: The input scalar constant that can be any of the following:
  1. the constant uppDia or an object of type uppDia_type implying that the upper-diagonal triangular block of the input chol must be used while the lower subset is not referenced.
  2. the constant lowDia or an object of type lowDia_type implying that the lower-diagonal triangular block of the input chol must be used while the upper subset is not referenced.
This argument is merely a convenience to differentiate the different procedure functionalities within this generic interface.
(optional. It must be present if and only if the input argument chol is present.)
[in]invGram: The input array of shape (1:ndim, 1:ndim, 1:nell) of the same type and kind as the output argument rand, containing the collection of square matrices of full inverse representative Gramian matrix of the \(\ndim\)-dimensional ellipsoids in the distribution.
This argument is needed to determine the membership of the output random vector(s) in the ellipsoids.
(optional. It must be present if and only if the input argument chol is present.)
[in]cumPropVol: The input vector of shape (0:nell) of the same type and kind as the output argument rand, the subset (1:nell) contains the cumulative proportions of the multiplicative traces of the input Cholesky factors chol of the ellipsoids in the distribution.
This argument can be computed as,
use pm_mathCumPropExp, only: setCumPropExp, sequence
real(typeof(rand)) :: cumPropVol(size(mean, 2))
do iell = 1, size(mean, 2, IK)
cumPropVol(iell) = getMatMulTraceLog(chol(:, :, iell))
end do
call setCumPropExp(cumPropVol, maxArray = maxval(cumPropVol), control = sequence)
Return the cumulative sum of the proportions of the exponential of the input array,...
Generate and return the natural logarithm of the multiplicative trace of an input square matrix of ty...
This module contains the procedures and interfaces for computing the cumulative sum of the exponentia...
This module contains procedures and generic interfaces for computing the additive or multiplicative t...
The presence of this argument significantly boosts the algorithm performance by avoiding redundant, repetitive computations.
(optional. It must be present if and only if the input arguments chol and invGram are present and the output argument rand is of rank 1.)


Possible calling interfaces

use pm_distUnifElls, only: setUnifEllsRand, uppDia, lowDia
! 1D output random vector.
call setUnifEllsRand(rng, rand(1:ndim), mahalSq(1:nell), invmul, membership, mean(1:ndim, 1:nell))
call setUnifEllsRand(rng, rand(1:ndim), mahalSq(1:nell), invmul, membership, mean(1:ndim, 1:nell), chol(1:ndim, 1:ndim, 1:nell), subset, invGram(1:ndim, 1:ndim, 1:nell), cumPropVol(1:nell))
! 2D output random vector.
call setUnifEllsRand(rng, rand(1:ndim, 1:nsam), mahalSq(1:nell, 1:nsam), invmul(1:nsam), membership(1:nsam), mean(1:ndim, 1:nell))
call setUnifEllsRand(rng, rand(1:ndim, 1:nsam), mahalSq(1:nell, 1:nsam), invmul(1:nsam), membership(1:nsam), mean(1:ndim, 1:nell), chol(1:ndim, 1:ndim, 1:nell), subset, invGram(1:ndim, 1:ndim, 1:nell))
Return a collection of random vectors of size ndim from the ndim-dimensional Multiple MultiVariate Un...
This module contains classes and procedures for computing various statistical quantities related to t...
Warning
The condition size(mean, 1) == size(rand, 1) must hold for the corresponding input arguments.
The condition size(mean, 2) == size(cumPropVol) must hold for the corresponding input arguments.
The condition rank(rand) == 2 .and. size(invmul) == size(rand, 2) must hold for the corresponding input arguments.
The condition rank(rand) == 2 .and. size(invmul) == size(membership) must hold for the corresponding input arguments.
The condition rank(rand) == 2 .and. all(shape(mahalSq) == [size(mean, 2), size(rand, 2)]) must hold for the corresponding input arguments.
The condition size(chol, 1) <= size(chol, 2) must hold for the corresponding input arguments (to ensure that arguments can be passed without data copy).
The condition all(shape(invGram) == [size(mean, 1), size(mean, 1), size(mean, 2)) must hold for the corresponding input arguments.
The condition all([size(chol, 1), size(chol, 3)] == shape(mean)) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure. The procedures of this generic interface are pure when the input argument rng is set to xoshiro256ssw_type and the compile-time macro CHECK_ENABLED is set to 0 or is undefined.
See also
pm_distNorm
pm_distUnif
pm_distUnifEll
pm_distUnifPar
pm_distUnifElls
pm_distMultiNorm
getUnifEllsLogPDF


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_kind, only: IK, LK
5 use pm_io, only: display_type
8 use pm_matrixChol, only: getMatChol, uppDia, lowDia
10 use pm_distChol, only: getCholRand, uppDia
11 use pm_matrixInv, only: getMatInv, choUpp
12 use pm_io, only: getErrTableWrite, trans
13
14 implicit none
15
16 integer(IK) :: ndim, iell, nell, nsam
17
18 type(display_type) :: disp
19 disp = display_type(file = "main.out.F90")
20
21 block
22 use pm_kind, only: TKG => RKD
23 real(TKG), allocatable :: mean(:,:), invGram(:, :, :), chol(:, :, :), rand(:, :), mahalSq(:, :), invmul(:)
24 integer(IK), allocatable :: membership(:)
25
26 call disp%skip()
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show("! Random vectors from Multiple Multivariate Uniform Balls with particular means.")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%skip()
31
32 call disp%skip()
33 call disp%show("ndim = 2; nell = getUnifRand(5, 20); nsam = nell * 1000")
34 ndim = 2; nell = getUnifRand(5, 20); nsam = nell * 1000
35 call disp%show("[ndim, nell, nsam]")
36 call disp%show( [ndim, nell, nsam] )
37 call disp%show("mean = getUnifRand(-5._TKG, +5._TKG, ndim, nell)")
38 mean = getUnifRand(-5._TKG, +5._TKG, ndim, nell)
39 call disp%show("mean")
40 call disp%show( mean )
41 call disp%show("call setResized(mahalSq, [nell, nsam])")
42 call setResized(mahalSq, [nell, nsam])
43 call disp%show("call setResized(rand, [ndim, nsam])")
44 call setResized(rand, [ndim, nsam])
45 call disp%show("call setResized(membership, nsam)")
46 call setResized(membership, nsam)
47 call disp%show("call setResized(invmul, nsam)")
48 call setResized(invmul, nsam)
49 call disp%show("call setUnifEllsRand(rngf_type(), rand, mahalSq, invmul, membership, mean)")
50 call setUnifEllsRand(rngf_type(), rand, mahalSq, invmul, membership, mean)
51 call disp%show("if (0 /= getErrTableWrite(SK_'setUnifEllsRandMean.txt', reshape([transpose(rand), real(membership, TKG)], [nsam, ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'")
52 if (0 /= getErrTableWrite(SK_'setUnifEllsRandMean.txt', reshape([transpose(rand), real(membership, TKG)], [nsam, ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
53 call disp%skip()
54
55 call disp%skip()
56 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
57 call disp%show("! Random vectors from Multiple Multivariate Uniform Ellipsoids with particular means.")
58 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%skip()
60
61 call disp%show("call setResized(chol, [ndim, ndim, nell])")
62 call setResized(chol, [ndim, ndim, nell])
63 call disp%show("call setResized(invGram, [ndim, ndim, nell])")
64 call setResized(invGram, [ndim, ndim, nell])
65 call disp%show("do iell = 1, nell; chol(:, :, iell) = getCholRand(0._TKG, ndim, uppDia); invGram(:, :, iell) = getMatInv(chol(:, :, iell), choUpp); end do")
66 do iell = 1, nell; chol(:, :, iell) = getCholRand(0._TKG, ndim, uppDia); invGram(:, :, iell) = getMatInv(chol(:, :, iell), choUpp); end do
67 call disp%show("chol")
68 call disp%show( chol )
69 call disp%show("call setUnifEllsRand(rngf_type(), rand, mahalSq, invmul, membership, mean, chol, uppDia, invGram)")
70 call setUnifEllsRand(rngf_type(), rand, mahalSq, invmul, membership, mean, chol, uppDia, invGram)
71 call disp%show("if (0 /= getErrTableWrite(SK_'setUnifEllsRandMeanChol.txt', reshape([transpose(rand), real(membership, TKG)], [nsam, ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'")
72 if (0 /= getErrTableWrite(SK_'setUnifEllsRandMeanChol.txt', reshape([transpose(rand), real(membership, TKG)], [nsam, ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
73 call disp%skip()
74
75 end block
76
77end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return a random upper and lower Cholesky factorization.
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
Definition: pm_io.F90:5940
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
Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-def...
Generate and return the full inverse of an input matrix of general or triangular form directly or thr...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains classes and procedures for generating random upper or lower Cholesky factor tria...
Definition: pm_distChol.F90:43
This module contains classes and procedures for computing various statistical quantities related to t...
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 RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
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
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
This module contains abstract and concrete derived types and procedures related to the inversion of s...
This is a concrete derived type whose instances can be used to define/request the default uniform ran...
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! Random vectors from Multiple Multivariate Uniform Balls with particular means.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = 2; nell = getUnifRand(5, 20); nsam = nell * 1000
8[ndim, nell, nsam]
9+2, +8, +8000
10mean = getUnifRand(-5._TKG, +5._TKG, ndim, nell)
11mean
12-2.9804421326474309, -1.3482927006721832, -1.8110630665431691, -1.2959631379537750, -2.3023795164443817, -3.2850247288326160, -3.7832705104089750, +2.5396567470788884
13-1.5727114409381304, +2.1622703623306281, -1.4907322130678793, -0.43266049777419990, +0.46715356928837792, -1.5122260044751368, +3.7114548256281150, -3.1941187613007731
14call setResized(mahalSq, [nell, nsam])
15call setResized(rand, [ndim, nsam])
16call setResized(membership, nsam)
17call setResized(invmul, nsam)
18call setUnifEllsRand(rngf_type(), rand, mahalSq, invmul, membership, mean)
19if (0 /= getErrTableWrite(SK_'setUnifEllsRandMean.txt', reshape([transpose(rand), real(membership, TKG)], [nsam, ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
20
21
22!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
23! Random vectors from Multiple Multivariate Uniform Ellipsoids with particular means.
24!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
25
26call setResized(chol, [ndim, ndim, nell])
27call setResized(invGram, [ndim, ndim, nell])
28do iell = 1, nell; chol(:, :, iell) = getCholRand(0._TKG, ndim, uppDia); invGram(:, :, iell) = getMatInv(chol(:, :, iell), choUpp); end do
29chol
30slice(:,:,1) =
31+1.0000000000000000, -0.95884781933593966
32+0.0000000000000000, +0.28392051591019835
33slice(:,:,2) =
34+1.0000000000000000, -0.95291922971898457
35+0.0000000000000000, +0.30322424314651525
36slice(:,:,3) =
37+1.0000000000000000, +0.45418607629981694
38+0.0000000000000000, +0.89090684591340796
39slice(:,:,4) =
40+1.0000000000000000, -0.95610742758034245E-1
41+0.0000000000000000, +0.99541879923440113
42slice(:,:,5) =
43+1.0000000000000000, +0.97555519530861168
44+0.0000000000000000, +0.21975454695267746
45slice(:,:,6) =
46+1.0000000000000000, +0.23659867030462828
47+0.0000000000000000, +0.97160746662944175
48slice(:,:,7) =
49+1.0000000000000000, -0.59234887991871332
50+0.0000000000000000, +0.80568157758449821
51slice(:,:,8) =
52+0.99999999999999989, -0.34683977841435364E-1
53+0.0000000000000000, +0.99939832983705490
54call setUnifEllsRand(rngf_type(), rand, mahalSq, invmul, membership, mean, chol, uppDia, invGram)
55if (0 /= getErrTableWrite(SK_'setUnifEllsRandMeanChol.txt', reshape([transpose(rand), real(membership, TKG)], [nsam, ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
56
57

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6import glob
7import sys
8import os
9procname = os.path.basename(os.getcwd())
10
11linewidth = 2
12fontsize = 17
13
14pattern = procname + "*.txt"
15fileList = glob.glob(pattern)
16
17for file in fileList:
18
19 df = pd.read_csv(file, delimiter = ",", header = None)
20
21 # definitions for the axes
22 left, width = 0.1, 0.65
23 bottom, height = 0.1, 0.65
24 spacing = 0.015
25
26 # start with a square Figure
27 fig = plt.figure(figsize = (8, 8), dpi = 200)
28
29 plt.rcParams.update({'font.size': fontsize - 2})
30 ax = fig.add_axes([left, bottom, width, height]) # scatter plot
31 ax_histx = fig.add_axes([left, bottom + height + spacing, width, 0.2], sharex = ax) # histx
32 ax_histy = fig.add_axes([left + width + spacing, bottom, 0.2, height], sharey = ax) # histy
33
34 for axes in [ax, ax_histx, ax_histy]:
35 axes.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
36 axes.tick_params(axis = "y", which = "minor")
37 axes.tick_params(axis = "x", which = "minor")
38
39 # no labels
40 ax_histy.tick_params(axis = "y", labelleft = False)
41 ax_histx.tick_params(axis = "x", labelbottom = False)
42
43 # the scatter plot:
44 ax.scatter ( df.values[:, 0]
45 , df.values[:, 1]
46 , c = df.values[:, -1]
47 , zorder = 1000
48 , marker = "o"
49 , alpha = 1
50 , s = 0.3
51 )
52
53 ax_histx.hist(df.values[:, 0], bins = 50, zorder = 1000)
54 ax_histy.hist(df.values[:, 1], bins = 50, orientation = "horizontal", zorder = 1000)
55
56 ax.set_xlabel("X", fontsize = 17)
57 ax.set_ylabel("Y", fontsize = 17)
58
59 plt.savefig(file.replace(".txt",".svg"))

Visualization of the example output
Test:
test_pm_distUnifElls


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, April 23, 2017, 12:36 AM, Institute for Computational Engineering and Sciences (ICES), University of Texas at Austin

Definition at line 915 of file pm_distUnifElls.F90.


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