File: stk_sampling_olhs.m

package info (click to toggle)
octave-stk 2.5.1-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 3,056 kB
  • sloc: ansic: 4,483; makefile: 32
file content (279 lines) | stat: -rw-r--r-- 8,139 bytes parent folder | download
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
% STK_SAMPLING_OLHS generates a random Orthogonal Latin Hypercube (OLH) sample
%
% CALL: X = stk_sampling_olhs (N)
%
%    generates a random Orthogonal Latin Hypercube (OLH) sample X, using the
%    construction of Ye (1998). The algorithm only works for sample sizes N
%    of the form 2^(R+1)+1, with R >= 1. Trying to generate an OLHS with a
%    value of N that is not of this form generates an error. The number of
%    factors is D = 2*R, and the OLHS is defined on [-1; 1]^D.
%
% CALL: X = stk_sampling_olhs (N, D)
%
%    does exactly the same thing, provided that there exists an integer R
%    such that N = 2^(R+1)+1 and D = 2*R (or D is empty).
%
% CALL: X = stk_sampling_olhs (N, D, BOX)
%
%    generates an OLHS on BOX. Again, D can be empty since the number of
%    factors can be deduced from N.
%
% CALL: X = stk_sampling_olhs (N, D, BOX, PERMUT)
%
%    uses a given permutation PERMUT, instead of a random permutation, to
%    initialize the construction of Ye (1998). As a result, the generated
%    OLHS is not random anymore. PERMUT must be a permutation of 1:2^R. If
%    BOX is empty, then the default domain [-1, 1]^D is used.
%
% NOTE: orthogonality
%
%    The samples generated by this functions are only orthogonal, stricty-
%    speaking, if BOX is a symmetric domain (e.g., [-1, 1] ^ D). Otherwise,
%    the generated samples should be called "uncorrelated".
%
% REFERENCE
%
%    Kenny Q. Ye, "Orthogonal Column Latin Hypercubes and Their
%    Application in Computer Experiments", Journal of the American
%    Statistical Association, 93(444), 1430-1439, 1998.
%    http://dx.doi.org/10.1080/01621459.1998.10473803
%
% See also: stk_sampling_randomlhs, stk_sampling_maximinlhs

% Copyright Notice
%
%    Copyright (C) 2017 CentraleSupelec
%    Copyright (C) 2012-2014 SUPELEC
%
%    Author:  Julien Bect  <julien.bect@centralesupelec.fr>

