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

Generate and return a (collection) of random vector(s) of size ndim from the ndim-dimensional MultiVariate Uniform Ellipsoidal (MVUE) distribution, optionally with the specified input mean(1:ndim) and the specified subset of the Cholesky Factorization of the Gramian matrix of the MVUE distribution. More...

Detailed Description

Generate and return a (collection) of random vector(s) of size ndim from the ndim-dimensional MultiVariate Uniform Ellipsoidal (MVUE) distribution, optionally with the specified input mean(1:ndim) and the specified subset of the Cholesky Factorization of the Gramian matrix of the MVUE distribution.

The procedures of this generic interface are merely wrappers around the subroutine interface setUnifRand.

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.
(optional, default = rngf_type.)
[in]mean: The input contiguous vector of the same type, kind, rank, and size as the output rand, representing the mean of the MVUE distribution.
(optional, default = [(0., i = 1, size(rand))]. It must be present if the input argument chol is missing.)
[in]chol: The input contiguous matrix of shape (ndim, ndim) whose specified triangular subset contains the Cholesky Factorization of the Gramian matrix of the MVUE distribution.
(optional, the default is the Identity matrix of rank ndim. It must be present if and only if the input argument subset is 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]nsam: The input scalar integer of default kind IK containing the number of random MVUE vectors to generate.
(optional. If present, the output rand is of rank 2, otherwise is of rank 1.)
Returns
rand : The output vector of
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the MVUE distributed random output vector(s):
  1. If the input argument nsam is missing, then rand shall be of shape (1:ndim).
  2. If the input argument nsam is present, then rand shall be of shape (1:ndim, 1:nsam).


Possible calling interfaces

