ParaMonte MATLAB 3.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
getRand.m
Go to the documentation of this file.
1%> \brief
2%> Return a (set of) multivariate Normal random vector(s).<br>
3%>
4%> \details
5%> This RNG function can be potentially faster than the intrinsic MATLAB equivalent
6%> because it uses the Cholesky factorization of the distribution covariance matrix.<br>
7%>
8%> \param[in] mean : The input vector of MATLAB ``real``,
9%> representing the mean of a Multivariate Normal
10%> distribution in ``size(mean)`` dimensional space.<br>
11%> (**optional**. default = ``[]``. It must be present if ``cholow`` is missing.)
12%> \param[in] cholow : The input square matrix of MATLAB ``real``,
13%> representing the lower-triangle of the Cholesky
14%> factorization of the covariance matrix of the target
15%> Multivariate Normal distribution in ``numel(mean)`` dimensional space.<br>
16%> This argument can be obtained by passing the covariance matrix ``covmat``
17%> of the distribution to the MATLAB intrinsic function ``chol(covmat, "lower")``.<br>
18%> (**optional**. default = ``[]``. It must be present if ``mean`` is missing.)
19%> \param[in] s1 : The input vector of MATLAB ``real``,
20%> representing the mean of a Multivariate Normal
21%> distribution in ``size(mean)`` dimensional space.<br>
22%> (**optional**. default = ``1``)
23%>
24%> \return
25%> ``rand`` : The output vector of MATLAB ``real`` of
26%> shape ``(numel(mean), 1)`` containing a random vector
27%> from the specified Multivariate Normal distribution.<br>
28%>
29%> \interface{getRand}
30%> \code{.m}
31%>
32%> rand = pm.stats.dist.mvn.getRand(mean)
33%> rand = pm.stats.dist.mvn.getRand([], cholow)
34%> rand = pm.stats.dist.mvn.getRand(mean, cholow)
35%> rand = pm.stats.dist.mvn.getRand([], cholow, s1)
36%> rand = pm.stats.dist.mvn.getRand(mean, cholow, s1)
37%>
38%> \endcode
39%>
40%> \example{getRand}
41%> \include{lineno} example/stats/dist/mvn/getRand/main.m
42%> \vis{getRand}
43%> \image html example/stats/dist/mvn/getRand/mvn.getRand.hist1.png width=700
44%> \image html example/stats/dist/mvn/getRand/mvn.getRand.hist2.png width=700
45%> \image html example/stats/dist/mvn/getRand/mvn.getRand.scatter1.png width=700
46%>
47%> \final{getRand}
48%>
49%> \author
50%> \JoshuaOsborne, May 21 2024, 12:06 AM, University of Texas at Arlington<br>
51%> \FatemehBagheri, May 20 2024, 1:25 PM, NASA Goddard Space Flight Center (GSFC), Washington, D.C.<br>
52%> \AmirShahmoradi, May 16 2016, 9:03 AM, Oden Institute for Computational Engineering and Sciences (ICES), UT Austin<br>
53function rand = getRand(mean, cholow, s1)
54 if nargin < 3
55 s1 = 1;
56 end
57 if 1 < nargin
58 if ~isempty(cholow)
59 if 1 < s1
60 rand = randn(size(cholow, 1), s1);
61 for irand = 1 : s1
62 rand(:, irand) = cholow * rand(:, irand);
63 end
64 else
65 rand = cholow * randn(size(cholow, 1), s1);
66 end
67 else
68 rand = randn(numel(mean), s1);
69 end
70 else
71 rand = randn(numel(mean), s1);
72 end
73 if ~isempty(mean)
74 rand = bsxfun(@plus, mean(:), rand);
75 end
76end
function getRand(in ndim, in scale)
Generate and return a random positive-definite (correlation or covariance) matrix using the Gram meth...