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,
-
type rngf_type, implying the use of intrinsic Fortran uniform RNG.
-
type xoshiro256ssw_type, implying the use of xoshiro256** uniform RNG.
(optional, default = rngf_type.) |
[in] | mean | : The input contiguous vector of shape (1:ndim) , of the same type and kind 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:
-
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.
-
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
-
type
real
of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the MVUE distributed random output vector(s):
-
If the input argument
nsam
is missing, then rand
shall be of shape (1:ndim)
.
-
If the input argument
nsam
is present, then rand
shall be of shape (1:ndim, 1:nsam)
.
Possible calling interfaces ⛓
rand(
1:ndim)
= getUnifEllRand(mean(
1:ndim), chol(
1:ndim,
1:ndim), subset)
rand(
1:ndim)
= getUnifEllRand(rng, mean(
1:ndim), chol(
1:ndim,
1:ndim), subset)
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)
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...
- See also
- pm_ellipsoid
getNormRand
setNormRand
getUnifRand
setUnifRand
Example usage ⛓
11 integer(IK),
parameter :: NDIM
= 2_IK
12 real(RK) :: mean(NDIM), gramian(NDIM, NDIM), choLow(NDIM, NDIM), rand(NDIM)
15 type(display_type) :: disp
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
35 call disp%show(
"rand = getUnifEllRand(mean)")
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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)")
54 call disp%show(
"rand = getUnifEllRand(choLow, lowDia)")
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
71 call disp%show(
"rand = getUnifEllRand(mean, choLow, lowDia)")
82 integer(IK) :: fileUnit, i
83 integer(IK),
parameter :: NVEC
= 5000_IK
84 open(newunit
= fileUnit, file
= "getUnifEllRandMean.RK.txt")
89 open(newunit
= fileUnit, file
= "getUnifEllRandChol.RK.txt")
94 open(newunit
= fileUnit, file
= "getUnifEllRandMeanChol.RK.txt")
96 write(fileUnit,
"(*(g0,:,','))")
getUnifEllRand(mean, choLow, lowDia)
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
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...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
17-5.7764870484536175,
-5.5893097670771903
25gramian
= reshape([
1.,
1.,
1.,
4.], shape
= [NDIM, NDIM])
28+1.0000000000000000,
+0.0000000000000000
29+1.0000000000000000,
+1.7320508075688772
32+0.41392759853468242,
+1.8316830882893160
41-5.0000000000000000,
-5.0000000000000000
43+1.0000000000000000,
+0.0000000000000000
44+1.0000000000000000,
+1.7320508075688772
47-4.9870789622046052,
-3.9674727451642324
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
14 pattern =
"*." + kind +
".txt"
15 fileList = glob.glob(pattern)
19 df = pd.read_csv(file, delimiter =
",", header =
None)
22 left, width = 0.1, 0.65
23 bottom, height = 0.1, 0.65
27 fig = plt.figure(figsize = (8, 8))
29 plt.rcParams.update({
'font.size': fontsize - 2})
30 ax = fig.add_axes([left, bottom, width, height])
31 ax_histx = fig.add_axes([left, bottom + height + spacing, width, 0.2], sharex = ax)
32 ax_histy = fig.add_axes([left + width + spacing, bottom, 0.2, height], sharey = ax)
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")
40 ax_histy.tick_params(axis =
"y", labelleft =
False)
41 ax_histx.tick_params(axis =
"x", labelbottom =
False)
44 ax.scatter ( df.values[:, 0]
50 ax_histx.hist(df.values[:, 0], bins = 50, zorder = 1000)
51 ax_histy.hist(df.values[:, 1], bins = 50, orientation =
"horizontal", zorder = 1000)
53 ax.set_xlabel(
"X", fontsize = 17)
54 ax.set_ylabel(
"Y", fontsize = 17)
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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- 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.