! single vector, using default rng
rand(1:ndim) = getUnifEllRand(mean(1:ndim))
rand(1:ndim) = getUnifEllRand(chol(1:ndim, 1:ndim), subset)
rand(1:ndim) = getUnifEllRand(mean(1:ndim), chol(1:ndim, 1:ndim), subset)
! single vector, using custom rng
rand(1:ndim) = getUnifEllRand(rng, mean(1:ndim))
rand(1:ndim) = getUnifEllRand(rng, chol(1:ndim, 1:ndim), subset)
rand(1:ndim) = getUnifEllRand(rng, mean(1:ndim), chol(1:ndim, 1:ndim), subset)
! collection of `nsam` vectors, using default rng
rand(1:ndim, 1:nsam) = getUnifEllRand(mean(1:ndim), nsam)
rand(1:ndim, 1:nsam) = getUnifEllRand(chol(1:ndim, 1:ndim), subset, nsam)
rand(1:ndim, 1:nsam) = getUnifEllRand(mean(1:ndim), chol(1:ndim, 1:ndim), subset, nsam)
! collection of `nsam` vectors, using custom rng
rand(1:ndim, 1:nsam) = getUnifEllRand(rng, mean(1:ndim), nsam)
rand(1:ndim, 1:nsam) = getUnifEllRand(rng, chol(1:ndim, 1:ndim), subset, nsam)
rand(1:ndim, 1:nsam) = getUnifEllRand(rng, mean(1:ndim), chol(1:ndim, 1:ndim), subset, nsam)
Generate and return a (collection) of random vector(s) of size ndim from the ndim-dimensional MultiVa...
This module contains classes and procedures for computing various statistical quantities related to t...
Remarks
The procedures under discussion are impure.
See also
pm_ellipsoid
getNormRand
setNormRand
getUnifRand
setUnifRand


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_kind, only: IK, LK, RK ! all real kinds are supported.
5 use pm_io, only: display_type
6 use pm_matrixChol, only: getMatChol, uppDia, lowDia
8
9 implicit none
10
11 integer(IK), parameter :: NDIM = 2_IK ! The number of dimensions of the MVN distribution.
12 real(RK) :: mean(NDIM), gramian(NDIM, NDIM), choLow(NDIM, NDIM), rand(NDIM)
13 integer(IK) :: info
14
15 type(display_type) :: disp
16 disp = display_type(file = "main.out.F90")
17
18 call disp%skip()
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%show("! Generate random numbers from the (Standard) Multivariate Uniform Ellipsoidal distribution.")
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
24 call disp%skip()
25
26 call disp%skip()
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show("! Multivariate Uniform Ellipsoidal random vector with a particular mean and Identity Gramian matrix.")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%skip()
31
32 call disp%skip()
33 call disp%show("mean = [-5., -5.]")
34 mean = [-5., -5.]
35 call disp%show("rand = getUnifEllRand(mean)")
36 rand = getUnifEllRand(mean)
37 call disp%show("rand")
38 call disp%show( rand )
39 call disp%skip()
40
41 call disp%skip()
42 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
43 call disp%show("! Multivariate Uniform Ellipsoidal random vector with zero mean and Gramian matrix specified via the Cholesky Lower Triangle.")
44 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
45 call disp%skip()
46
47 call disp%skip()
48 call disp%show("gramian = reshape([1., 1., 1., 4.], shape = [NDIM, NDIM])")
49 gramian = reshape([1., 1., 1., 4.], shape = [NDIM, NDIM])
50 call disp%show("choLow = getMatChol(gramian, lowDia)")
51 choLow = getMatChol(gramian, lowDia)
52 call disp%show("choLow")
53 call disp%show( choLow )
54 call disp%show("rand = getUnifEllRand(choLow, lowDia)")
55 rand = getUnifEllRand(choLow, lowDia)
56 call disp%show("rand")
57 call disp%show( rand )
58 call disp%skip()
59
60 call disp%skip()
61 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
62 call disp%show("! Multivariate Uniform Ellipsoidal random vector with given mean and Gramian matrix specified via the Cholesky Lower Triangle.")
63 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
64 call disp%skip()
65
66 call disp%skip()
67 call disp%show("mean")
68 call disp%show( mean )
69 call disp%show("choLow")
70 call disp%show( choLow )
71 call disp%show("rand = getUnifEllRand(mean, choLow, lowDia)")
72 rand = getUnifEllRand(mean, choLow, lowDia)
73 call disp%show("rand")
74 call disp%show( rand )
75 call disp%skip()
76
77 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
78 ! Output an example rand array for visualization.
79 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
80
81 block
82 integer(IK) :: fileUnit, i
83 integer(IK), parameter :: NVEC = 5000_IK
84 open(newunit = fileUnit, file = "getUnifEllRandMean.RK.txt")
85 do i = 1, NVEC
86 write(fileUnit,"(*(g0,:,','))") getUnifEllRand(mean)
87 end do
88 close(fileUnit)
89 open(newunit = fileUnit, file = "getUnifEllRandChol.RK.txt")
90 do i = 1, NVEC
91 write(fileUnit,"(*(g0,:,','))") getUnifEllRand(choLow, lowDia)
92 end do
93 close(fileUnit)
94 open(newunit = fileUnit, file = "getUnifEllRandMeanChol.RK.txt")
95 do i = 1, NVEC
96 write(fileUnit,"(*(g0,:,','))") getUnifEllRand(mean, choLow, lowDia)
97 end do
98 close(fileUnit)
99 end block
100
101end program example
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...
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 RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
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
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
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!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4! Generate random numbers from the (Standard) Multivariate Uniform Ellipsoidal distribution.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Multivariate Uniform Ellipsoidal random vector with a particular mean and Identity Gramian matrix.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14mean = [-5., -5.]
15rand = getUnifEllRand(mean)
16rand
17-0.78040564368514631, +0.25372573286899741
18
19
20!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
21! Multivariate Uniform Ellipsoidal random vector with zero mean and Gramian matrix specified via the Cholesky Lower Triangle.
22!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
23
24
25gramian = reshape([1., 1., 1., 4.], shape = [NDIM, NDIM])
26choLow = getMatChol(gramian, lowDia)
27choLow
28+1.0000000000000000, +0.0000000000000000
29+1.0000000000000000, +1.7320508075688772
30rand = getUnifEllRand(choLow, lowDia)
31rand
32-0.48470720199361156, +0.40355619489689837
33
34
35!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
36! Multivariate Uniform Ellipsoidal random vector with given mean and Gramian matrix specified via the Cholesky Lower Triangle.
37!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
38
39
40mean
41-5.0000000000000000, -5.0000000000000000
42choLow
43+1.0000000000000000, +0.0000000000000000
44+1.0000000000000000, +1.7320508075688772
45rand = getUnifEllRand(mean, choLow, lowDia)
46rand
47-5.6540525514129572, -5.8891770890215769
48
49

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
8
9linewidth = 2
10fontsize = 17
11
12for kind in ["RK"]:
13
14 pattern = "*." + kind + ".txt"
15 fileList = glob.glob(pattern)
16
17 for 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))
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 , s = 8
47 , zorder = 1000
48 )
49
50 ax_histx.hist(df.values[:, 0], bins = 50, zorder = 1000)
51 ax_histy.hist(df.values[:,1], bins = 50, orientation = "horizontal", zorder= 1000)
52
53 ax.set_xlabel("X", fontsize = 17)
54 ax.set_ylabel("Y", fontsize = 17)
55
56 plt.savefig(file.replace(".txt",".png"))

Visualization of the example output
Test:
test_pm_distUnifEll


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 559 of file pm_distUnifEll.F90.


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