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

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

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.

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.)
[out]rand: The output contiguous vector of shape (1:ndim) or matrix of shape (1:ndim, 1:nsam) of
  • type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the random output vector(s).
[in]mean: The input contiguous vector of shape (1:ndim), of the same type and kind as the output rand, representing the center of the distribution.
(optional, default = [(0., i = 1, size(rand))].)
[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.)


Possible calling interfaces

! single vector, using default rng
call setUnifEllRand(rand(1:ndim), mean(1:ndim))
call setUnifEllRand(rand(1:ndim), chol(1:ndim, 1:ndim), subset)
call setUnifEllRand(rand(1:ndim), mean(1:ndim), chol(1:ndim, 1:ndim), subset)
! single vector, using custom rng
call setUnifEllRand(rng, rand(1:ndim), mean(1:ndim))
call setUnifEllRand(rng, rand(1:ndim), chol(1:ndim, 1:ndim), subset)
call setUnifEllRand(rng, rand(1:ndim), mean(1:ndim), chol(1:ndim, 1:ndim), subset)
! collection of `nsam` vectors, using default rng
call setUnifEllRand(rand(1:ndim, 1:nsam), mean(1:ndim))
call setUnifEllRand(rand(1:ndim, 1:nsam), chol(1:ndim, 1:ndim), subset)
call setUnifEllRand(rand(1:ndim, 1:nsam), mean(1:ndim), chol(1:ndim, 1:ndim), subset)
! collection of `nsam` vectors, using custom rng
call setUnifEllRand(rng, rand(1:ndim, 1:nsam), mean(1:ndim))
call setUnifEllRand(rng, rand(1:ndim, 1:nsam), chol(1:ndim, 1:ndim), subset)
call setUnifEllRand(rng, rand(1:ndim, 1:nsam), mean(1:ndim), chol(1:ndim, 1:ndim), subset)
Return a (collection) of random vector(s) of size ndim from the ndim-dimensional MultiVariate Uniform...
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 all(shape(chol) == size(rand, 1)) 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.
The procedures under discussion are impure.
See also
getNormRand
setNormRand
getMultiNormLogPDF


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 from a Standard distribution.")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%skip()
31
32 call disp%skip()
33 call disp%show("call setUnifEllRand(rand)")
34 call setUnifEllRand(rand)
35 call disp%show("rand")
36 call disp%show( rand )
37 call disp%skip()
38
39 call disp%skip()
40 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
41 call disp%show("! Multivariate Uniform Ellipsoidal random vector with a particular mean and Identity Gramian matrix.")
42 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
43 call disp%skip()
44
45 call disp%skip()
46 call disp%show("mean = [-5., -5.]")
47 mean = [-5., -5.]
48 call disp%show("call setUnifEllRand(rand, mean)")
49 call setUnifEllRand(rand, mean)
50 call disp%show("rand")
51 call disp%show( rand )
52 call disp%skip()
53
54 call disp%skip()
55 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
56 call disp%show("! Multivariate Uniform Ellipsoidal random vector with zero mean and Gramian matrix specified via the Cholesky Lower Triangle.")
57 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
58 call disp%skip()
59
60 call disp%skip()
61 call disp%show("gramian = reshape([1., 1., 1., 4.], shape = [NDIM, NDIM])")
62 gramian = reshape([1., 1., 1., 4.], shape = [NDIM, NDIM])
63 call disp%show("choLow = getMatChol(gramian, lowDia)")
64 choLow = getMatChol(gramian, lowDia)
65 call disp%show("choLow")
66 call disp%show( choLow )
67 call disp%show("call setUnifEllRand(rand, choLow, lowDia)")
68 call setUnifEllRand(rand, choLow, lowDia)
69 call disp%show("rand")
70 call disp%show( rand )
71 call disp%skip()
72
73 call disp%skip()
74 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
75 call disp%show("! Multivariate Uniform Ellipsoidal random vector with given mean and Gramian matrix specified via the Cholesky Lower Triangle.")
76 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
77 call disp%skip()
78
79 call disp%skip()
80 call disp%show("mean")
81 call disp%show( mean )
82 call disp%show("choLow")
83 call disp%show( choLow )
84 call disp%show("call setUnifEllRand(rand, mean, choLow, lowDia)")
85 call setUnifEllRand(rand, mean, choLow, lowDia)
86 call disp%show("rand")
87 call disp%show( rand )
88 call disp%skip()
89
90 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
91 ! Output an example rand array for visualization.
92 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
93
94 block
95 integer(IK) :: fileUnit, i
96 integer(IK), parameter :: NVEC = 5000_IK
97 open(newunit = fileUnit, file = "setUnifEllRand.RK.txt")
98 do i = 1, NVEC
99 call setUnifEllRand(rand)
100 write(fileUnit,"(*(g0,:,','))") rand
101 end do
102 close(fileUnit)
103 open(newunit = fileUnit, file = "setUnifEllRandMean.RK.txt")
104 do i = 1, NVEC
105 call setUnifEllRand(rand, mean)
106 write(fileUnit,"(*(g0,:,','))") rand
107 end do
108 close(fileUnit)
109 open(newunit = fileUnit, file = "setUnifEllRandChol.RK.txt")
110 do i = 1, NVEC
111 call setUnifEllRand(rand, choLow, lowDia)
112 write(fileUnit,"(*(g0,:,','))") rand
113 end do
114 close(fileUnit)
115 open(newunit = fileUnit, file = "setUnifEllRandMeanChol.RK.txt")
116 do i = 1, NVEC
117 call setUnifEllRand(rand, mean, choLow, lowDia)
118 write(fileUnit,"(*(g0,:,','))") rand
119 end do
120 close(fileUnit)
121 end block
122
123end 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 from a Standard distribution.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14call setUnifEllRand(rand)
15rand
16-0.66924052351860452, -0.39463295180884450
17
18
19!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
20! Multivariate Uniform Ellipsoidal random vector with a particular mean and Identity Gramian matrix.
21!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
22
23
24mean = [-5., -5.]
25call setUnifEllRand(rand, mean)
26rand
27-5.0544090255341958, -4.0071909973918673
28
29
30!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
31! Multivariate Uniform Ellipsoidal random vector with zero mean and Gramian matrix specified via the Cholesky Lower Triangle.
32!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
33
34
35gramian = reshape([1., 1., 1., 4.], shape = [NDIM, NDIM])
36choLow = getMatChol(gramian, lowDia)
37choLow
38+1.0000000000000000, +0.0000000000000000
39+1.0000000000000000, +1.7320508075688772
40call setUnifEllRand(rand, choLow, lowDia)
41rand
42-0.28254481136513249, -1.6790197729545753
43
44
45!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
46! Multivariate Uniform Ellipsoidal random vector with given mean and Gramian matrix specified via the Cholesky Lower Triangle.
47!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
48
49
50mean
51-5.0000000000000000, -5.0000000000000000
52choLow
53+1.0000000000000000, +0.0000000000000000
54+1.0000000000000000, +1.7320508075688772
55call setUnifEllRand(rand, mean, choLow, lowDia)
56rand
57-5.8651233976746084, -6.7217395072616135
58
59

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


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