% Copying Permission Statement
%
%    This file is part of
%
%            STK: a Small (Matlab/Octave) Toolbox for Kriging
%               (http://sourceforge.net/projects/kriging)
%
%    STK is free software: you can redistribute it and/or modify it under
%    the terms of the GNU General Public License as published by the Free
%    Software Foundation,  either version 3  of the License, or  (at your
%    option) any later version.
%
%    STK is distributed  in the hope that it will  be useful, but WITHOUT
%    ANY WARRANTY;  without even the implied  warranty of MERCHANTABILITY
%    or FITNESS  FOR A  PARTICULAR PURPOSE.  See  the GNU  General Public
%    License for more details.
%
%    You should  have received a copy  of the GNU  General Public License
%    along with STK.  If not, see <http://www.gnu.org/licenses/>.

function [x, aux] = stk_sampling_olhs (n, d, box, permut, extended)

if nargin > 5,
    stk_error ('Too many input arguments.', 'TooManyInputArgs');
end

% Read argument dim
if (nargin < 2) || ((nargin < 3) && (isempty (d)))
    d = 1;  % Default dimension
elseif (nargin > 2) && (~ isempty (box))
    d = size (box, 2);
end

if nargin < 5, extended = false; end

%%% PROCESS INPUT ARGUMENTS %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

% recover the "order" r from the value of n
r = floor(log2(n - 1) - 1);
n_ = 2^(r + 1) + 1;
if (r <= 0) || (abs(n - n_) > eps),
    errmsg = 'n must be an integer of the form 2^(r+1) + 1 with r > 0';
    stk_error(errmsg, 'InvalidArgument');
end
n = n_;

if ~extended
    % check that d has the correct value for a "full" Ye98-OLHS
    % (other values of d can be reached by removing columns from a full OLHS)
    if (nargin > 1) && ~isempty(d),
        if d ~= 2 * r,
            errmsg = 'Incorrect value of d, please read the documentation...';
            stk_error(errmsg, 'InvalidArgument');
        end
    else
        d = 2 * r;
    end
else
    % check that d has the correct value for a "full" Cioppa-Lucas(2007) NOLHS
    % (other values of d can be reached by removing columns from a full NOLHS)
    if (nargin > 1) && ~isempty(d),
        if d ~= r + 1 + nchoosek(r, 2),
            errmsg = 'Incorrect value of d, please read the documentation...';
            stk_error(errmsg, 'InvalidArgument');
        end
    else
        d = r + 1 + nchoosek(r, 2);
    end
end

% number of "positive levels"
q = 2^r; % = (n - 1)/2

% read argument 'box'
if (nargin < 3) || isempty (box)
    box = stk_hrect (repmat ([-1; 1], 1, d));  % build a default box
else
    box = stk_hrect (box);  % convert input argument to a proper box
end

% permutation
if (nargin < 4) || isempty(permut),
    permut = randperm(q)';
else
    permut = permut(:);
    if ~isequal(sort(permut), (1:q)'),
        errmsg = sprintf('permut should be a permutation of 1:%d.', q);
        stk_error(errmsg, 'InvalidArgument');
    end
end

%%% CONSTRUCT AN "ORTHOGONAL" LHS FOLLOWING PROCESS INPUT ARGUMENTS %%%%%%%%%%%

% Construct permutation matrices A1, A2, ..., Ar
A = cell(1, r);
for i = 1:r,
    Ai = 1;
    for j = 1:i,
        Z  = zeros(size(Ai));
        Ai = [Z Ai; Ai Z]; %#ok<AGROW>
    end
    for j = (i+1):r,
        Z  = zeros(size(Ai));
        Ai = [Ai Z; Z Ai]; %#ok<AGROW>
    end
    A{i} = Ai;
end

% Construct the matrix M
M = [permut zeros(q, 2*r-1)];
for j = 1:r, % from column 2 to column r+1
    M(:, j+1) = A{j} * permut;
end
if ~extended, % OLHS / Ye (1998)
    for j = 1:(r-1), % from column r+2 to column 2*r
        M(:, j+r+1) = A{j} * A{r} * permut;
    end
else % NOLHS / Cioppa & Lucas(2007)
    col = r + 2;
    for j = 1:(r-1),
        for k = (j + 1):r,
            M(:, col) = A{j} * A{k} * permut;
            col = col + 1;
        end
    end
end

% Construct the matrix S
S = ones(q, 2*r);
for j = 1:r, % from column 2 to column r+1
    aj = 1;
    for l = r:(-1):1,
        if l == r - j + 1,
            aj = [-aj; aj]; %#ok<AGROW>
        else
            aj = [aj; aj]; %#ok<AGROW>
        end
    end
    S(:, j+1) = aj;
end
if ~extended, % OLHS / Ye (1998)
    for j = 1:(r-1), % from column r+2 to column 2*r
        S(:, r+1+j) = S(:, 2) .* S(:, j+2);
    end
else % NOLHS / Cioppa & Lucas(2007)
    col = r + 2;
    for j = 1:(r - 1),
        for k = (j + 1):r,
            S(:, col) = S(:, j+1) .* S(:, k+1);
            col = col + 1;
        end
    end
end

% Construct the matrix T
T = M .* S;

% Construct the OLHS X (with integer levels -q, ..., 0, ... +q)
x_integer_levels = [T; zeros(1, d); -T];

%%% CONVERT TO THE REQUESTED BOX %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

% Convert to positive integer levels (1 ... n)
x = x_integer_levels + q + 1;

% Convert to [0; 1]-valued levels
x = (2*x - 1) / (2*n);

% And, finally, convert to box
x = stk_dataframe (stk_rescale (x, [], box), box.colnames);

% Note: the results reported in Cioppa & Lucas correspond to the scaling
%   x = struct('a', stk_rescale(x, [min(x); max(x)], box));

%%% OUTPUT SOME AUXILIARY DATA IF REQUESTED %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

if nargout > 1,
    aux = struct('M', M, 'S', S, 'X', x_integer_levels);
end

end % function


%%
% Check error for incorrect number of input arguments

%!shared x, n, d, box, permut
%! n = 5; d = 2; box = [0 0; 1, 1]; permut = 1:2;

%!error x = stk_sampling_olhs();
%!test  x = stk_sampling_olhs(n);
%!test  x = stk_sampling_olhs(n, d);
%!test  x = stk_sampling_olhs(n, d, box);
%!test  x = stk_sampling_olhs(n, d, box, permut);
%!error x = stk_sampling_olhs(n, d, box, permut, pi);

%%
% Check that the output is a dataframe
% (all stk_sampling_* functions should behave similarly in this respect)

%!assert (isa (x, 'stk_dataframe'));

%%
% Check that column names are properly set, if available in box

%!assert (isequal (x.colnames, {}));

%!test
%! cn = {'W', 'H'};  box = stk_hrect (box, cn);
%! x = stk_sampling_olhs (n, d, box);
%! assert (isequal (x.colnames, cn));

%%
% Check output argument

%!test
%! for r = 1:5
%!
%!   n = 2 ^ (r + 1) + 1;  d = 2 * r;
%!   x = stk_sampling_olhs (n, d);
%!
%!   assert (isequal (size (x), [n d]));
%!
%!   box = repmat ([-1; 1], 1, d);
%!   assert (stk_is_lhs (x, n, d, box));
%!
%!   x = double (x);  w = x' * x;
%!   assert (stk_isequal_tolabs (w / w(1,1), eye (d)));
%!
%! end