-
Notifications
You must be signed in to change notification settings - Fork 240
/
multifi_cokriging.py
1047 lines (888 loc) · 39.8 KB
/
multifi_cokriging.py
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
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
Integrates the Multi-Fidelity Co-Kriging method described in [LeGratiet2013].
(Author: Remi Vauclin vauclin.remi@gmail.com)
This code was implemented using the package scikit-learn as basis.
(Author: Vincent Dubourg, vincent.dubourg@gmail.com)
OpenMDAO adaptation. Regression and correlation functions were directly copied
from scikit-learn package here to avoid scikit-learn dependency.
(Author: Remi Lafage, remi.lafage@onera.fr)
ISAE/DMSM - ONERA/DCPS
"""
import numpy as np
from numpy import atleast_2d as array2d
from scipy import linalg
from scipy.optimize import minimize
from scipy.spatial.distance import squareform
from openmdao.surrogate_models.surrogate_model import MultiFiSurrogateModel
import logging
_logger = logging.getLogger()
MACHINE_EPSILON = np.finfo(np.double).eps # machine precision
NUGGET = 10. * MACHINE_EPSILON # nugget for robustness
INITIAL_RANGE_DEFAULT = 0.3 # initial range for optimizer
TOLERANCE_DEFAULT = 1e-6 # stopping criterion for MLE optimization
THETA0_DEFAULT = 0.5
THETAL_DEFAULT = 1e-5
THETAU_DEFAULT = 50
if hasattr(linalg, 'solve_triangular'):
# only in scipy since 0.9
solve_triangular = linalg.solve_triangular
else:
# slower, but works
def solve_triangular(x, y, lower=True):
"""Solve triangular."""
return linalg.solve(x, y)
def constant_regression(x):
"""
Zero order polynomial (constant, p = 1) regression model.
x --> f(x) = 1
Parameters
----------
x : array_like
Input data.
Returns
-------
array_like
Constant regression output.
"""
x = np.asarray(x, dtype=np.float)
n_eval = x.shape[0]
f = np.ones([n_eval, 1])
return f
def linear_regression(x):
"""
First order polynomial (linear, p = n+1) regression model.
x --> f(x) = [ 1, x_1, ..., x_n ].T
Parameters
----------
x : array_like
Input data.
Returns
-------
array_like
Linear regression output.
"""
x = np.asarray(x, dtype=np.float)
n_eval = x.shape[0]
f = np.hstack([np.ones([n_eval, 1]), x])
return f
def squared_exponential_correlation(theta, d):
"""
Squared exponential correlation model (Radial Basis Function).
(Infinitely differentiable stochastic process, very smooth)::
n
theta, dx --> r(theta, dx) = exp( sum - theta_i * (dx_i)^2 )
i = 1
Parameters
----------
theta : array_like
An array with shape 1 (isotropic) or n (anisotropic) giving the
autocorrelation parameter(s).
d : array_like
An array with shape (n_eval, n_features) giving the componentwise
distances between locations x and x' at which the correlation model
should be evaluated.
Returns
-------
array_like
An array with shape (n_eval, ) containing the values of the
autocorrelation model.
"""
theta = np.asarray(theta, dtype=np.float)
d = np.asarray(d, dtype=np.float)
if d.ndim > 1:
n_features = d.shape[1]
else:
n_features = 1
if theta.size == 1:
return np.exp(-theta[0] * np.sum(d ** 2, axis=1))
elif theta.size != n_features:
raise ValueError("Length of theta must be 1 or %s" % n_features)
else:
return np.exp(-np.sum(theta.reshape(1, n_features) * d ** 2, axis=1))
def l1_cross_distances(X, Y=None):
"""
Compute the nonzero componentwise L1 cross-distances between the vectors in X and Y.
Parameters
----------
X : array_like
An array with shape (n_samples_X, n_features).
Y : array_like
An array with shape (n_samples_Y, n_features).
Returns
-------
array with shape (n_samples * (n_samples - 1) / 2, n_features)
The array of componentwise L1 cross-distances.
"""
X = array2d(X)
if Y is None:
n_samples, n_features = X.shape
n_nonzero_cross_dist = n_samples * (n_samples - 1) // 2
D = np.zeros((n_nonzero_cross_dist, n_features))
ll_1 = 0
for k in range(n_samples - 1):
ll_0 = ll_1
ll_1 = ll_0 + n_samples - k - 1
D[ll_0:ll_1] = np.abs(X[k] - X[(k + 1):])
else:
Y = array2d(Y)
n_samples_X, n_features_X = X.shape
n_samples_Y, n_features_Y = Y.shape
if n_features_X != n_features_Y:
raise ValueError("X and Y must have the same dimensions.")
n_features = n_features_X
n_nonzero_cross_dist = n_samples_X * n_samples_Y
D = np.zeros((n_nonzero_cross_dist, n_features))
ll_1 = 0
for k in range(n_samples_X):
ll_0 = ll_1
ll_1 = ll_0 + n_samples_Y # - k - 1
D[ll_0:ll_1] = np.abs(X[k] - Y)
return D
class MultiFiCoKriging(object):
"""
Integrate the Multi-Fidelity Co-Kriging method described in [LeGratiet2013].
Parameters
----------
regr : str or callable, optional
A regression function returning an array of outputs of the linear
regression functional basis for Universal Kriging purpose.
regr is assumed to be the same for all levels of code.
Default assumes a simple constant regression trend.
Available built-in regression models are:
'constant', 'linear'.
rho_regr : str or callable, optional
A regression function returning an array of outputs of the linear
regression functional basis. Defines the regression function for the
autoregressive parameter rho.
rho_regr is assumed to be the same for all levels of code.
Default assumes a simple constant regression trend.
Available built-in regression models are:
'constant', 'linear'.
normalize : bool, optional
When true, normalize X and Y so that the mean is at zero.
theta : double, array_like or list, optional
Value of correlation parameters if they are known; no optimization is run.
Default is None, so that optimization is run.
if double: value is replicated for all features and all levels.
if array_like: an array with shape (n_features, ) for
isotropic calculation. It is replicated for all levels.
if list: a list of nlevel arrays specifying value for each level.
theta0 : double, array_like or list, optional
Starting point for the maximum likelihood estimation of the
best set of parameters.
Default is None and meaning use of the default 0.5*np.ones(n_features)
if double: value is replicated for all features and all levels.
if array_like: an array with shape (n_features, ) for
isotropic calculation. It is replicated for all levels.
if list: a list of nlevel arrays specifying value for each level.
thetaL : double, array_like or list, optional
Lower bound on the autocorrelation parameters for maximum
likelihood estimation.
Default is None meaning use of the default 1e-5*np.ones(n_features).
if double: value is replicated for all features and all levels.
if array_like: An array with shape matching theta0's. It is replicated
for all levels of code.
if list: a list of nlevel arrays specifying value for each level.
thetaU : double, array_like or list, optional
Upper bound on the autocorrelation parameters for maximum
likelihood estimation.
Default is None meaning use of default value 50*np.ones(n_features).
if double: value is replicated for all features and all levels.
if array_like: An array with shape matching theta0's. It is replicated
for all levels of code.
if list: a list of nlevel arrays specifying value for each level.
Attributes
----------
corr : object
Correlation function to use, default is squared_exponential_correlation.
n_features : ndarray
Number of features for each fidelity level.
n_samples : ndarray
Number of samples for each fidelity level.
nlevel : int
Number of fidelity levels.
normalize : bool, optional
When true, normalize X and Y so that the mean is at zero.
regr : str or callable
A regression function returning an array of outputs of the linear
regression functional basis for Universal Kriging purpose.
regr is assumed to be the same for all levels of code.
Default assumes a simple constant regression trend.
Available built-in regression models are:
'constant', 'linear'
rho_regr : str or callable or None
A regression function returning an array of outputs of the linear
regression functional basis. Defines the regression function for the
autoregressive parameter rho.
rho_regr is assumed to be the same for all levels of code.
Default assumes a simple constant regression trend.
Available built-in regression models are:
'constant', 'linear'
theta : double, array_like or list or None
Value of correlation parameters if they are known; no optimization is run.
Default is None, so that optimization is run.
if double: value is replicated for all features and all levels.
if array_like: an array with shape (n_features, ) for
isotropic calculation. It is replicated for all levels.
if list: a list of nlevel arrays specifying value for each level
theta0 : double, array_like or list or None
Starting point for the maximum likelihood estimation of the
best set of parameters.
Default is None and meaning use of the default 0.5*np.ones(n_features)
if double: value is replicated for all features and all levels.
if array_like: an array with shape (n_features, ) for
isotropic calculation. It is replicated for all levels.
if list: a list of nlevel arrays specifying value for each level
thetaL : double, array_like or list or None
Lower bound on the autocorrelation parameters for maximum
likelihood estimation.
Default is None meaning use of the default 1e-5*np.ones(n_features).
if double: value is replicated for all features and all levels.
if array_like: An array with shape matching theta0's. It is replicated
for all levels of code.
if list: a list of nlevel arrays specifying value for each level
thetaU : double, array_like or list or None
Upper bound on the autocorrelation parameters for maximum
likelihood estimation.
Default is None meaning use of default value 50*np.ones(n_features).
if double: value is replicated for all features and all levels.
if array_like: An array with shape matching theta0's. It is replicated
for all levels of code.
if list: a list of nlevel arrays specifying value for each level
X_mean : float
Mean of the low fidelity training data for X.
X_std : float
Standard deviation of the low fidelity training data for X.
y_mean : float
Mean of the low fidelity training data for y.
y_std : float
Standard deviation of the low fidelity training data for y.
_nfev : int
Number of function evaluations.
Notes
-----
Implementation is based on the Package Scikit-Learn
(Author: Vincent Dubourg, vincent.dubourg@gmail.com) which translates
the DACE Matlab toolbox, see [NLNS2002]_.
References
----------
.. [NLNS2002] H. B. Nielsen, S. N. Lophaven, and J. Sondergaard.
`DACE - A MATLAB Kriging Toolbox.` (2002)
http://www2.imm.dtu.dk/~hbn/dace/dace.pdf
.. [WBSWM1992] W. J. Welch, R. J. Buck, J. Sacks, H. P. Wynn, T. J. Mitchell,
and M. D. Morris (1992). "Screening, predicting, and computer experiments."
`Technometrics,` 34(1) 15--25.
http://www.jstor.org/pss/1269548
.. [LeGratiet2013] L. Le Gratiet (2013). "Multi-fidelity Gaussian process
regression for computer experiments."
PhD thesis, Universite Paris-Diderot-Paris VII.
.. [TBKH2011] Toal, D. J., Bressloff, N. W., Keane, A. J., & Holden, C. M. E. (2011).
"The development of a hybridized particle swarm for kriging hyperparameter
tuning." `Engineering optimization`, 43(6), 675-699.
Examples
--------
>>> from openmdao.surrogate_models.multifi_cokriging import MultiFiCoKriging
>>> import numpy as np
>>> # Xe: DOE for expensive code (nested in Xc)
>>> # Xc: DOE for cheap code
>>> # ye: expensive response
>>> # yc: cheap response
>>> Xe = np.array([[0],[0.4],[1]])
>>> Xc = np.vstack((np.array([[0.1],[0.2],[0.3],[0.5],[0.6],[0.7],[0.8],[0.9]]),Xe))
>>> ye = ((Xe*6-2)**2)*np.sin((Xe*6-2)*2)
>>> yc = 0.5*((Xc*6-2)**2)*np.sin((Xc*6-2)*2)+(Xc-0.5)*10. - 5
>>> model = MultiFiCoKriging(theta0=1, thetaL=1e-5, thetaU=50.)
>>> model.fit([Xc, Xe], [yc, ye])
>>> # Prediction on x=0.05
>>> np.abs(float(model.predict([0.05])[0])- ((0.05*6-2)**2)*np.sin((0.05*6-2)*2)) < 0.05
True
"""
_regression_types = {
'constant': constant_regression,
'linear': linear_regression
}
def __init__(self, regr='constant', rho_regr='constant', normalize=True,
theta=None, theta0=None, thetaL=None, thetaU=None):
"""
Initialize all attributes.
"""
self.corr = squared_exponential_correlation
self.regr = regr
self.rho_regr = rho_regr
self.theta = theta
self.theta0 = theta0
self.thetaL = thetaL
self.thetaU = thetaU
self.normalize = normalize
self.X_mean = 0
self.X_std = 1
self.y_mean = 0
self.y_std = 1
self.n_features = None
self.n_samples = None
self.nlevel = None
self._nfev = 0
def _build_R(self, lvl, theta):
"""
Build the correlation matrix with given theta for the specified level.
Parameters
----------
lvl : int
Level of fidelity
theta : array_like
An array containing the autocorrelation parameters at which the
Gaussian Process model parameters should be determined.
Default uses the built-in autocorrelation parameters
(ie ``theta = self.theta``).
Returns
-------
ndarray
Correlation matrix.
"""
D = self.D[lvl]
n_samples = self.n_samples[lvl]
R = np.eye(n_samples) * (1. + NUGGET)
corr = squareform(self.corr(theta, D))
R = R + corr
return R
def fit(self, X, y, initial_range=INITIAL_RANGE_DEFAULT, tol=TOLERANCE_DEFAULT):
"""
Implement the Multi-Fidelity co-kriging model fitting method.
Parameters
----------
X : list of double array_like elements
A list of arrays with the input at which observations were made, from lowest
fidelity to highest fidelity. Designs must be nested
with X[i] = np.vstack([..., X[i+1]).
y : list of double array_like elements
A list of arrays with the observations of the scalar output to be predicted,
from lowest fidelity to highest fidelity.
initial_range : float
Initial range for the optimizer.
tol : float
Optimizer terminates when the tolerance tol is reached.
"""
# Run input checks
# Transforms floats and arrays in lists to have a multifidelity
# structure
self._check_list_structure(X, y)
# Checks if all parameters are structured as required
self._check_params()
X = self.X
y = self.y
nlevel = self.nlevel
n_samples = self.n_samples
# initialize lists
self.beta = nlevel * [0]
self.beta_rho = nlevel * [None]
self.beta_regr = nlevel * [None]
self.C = nlevel * [0]
self.D = nlevel * [0]
self.F = nlevel * [0]
self.p = nlevel * [0]
self.q = nlevel * [0]
self.G = nlevel * [0]
self.sigma2 = nlevel * [0]
self._R_adj = nlevel * [None]
# Training data will be normalized using statistical quantities from the low fidelity set.
if self.normalize:
self.X_mean = X_mean = np.mean(X[0], axis=0)
self.X_std = X_std = np.std(X[0], axis=0)
self.y_mean = y_mean = np.mean(y[0], axis=0)
self.y_std = y_std = np.std(y[0], axis=0)
X_std[X_std == 0.] = 1.
y_std[y_std == 0.] = 1.
for lvl in range(nlevel):
if self.normalize:
X[lvl] = (X[lvl] - X_mean) / X_std
y[lvl] = (y[lvl] - y_mean) / y_std
# Calculate matrix of distances D between samples
self.D[lvl] = l1_cross_distances(X[lvl])
if (np.min(np.sum(self.D[lvl], axis=1)) == 0.):
raise ValueError("Multiple input features cannot have the same value.")
# Regression matrix and parameters
self.F[lvl] = self.regr(X[lvl])
self.p[lvl] = self.F[lvl].shape[1]
# Concatenate the autoregressive part for levels > 0
if lvl > 0:
F_rho = self.rho_regr(X[lvl])
self.q[lvl] = F_rho.shape[1]
self.F[lvl] = np.hstack((F_rho * np.dot((self.y[lvl - 1])[-n_samples[lvl]:],
np.ones((1, self.q[lvl]))), self.F[lvl]))
else:
self.q[lvl] = 0
n_samples_F_i = self.F[lvl].shape[0]
if n_samples_F_i != n_samples[lvl]:
raise ValueError("Number of rows in F and X do not match. Most "
"likely something is going wrong with the "
"regression model.")
if int(self.p[lvl] + self.q[lvl]) >= n_samples_F_i:
raise ValueError(f"Ordinary least squares problem is undetermined "
f"n_samples={n_samples[lvl]} must be greater than the regression"
f" model size p+q={self.p[lvl] + self.q[lvl]}.")
# Set attributes
self.X = X
self.y = y
self.rlf_value = np.zeros(nlevel)
for lvl in range(nlevel):
# Determine Gaussian Process model parameters
if self.theta[lvl] is None:
# Maximum Likelihood Estimation of the parameters
sol = self._max_rlf(lvl=lvl, initial_range=initial_range, tol=tol)
self.theta[lvl] = sol['theta']
self.rlf_value[lvl] = sol['rlf_value']
if np.isinf(self.rlf_value[lvl]):
raise ValueError("Bad parameter region. Try increasing upper bound")
else:
self.rlf_value[lvl] = self.rlf(lvl=lvl)
if np.isinf(self.rlf_value[lvl]):
raise ValueError("Bad point. Try increasing theta0.")
return
def rlf(self, lvl, theta=None):
"""
Determine BLUP parameters and evaluate negative reduced likelihood function for theta.
Maximizing this function wrt the autocorrelation parameters theta is
equivalent to maximizing the likelihood of the assumed joint Gaussian
distribution of the observations y evaluated onto the design of
experiments X.
Parameters
----------
lvl : int
Level of fidelity.
theta : array_like, optional
An array containing the autocorrelation parameters at which the
Gaussian Process model parameters should be determined.
Default uses the built-in autocorrelation parameters
(ie ``theta = self.theta``).
Returns
-------
double
The value of the negative concentrated reduced likelihood function
associated to the given autocorrelation parameters theta.
"""
if theta is None:
# Use built-in autocorrelation parameters
theta = self.theta[lvl]
# Initialize output
rlf_value = 1e20
# Retrieve data
n_samples = self.n_samples[lvl]
y = self.y[lvl]
F = self.F[lvl]
p = self.p[lvl]
q = self.q[lvl]
R = self._build_R(lvl, theta)
try:
C = linalg.cholesky(R, lower=True)
except linalg.LinAlgError:
_logger.warning(('Cholesky decomposition of R at level %i failed' % lvl) +
' with theta=' + str(theta))
return rlf_value
# Get generalized least squares solution
Ft = solve_triangular(C, F, lower=True)
Yt = solve_triangular(C, y, lower=True)
try:
Q, G = linalg.qr(Ft, econ=True)
except TypeError: # qr() got an unexpected keyword argument 'econ'
# DeprecationWarning: qr econ argument will be removed after scipy
# 0.7. The economy transform will then be available through the
# mode='economic' argument.
Q, G = linalg.qr(Ft, mode='economic')
pass
# Universal Kriging
beta = solve_triangular(G, np.dot(Q.T, Yt))
err = Yt - np.dot(Ft, beta)
err2 = np.dot(err.T, err)[0, 0]
self._err = err
sigma2 = err2 / (n_samples - p - q)
detR = ((np.diag(C))**(2. / n_samples)).prod()
rlf_value = (n_samples - p - q) * np.log10(sigma2) \
+ n_samples * np.log10(detR)
self.beta_rho[lvl] = beta[:q]
self.beta_regr[lvl] = beta[q:]
self.beta[lvl] = beta
self.sigma2[lvl] = sigma2
self.C[lvl] = C
self.G[lvl] = G
return rlf_value
def _max_rlf(self, lvl, initial_range, tol):
"""
Estimate autocorrelation parameter theta as maximizer of the reduced likelihood function.
(Minimization of the negative reduced likelihood function is used for convenience.)
Parameters
----------
lvl : int
Level of fidelity
initial_range : float
Initial range of the optimizer
tol : float
Optimizer terminates when the tolerance tol is reached.
Returns
-------
array_like
The optimal hyperparameters.
double
The optimal negative reduced likelihood function value.
dict
res['theta']: optimal theta
res['rlf_value']: optimal value for likelihood
"""
# Initialize input
thetaL = self.thetaL[lvl]
thetaU = self.thetaU[lvl]
def rlf_transform(x):
return self.rlf(theta=10.**x, lvl=lvl)
# Use specified starting point as first guess
theta0 = self.theta0[lvl]
x0 = np.log10(theta0[0])
constraints = []
for i in range(theta0.size):
constraints.append({'type': 'ineq', 'fun': lambda log10t, i=i:
log10t[i] - np.log10(thetaL[0][i])})
constraints.append({'type': 'ineq', 'fun': lambda log10t, i=i:
np.log10(thetaU[0][i]) - log10t[i]})
constraints = tuple(constraints)
sol = minimize(rlf_transform, x0, method='COBYLA',
constraints=constraints,
options={'rhobeg': initial_range,
'tol': tol, 'disp': 0})
log10_optimal_x = sol['x']
optimal_rlf_value = sol['fun']
self._nfev += sol['nfev']
optimal_theta = 10. ** log10_optimal_x
res = {}
res['theta'] = optimal_theta
res['rlf_value'] = optimal_rlf_value
return res
def predict(self, X, eval_MSE=True):
"""
Perform the predictions of the kriging model on X.
Parameters
----------
X : array_like
An array with shape (n_eval, n_features) giving the point(s) at
which the prediction(s) should be made.
eval_MSE : bool, optional
A boolean specifying whether the Mean Squared Error should be
evaluated or not. Default assumes evalMSE is True.
Returns
-------
array_like
An array with shape (n_eval, ) with the Best Linear Unbiased
Prediction at X. If all_levels is set to True, an array
with shape (n_eval, nlevel) giving the BLUP for all levels.
array_like, optional (if eval_MSE is True)
An array with shape (n_eval, ) with the Mean Squared Error at X.
If all_levels is set to True, an array with shape (n_eval, nlevel)
giving the MSE for all levels.
"""
X = array2d(X)
nlevel = self.nlevel
n_eval, n_features_X = X.shape
# Normalize
if self.normalize:
X = (X - self.X_mean) / self.X_std
# Calculate kriging mean and variance at level 0
mu = np.zeros((n_eval, nlevel))
f = self.regr(X)
f0 = self.regr(X)
dx = l1_cross_distances(X, Y=self.X[0])
# Get regression function and correlation
F = self.F[0]
C = self.C[0]
beta = self.beta[0]
Ft = solve_triangular(C, F, lower=True)
yt = solve_triangular(C, self.y[0], lower=True)
r_ = self.corr(self.theta[0], dx).reshape(n_eval, self.n_samples[0])
gamma = solve_triangular(C.T, yt - np.dot(Ft, beta), lower=False)
# Scaled predictor
mu[:, 0] = (np.dot(f, beta) + np.dot(r_, gamma)).ravel()
if eval_MSE:
MSE = np.zeros((n_eval, nlevel))
r_t = solve_triangular(C, r_.T, lower=True)
G = self.G[0]
u_ = solve_triangular(G.T, f.T - np.dot(Ft.T, r_t), lower=True)
MSE[:, 0] = self.sigma2[0] * \
(1 - (r_t**2).sum(axis=0) + (u_**2).sum(axis=0))
# Calculate recursively kriging mean and variance at level i
for i in range(1, nlevel):
C = self.C[i]
F = self.F[i]
g = self.rho_regr(X)
dx = l1_cross_distances(X, Y=self.X[i])
r_ = self.corr(self.theta[i], dx).reshape(
n_eval, self.n_samples[i])
f = np.vstack((g.T * mu[:, i - 1], f0.T))
Ft = solve_triangular(C, F, lower=True)
yt = solve_triangular(C, self.y[i], lower=True)
r_t = solve_triangular(C, r_.T, lower=True)
G = self.G[i]
beta = self.beta[i]
# scaled predictor
mu[:, i] = (np.dot(f.T, beta)
+ np.dot(r_t.T, yt - np.dot(Ft, beta))).ravel()
if eval_MSE:
Q_ = (np.dot((yt - np.dot(Ft, beta)).T,
yt - np.dot(Ft, beta)))[0, 0]
u_ = solve_triangular(G.T, f - np.dot(Ft.T, r_t), lower=True)
sigma2_rho = np.dot(g,
self.sigma2[
i] * linalg.inv(np.dot(G.T, G))[:self.q[i], :self.q[i]]
+ np.dot(beta[:self.q[i]], beta[:self.q[i]].T))
sigma2_rho = (sigma2_rho * g).sum(axis=1)
MSE[:, i] = sigma2_rho * MSE[:, i - 1] \
+ Q_ / (2 * (self.n_samples[i] - self.p[i] - self.q[i])) \
* (1 - (r_t**2).sum(axis=0)) \
+ self.sigma2[i] * (u_**2).sum(axis=0)
# scaled predictor
for i in range(nlevel): # Predictor
mu[:, i] = self.y_mean + self.y_std * mu[:, i]
if eval_MSE:
MSE[:, i] = self.y_std**2 * MSE[:, i]
if eval_MSE:
return mu[:, -1].reshape((n_eval, 1)), MSE[:, -1].reshape((n_eval, 1))
else:
return mu[:, -1].reshape((n_eval, 1))
def _check_list_structure(self, X, y):
"""
Transform floats and arrays in the training data lists to have a multifidelity structure.
Parameters
----------
X : list of double array_like elements
A list of arrays with the input at which observations were made, from lowest
fidelity to highest fidelity. Designs must be nested
with X[i] = np.vstack([..., X[i+1])
y : list of double array_like elements
A list of arrays with the observations of the scalar output to be predicted,
from lowest fidelity to highest fidelity.
"""
if type(X) is not list:
nlevel = 1
X = [X]
else:
nlevel = len(X)
if type(y) is not list:
y = [y]
if len(X) != len(y):
raise ValueError("X and y must have the same length.")
n_samples = np.zeros(nlevel, dtype=int)
n_features = np.zeros(nlevel, dtype=int)
n_samples_y = np.zeros(nlevel, dtype=int)
for i in range(nlevel):
n_samples[i], n_features[i] = X[i].shape
if i > 0 and n_features[i] != n_features[i - 1]:
raise ValueError("All X must have the same number of columns.")
y[i] = np.asarray(y[i]).ravel()[:, np.newaxis]
n_samples_y[i] = y[i].shape[0]
if n_samples[i] != n_samples_y[i]:
raise ValueError("X and y must have the same number of rows.")
self.n_features = n_features[0]
if type(self.theta) is not list:
self.theta = nlevel * [self.theta]
elif len(self.theta) != nlevel:
raise ValueError(f"theta must be a list of {nlevel} element(s).")
if type(self.theta0) is not list:
self.theta0 = nlevel * [self.theta0]
elif len(self.theta0) != nlevel:
raise ValueError(f"theta0 must be a list of {nlevel} element(s).")
if type(self.thetaL) is not list:
self.thetaL = nlevel * [self.thetaL]
elif len(self.thetaL) != nlevel:
raise ValueError(f"thetaL must be a list of {nlevel} element(s).")
if type(self.thetaU) is not list:
self.thetaU = nlevel * [self.thetaU]
elif len(self.thetaU) != nlevel:
raise ValueError(f"thetaU must be a list of {nlevel} element(s).")
self.nlevel = nlevel
self.X = X[:]
self.y = y[:]
self.n_samples = n_samples
return
def _check_params(self):
"""
Perform sanity checks on all parameters.
"""
# Check regression model
if not callable(self.regr):
if self.regr in self._regression_types:
self.regr = self._regression_types[self.regr]
else:
raise ValueError(f"regr should be one of {self._regression_types.keys()} or "
f"callable, {self.regr} was given.")
# Check rho regression model
if not callable(self.rho_regr):
if self.rho_regr in self._regression_types:
self.rho_regr = self._regression_types[self.rho_regr]
else:
raise ValueError(f"regr should be one of {self._regression_types.keys()} or "
f"callable, {self.rho_regr} was given.")
for i in range(self.nlevel):
# Check correlation parameters
if self.theta[i] is not None:
self.theta[i] = array2d(self.theta[i])
if np.any(self.theta[i] <= 0):
raise ValueError("theta must be strictly positive.")
if self.theta0[i] is not None:
self.theta0[i] = array2d(self.theta0[i])
if np.any(self.theta0[i] <= 0):
raise ValueError("theta0 must be strictly positive.")
else:
self.theta0[i] = array2d(self.n_features * [THETA0_DEFAULT])
lth = self.theta0[i].size
if self.thetaL[i] is not None:
self.thetaL[i] = array2d(self.thetaL[i])
if self.thetaL[i].size != lth:
raise ValueError("theta0 and thetaL must have the same length.")
else:
self.thetaL[i] = array2d(self.n_features * [THETAL_DEFAULT])
if self.thetaU[i] is not None:
self.thetaU[i] = array2d(self.thetaU[i])
if self.thetaU[i].size != lth:
raise ValueError("theta0 and thetaU must have the same length.")
else:
self.thetaU[i] = array2d(self.n_features * [THETAU_DEFAULT])
if np.any(self.thetaL[i] <= 0) or np.any(self.thetaU[i] < self.thetaL[i]):
raise ValueError("The bounds must satisfy O < thetaL <= thetaU.")
return
class MultiFiCoKrigingSurrogate(MultiFiSurrogateModel):
"""
OpenMDAO adapter of multi-fidelity recursive cokriging method described in [LeGratiet2013].
See MultiFiCoKriging class.
Parameters
----------
**kwargs : keyword args
Some implementations of record_derivatives need additional args.
Attributes
----------
model : MultiFiCoKriging
Contains MultiFiCoKriging surrogate.
"""
def __init__(self, **kwargs):
"""
Initialize all attributes.
"""
super().__init__(**kwargs)
self.model = None
def _declare_options(self):
"""
Declare options before kwargs are processed in the init method.
"""
opt = self.options
opt.declare('normalize', default=True, types=bool,
desc="When true, normalize X and Y so that the mean is at zero.")
opt.declare('regr', default='constant', types=(object, ),
desc="A regression function returning an array of outputs of the linear "
"regression functional basis for Universal Kriging purpose. regr is assumed "
"to be the same for all levels of code. Default assumes a simple constant "
"regression trend. Available built-in regression models can be accessed by "
"setting this option to the strings 'constant' or 'linear'")
opt.declare('rho_regr', default='constant', types=(object, ),
desc="A regression function returning an array of outputs of the linear "
"regression functional basis. Defines the regression function for the "
"autoregressive parameter rho.. regr is assumed to be the same for all levels "
"of code. Default assumes a simple constant regression trend. Available "
"built-in regression models can be accessed by setting this option to the "
"strings 'constant' or 'linear'")
opt.declare('theta', default=None, allow_none=True,
desc="Value of correlation parameters. If they are known, then no "
"optimization is run. Default is None, so that optimization is run. if double, "
"then value is replicated for all features and all levels. if array_like, "
"then an array with shape (n_features, ) for isotropic calculation. It is "
"replicated for all levels. if list, then a list of nlevel arrays specifying "
"value for each level")
opt.declare('theta0', default=None, allow_none=True,
desc="Starting point for the maximum likelihood estimation of the best set "
"of parameters. "
"Default is None and meaning use of the default 0.5*np.ones(n_features) "
"if double: value is replicated for all features and all levels. "
"if array_like: an array with shape (n_features, ) for "
"isotropic calculation. It is replicated for all levels. "
"if list: a list of nlevel arrays specifying value for each level")
opt.declare('thetaL', default=None, allow_none=True,
desc="Lower bound on the autocorrelation parameters for maximum "
"likelihood estimation."
"Default is None meaning use of the default 1e-5*np.ones(n_features). "
"if double: value is replicated for all features and all levels. "
"if array_like: An array with shape matching theta0s. It is replicate "
"for all levels of code. "
"if list: a list of nlevel arrays specifying value for each level")
opt.declare('thetaU', default=None, allow_none=True,
desc="Upper bound on the autocorrelation parameters for maximum "
"likelihood estimation. "
"Default is None meaning use of default value 50*np.ones(n_features). "
"if double: value is replicated for all features and all levels. "
"if array_like: An array with shape matching theta0's. It is replicated "
"for all levels of code. "
"if list: a list of nlevel arrays specifying value for each level")
opt.declare('tolerance', default=TOLERANCE_DEFAULT,
desc='Optimizer terminates when the tolerance tol is reached.')
opt.declare('initial_range', default=INITIAL_RANGE_DEFAULT,
desc='Initial range for the optimizer.')
def predict(self, new_x):
"""
Calculate a predicted value of the response based on the current trained model.
Parameters
----------
new_x : array_like
An array with shape (n_eval, n_features) giving the point(s) at
which the prediction(s) should be made.
Returns
-------
array_like
An array with shape (n_eval, ) with the Best Linear Unbiased
Prediction at X. If all_levels is set to True, an array
with shape (n_eval, nlevel) giving the BLUP for all levels.
array_like
An array with shape (n_eval, ) with the square root of the Mean Squared Error at X.
"""
Y_pred, MSE = self.model.predict([new_x])
return Y_pred, np.sqrt(np.abs(MSE))
def train_multifi(self, X, Y):
"""
Train the surrogate model with the given set of inputs and outputs.
Parameters
----------
X : list of double array_like elements
A list of arrays with the input at which observations were made, from highest
fidelity to lowest fidelity. Designs must be nested