-
Notifications
You must be signed in to change notification settings - Fork 0
/
mvrandn.m
295 lines (291 loc) · 10.3 KB
/
mvrandn.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
function rv=mvrandn(l,u,Sig,n)
%% truncated multivariate normal generator
% simulates 'n' random vectors exactly/perfectly distributed
% from the d-dimensional N(0,Sig) distribution (zero-mean normal
% with covariance 'Sig') conditional on l<X<u;
% infinite values for 'l' and 'u' are accepted;
% output: 'd' times 'n' array 'rv' storing random vectors;
%
% * Example:
% d=60;n=10^3;Sig=0.9*ones(d,d)+0.1*eye(d);l=(1:d)/d*4;u=l+2;
% X=mvrandn(l,u,Sig,n);boxplot(X','plotstyle','compact') % plot marginals
%
% * Notes: Algorithm may not work if 'Sig' is close to being rank deficient.
% Reference:
% Z. I. Botev (2015), "The Normal Law Under Linear Restrictions:
% Simulation and Estimation via Minimax Tilting", submitted to JRSS(B)
l=l(:); u=u(:); % set to column vectors
d=length(l); % basic input check
if (length(u)~=d)|(d~=sqrt(prod(size(Sig)))|any(l>u))
error('l, u, and Sig have to match in dimension with u>l')
end
% Cholesky decomposition of matrix with permuation
[Lfull,l,u,perm]=cholperm(Sig,l,u); % outputs the permutation
D=diag(Lfull);
if any(D<eps)
warning('Method may fail as covariance matrix is singular!')
end
L=Lfull./repmat(D,1,d);u=u./D; l=l./D; % rescale
L=L-eye(d); % remove diagonal
% find optimal tilting parameter non-linear equation solver
options=optimset('Diagnostics','off','Display','off',...
'Algorithm','trust-region-dogleg','Jacobian','on');
[soln,fval,exitflag] = fsolve(@(x)gradpsi(x,L,l,u),zeros(2*(d-1),1),options);
if exitflag~=1
warning('Method may fail as covariance matrix is close to singular!')
end
x=soln(1:(d-1));mu=soln(d:(2*d-2));
% compute psi star
psistar=psy(x,L,l,u,mu);
% start acceptance rejection sampling
rv=[]; accept=0; iter=0;
while accept<n % while # of accepted is less than n
[logpr,Z]=mvnrnd(n,L,l,u,mu); % simulate n proposals
idx=-log(rand(1,n))>(psistar-logpr); % acceptance tests
rv=[rv,Z(:,idx)]; % accumulate accepted
accept=size(rv,2); % keep track of # of accepted
iter=iter+1; % keep track of while loop iterations
if iter==10^3 % if iterations are getting large, give warning
warning('Acceptance prob. smaller than 0.001')
elseif iter>10^4 % if iterations too large, seek approximation only
accept=n;rv=[rv,Z]; % add the approximate samples
warning('Sample is only approximately distributed.')
end
end
% finish sampling; postprocessing
[dum,order]=sort(perm,'ascend');
rv=rv(:,1:n); % cut-down the array to desired n samples
rv=Lfull*rv; % reverse scaling of L
rv=rv(order,:); % reverse the Cholesky permutation
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function [p,Z]=mvnrnd(n,L,l,u,mu)
% generates the proposals from the exponentially tilted
% sequential importance sampling pdf;
% output: 'p', log-likelihood of sample
% Z, random sample
d=length(l); % Initialization
mu(d)=0;
Z=zeros(d,n); % create array for variables
p=0;
for k=1:d
% compute matrix multiplication L*Z
col=L(k,1:k)*Z(1:k,:);
% compute limits of truncation
tl=l(k)-mu(k)-col;
tu=u(k)-mu(k)-col;
%simulate N(mu,1) conditional on [tl,tu]
Z(k,:)=mu(k)+trandn(tl',tu');
% update likelihood ratio
p = p+lnNpr(tl,tu)+.5*mu(k)^2-mu(k)*Z(k,:);
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function [grad,J]=gradpsi(y,L,l,u)
% implements gradient of psi(x) to find optimal exponential twisting;
% assumes scaled 'L' with zero diagonal;
d=length(u);c=zeros(d,1);x=c;mu=c;
x(1:(d-1))=y(1:(d-1));mu(1:(d-1))=y(d:end);
% compute now ~l and ~u
c(2:d)=L(2:d,:)*x;lt=l-mu-c;ut=u-mu-c;
% compute gradients avoiding catastrophic cancellation
w=lnNpr(lt,ut);
pl=exp(-0.5*lt.^2-w)/sqrt(2*pi);
pu=exp(-0.5*ut.^2-w)/sqrt(2*pi);
P=pl-pu;
% output the gradient
dfdx=-mu(1:(d-1))+(P'*L(:,1:(d-1)))';
dfdm= mu-x+P;
grad=[dfdx;dfdm(1:(d-1))];
if nargout>1 % here compute Jacobian matrix
lt(isinf(lt))=0; ut(isinf(ut))=0;
dP=-P.^2+lt.*pl-ut.*pu; % dPdm
DL=repmat(dP,1,d).*L;
mx=-eye(d)+DL;
xx=L'*DL;
mx=mx(1:d-1,1:d-1);
xx=xx(1:d-1,1:d-1);
J=[xx,mx';
mx,diag(1+dP(1:d-1))];
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function p=psy(x,L,l,u,mu)
% implements psi(x,mu); assumes scaled 'L' without diagonal;
d=length(u);x(d)=0;mu(d)=0;x=x(:);mu=mu(:);
% compute now ~l and ~u
c=L*x;l=l-mu-c;u=u-mu-c;
p=sum(lnNpr(l,u)+.5*mu.^2-x.*mu);
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function p=lnNpr(a,b)
% computes ln(P(a<Z<b))
% where Z~N(0,1) very accurately for any 'a', 'b'
p=zeros(size(a));
% case b>a>0
I=a>0;
if any(I)
pa=lnPhi(a(I)); % log of upper tail
pb=lnPhi(b(I));
p(I)=pa+log1p(-exp(pb-pa));
end
% case a<b<0
idx=b<0;
if any(idx)
pa=lnPhi(-a(idx)); % log of lower tail
pb=lnPhi(-b(idx));
p(idx)=pb+log1p(-exp(pa-pb));
end
% case a<0<b
I=(~I)&(~idx);
if any(I)
pa=erfc(-a(I)/sqrt(2))/2; % lower tail
pb=erfc(b(I)/sqrt(2))/2; % upper tail
p(I)=log1p(-pa-pb);
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function p=lnPhi(x)
% computes logarithm of tail of Z~N(0,1) mitigating
% numerical roundoff errors;
p=-0.5*x.^2-log(2)+reallog(erfcx(x/sqrt(2)));
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function x=trandn(l,u)
%% truncated normal generator
% * efficient generator of a vector of length(l)=length(u)
% from the standard multivariate normal distribution,
% truncated over the region [l,u];
% infinite values for 'u' and 'l' are accepted;
% * Remark:
% If you wish to simulate a random variable
% 'Z' from the non-standard Gaussian N(m,s^2)
% conditional on l<Z<u, then first simulate
% X=trandn((l-m)/s,(u-m)/s) and set Z=m+s*X;
l=l(:);u=u(:); % make 'l' and 'u' column vectors
if length(l)~=length(u)
error('Truncation limits have to be vectors of the same length')
end
x=nan(size(l));
a=.66; % treshold for switching between methods
% threshold can be tuned for maximum speed for each Matlab version
% three cases to consider:
% case 1: a<l<u
I=l>a;
if any(I)
tl=l(I); tu=u(I); x(I)=ntail(tl,tu);
end
% case 2: l<u<-a
J=u<-a;
if any(J)
tl=-u(J); tu=-l(J); x(J)=-ntail(tl,tu);
end
% case 3: otherwise use inverse transform or accept-reject
I=~(I|J);
if any(I)
tl=l(I); tu=u(I); x(I)=tn(tl,tu);
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function x=trnd(l,u)
% uses acceptance rejection to simulate from truncated normal
x=randn(size(l)); % sample normal
% keep list of rejected
I=find(x<l|x>u); d=length(I);
while d>0 % while there are rejections
ly=l(I); % find the thresholds of rejected
uy=u(I);
y=randn(size(ly));
idx=y>ly&y<uy; % accepted
x(I(idx))=y(idx); % store the accepted
I=I(~idx); % remove accepted from list
d=length(I); % number of rejected
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function x=tn(l,u)
% samples a column vector of length=length(l)=length(u)
% from the standard multivariate normal distribution,
% truncated over the region [l,u], where -a<l<u<a for some
% 'a' and l and u are column vectors;
% uses acceptance rejection and inverse-transform method;
tol=2; % controls switch between methods
% threshold can be tuned for maximum speed for each platform
% case: abs(u-l)>tol, uses accept-reject from randn
I=abs(u-l)>tol; x=l;
if any(I)
tl=l(I); tu=u(I); x(I)=trnd(tl,tu);
end
% case: abs(u-l)<tol, uses inverse-transform
I=~I;
if any(I)
tl=l(I); tu=u(I); pl=erfc(tl/sqrt(2))/2; pu=erfc(tu/sqrt(2))/2;
x(I)=sqrt(2)*erfcinv(2*(pl-(pl-pu).*rand(size(tl))));
end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function x=ntail(l,u)
% samples a column vector of length=length(l)=length(u)
% from the standard multivariate normal distribution,
% truncated over the region [l,u], where l>0 and
% l and u are column vectors;
% uses acceptance-rejection from Rayleigh distr.
% similar to Marsaglia (1964);
c=l.^2/2; n=length(l); f=expm1(c-u.^2/2);
x=c-reallog(1+rand(n,1).*f); % sample using Rayleigh
% keep list of rejected
I=find(rand(n,1).^2.*x>c); d=length(I);
while d>0 % while there are rejections
cy=c(I); % find the thresholds of rejected
y=cy-reallog(1+rand(d,1).*f(I));
idx=rand(d,1).^2.*y<cy; % accepted
x(I(idx))=y(idx); % store the accepted
I=I(~idx); % remove accepted from list
d=length(I); % number of rejected
end
x=sqrt(2*x); % this Rayleigh transform can be delayed till the end
end
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
function [ L, l, u, perm ] = cholperm( Sig, l, u )
% Computes permuted lower Cholesky factor L for Sig
% by permuting integration limit vectors l and u.
% Outputs perm, such that Sig(perm,perm)=L*L'.
%
% Reference:
% Gibson G. J., Glasbey C. A., Elston D. A. (1994),
% "Monte Carlo evaluation of multivariate normal integrals and
% sensitivity to variate ordering",
% In: Advances in Numerical Methods and Applications, pages 120--126
d=length(l);perm=1:d; % keep track of permutation
L=zeros(d,d);z=zeros(d,1);
for j=1:d
pr=Inf(d,1); % compute marginal prob.
I=j:d; % search remaining dimensions
D=diag(Sig);
s=D(I)-sum(L(I,1:j-1).^2,2);
s(s<0)=eps;s=sqrt(s);
tl=(l(I)-L(I,1:j-1)*z(1:j-1))./s;
tu=(u(I)-L(I,1:j-1)*z(1:j-1))./s;
pr(I)=lnNpr(tl,tu);
% find smallest marginal dimension
[dummy,k]=min(pr);
% flip dimensions k-->j
jk=[j,k];kj=[k,j];
Sig(jk,:)=Sig(kj,:);Sig(:,jk)=Sig(:,kj); % update rows and cols of Sig
L(jk,:)=L(kj,:); % update only rows of L
l(jk)=l(kj);u(jk)=u(kj); % update integration limits
perm(jk)=perm(kj); % keep track of permutation
% construct L sequentially via Cholesky computation
s=Sig(j,j)-sum(L(j,1:j-1).^2);
if s<-0.01
error('Sigma is not positive semi-definite')
end
s(s<0)=eps;L(j,j)=sqrt(s);
L(j+1:d,j)=(Sig(j+1:d,j)-L(j+1:d,1:j-1)*(L(j,1:j-1))')/L(j,j);
% find mean value, z(j), of truncated normal:
tl=(l(j)-L(j,1:j-1)*z(1:j-1))/L(j,j);
tu=(u(j)-L(j,1:j-1)*z(1:j-1))/L(j,j);
w=lnNpr(tl,tu); % aids in computing expected value of trunc. normal
z(j)=(exp(-.5*tl.^2-w)-exp(-.5*tu.^2-w))/sqrt(2*pi);
end
end