-
Notifications
You must be signed in to change notification settings - Fork 5
/
EnsembleKalmanSchemes.jl
executable file
·2872 lines (2565 loc) · 119 KB
/
EnsembleKalmanSchemes.jl
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
##############################################################################################
module EnsembleKalmanSchemes
##############################################################################################
# imports and exports
using Random, Distributions, Statistics, StatsBase
using LinearAlgebra, SparseArrays
using ..DataAssimilationBenchmarks
using Optim, LineSearches, LinearAlgebra
export analyze_ens, analyze_ens_param, rand_orth, inflate_state!, inflate_param!,
transform_R, ens_gauss_newton,
square_root, square_root_inv,
ensemble_filter, ls_smoother_classic,
ls_smoother_single_iteration, ls_smoother_gauss_newton
##############################################################################################
# Main methods, debugged and validated
##############################################################################################
"""
analyze_ens(ens::ArView(T), truth::VecA(T)) where T <: Float64
Computes the ensemble state RMSE as compared with truth twin, and the ensemble spread.
```
return rmse, spread
```
Note: the ensemble `ens` should only include the state vector components to compare with the
truth twin state vector `truth`, without replicates of the model parameters. These can be
passed as an [`ArView`](@ref) for efficient memory usage.
"""
function analyze_ens(ens::ArView(T), truth::VecA(T)) where T <: Float64
# infer the shapes
sys_dim, N_ens = size(ens)
# compute the ensemble mean
x_bar = mean(ens, dims=2)
# compute the RMSE of the ensemble mean
rmse = rmsd(truth, x_bar)
# compute the spread as in whitaker & louge 98 by the standard deviation
# of the mean square deviation of the ensemble from its mean
spread = sqrt( ( 1.0 / (N_ens - 1.0) ) * sum(mean((ens .- x_bar).^2.0, dims=1)))
# return the tuple pair
rmse, spread
end
##############################################################################################
"""
analyze_ens_param(ens::ArView(T), truth::VecA(T)) where T <: Float64
Computes the ensemble parameter RMSE as compared with truth twin, and the ensemble spread.
```
return rmse, spread
```
Note: the ensemble `ens` should only include the extended state vector components
consisting of model parameter replicates to compare with the truth twin's governing
model parameters `truth`. These can be passed as an [`ArView`](@ref) for
efficient memory usage.
"""
function analyze_ens_param(ens::ArView(T), truth::VecA(T)) where T <: Float64
# infer the shapes
param_dim, N_ens = size(ens)
# compute the ensemble mean
x_bar = mean(ens, dims=2)
# compute the RMSE of relative to the magnitude of the parameter
rmse = sqrt( mean( (truth - x_bar).^2.0 ./ truth.^2.0 ) )
# compute the spread as in whitaker & louge 98 by the standard deviation
# of the mean square deviation of the ensemble from its mean,
# with the weight by the size of the parameter square
spread = sqrt( ( 1.0 / (N_ens - 1.0) ) *
sum(mean( (ens .- x_bar).^2.0 ./
(ones(param_dim, N_ens) .* truth.^2.0), dims=1)))
# return the tuple pair
rmse, spread
end
##############################################################################################
"""
rand_orth(N_ens::Int64)
This generates a random, mean-preserving, orthogonal matrix as in [Sakov et al.
2008](https://journals.ametsoc.org/view/journals/mwre/136/3/2007mwr2021.1.xml), depending on
the esemble size `N_ens`.
```
return U
```
"""
function rand_orth(N_ens::Int64)
# generate the random, mean preserving orthogonal transformation within the
# basis given by the B matrix
Q = rand(Normal(), N_ens - 1, N_ens - 1)
Q, R = qr!(Q)
U_p = zeros(N_ens, N_ens)
U_p[1, 1] = 1.0
U_p[2:end, 2:end] = Q
# generate the B basis for which the first basis vector is the vector of 1/sqrt(N)
b_1 = ones(N_ens) / sqrt(N_ens)
B = zeros(N_ens, N_ens)
B[:, 1] = b_1
# note, this uses the "full" QR decomposition so that the singularity is encoded in R
# and B is a full-size orthogonal matrix
B, R = qr!(B)
U = B * U_p * transpose(B)
U
end
##############################################################################################
"""
inflate_state!(ens::ArView(T), inflation::Float64, sys_dim::Int64,
state_dim::Int64) where T <: Float64
Applies multiplicative covariance inflation to the state components of the ensemble matrix.
```
return ens
```
The first index of the ensemble matrix `ens` corresponds to the length `sys_dim` (extended)
state dimension while the second index corresponds to the ensemble dimension. Dynamic state
variables are assumed to be in the leading `state_dim` rows of `ens`, while extended state
parameter replicates are after. Multiplicative inflation is performed only in the leading
components of the ensemble anomalies from the ensemble mean, in-place in memory.
"""
function inflate_state!(ens::ArView(T), inflation::Float64, sys_dim::Int64,
state_dim::Int64) where T <: Float64
if inflation != 1.0
x_mean = mean(ens, dims=2)
X = ens .- x_mean
infl = Matrix(1.0I, sys_dim, sys_dim)
infl[1:state_dim, 1:state_dim] .*= inflation
ens .= x_mean .+ infl * X
return ens
end
end
##############################################################################################
"""
inflate_param!(ens::ArView(T), inflation::Float64, sys_dim::Int64,
state_dim::Int64) where T <: Float64
Applies multiplicative covariance inflation to parameter replicates in the ensemble matrix.
```
return ens
```
The first index of the ensemble matrix `ens` corresponds to the length `sys_dim` (extended)
state dimension while the second index corresponds to the ensemble dimension. Dynamic state
variables are assumed to be in the leading `state_dim` rows of `ens`, while extended state
parameter replicates are after. Multiplicative inflation is performed only in the trailing
`state_dim + 1: state_dim` components of the ensemble anomalies from the ensemble mean,
in-place in memory.
"""
function inflate_param!(ens::ArView(T), inflation::Float64, sys_dim::Int64,
state_dim::Int64) where T <: Float64
if inflation == 1.0
return ens
else
x_mean = mean(ens, dims=2)
X = ens .- x_mean
infl = Matrix(1.0I, sys_dim, sys_dim)
infl[state_dim+1: end, state_dim+1: end] .*= inflation
ens .= x_mean .+ infl * X
return ens
end
end
##############################################################################################
"""
square_root(M::CovM(T)) where T <: Real
Computes the square root of covariance matrices with parametric type.
Multiple dispatches for the method are defined according to the sub-type of [`CovM`](@ref),
where the square roots of `UniformScaling` and `Diagonal` covariance matrices are computed
directly, while the square roots of the more general class of `Symmetric` covariance
matrices are computed via the singular value decomposition, for stability and accuracy
for close-to-singular matrices.
```
return S
```
"""
function square_root(M::UniformScaling{T}) where T <: Real
S = M^0.5
end
function square_root(M::Diagonal{T, Vector{T}}) where T <: Real
S = sqrt(M)
end
function square_root(M::Symmetric{T, Matrix{T}}) where T <: Real
F = svd(M)
S = Symmetric(F.U * Diagonal(sqrt.(F.S)) * F.Vt)
end
##############################################################################################
"""
square_root_inv(M::CovM(T); sq_rt::Bool=false, inverse::Bool=false,
full::Bool=false) where T <: Real
Computes the square root inverse of covariance matrices with parametric type.
Multiple dispatches for the method are defined according to the sub-type of [`CovM`](@ref),
where the square root inverses of `UniformScaling` and `Diagonal` covariance matrices
are computed directly, while the square root inverses of the more general class of
`Symmetric` covariance matrices are computed via the singular value decomposition, for
stability and accuracy for close-to-singular matrices. This will optionally return a
computation of the inverse and the square root itself all as a byproduct of the singular
value decomposition for efficient numerical computation of ensemble
analysis / update routines.
Optional keyword arguments are specified as:
* `sq_rt=true` returns the matrix square root in addition to the square root inverse
* `inverse=true` returns the matrix inverse in addition to the square root inverse
* `full=true` returns the square root and the matrix inverse in addition to the square
root inverse
and are evaluated in the above order.
Output follows control flow:
```
if sq_rt
return S_inv, S
elseif inverse
return S_inv, M_inv
elseif full
return S_inv, S, M_inv
else
return S_inv
end
```
"""
function square_root_inv(M::UniformScaling{T}; sq_rt::Bool=false, inverse::Bool=false,
full::Bool=false) where T <: Real
if sq_rt
S = M^0.5
S_inv = S^(-1.0)
S_inv, S
elseif inverse
M_inv = M^(-1.0)
S_inv = M_inv^0.5
S_inv, M_inv
elseif full
M_inv = M^(-1.0)
S = M^0.5
S_inv = S^(-1.0)
S_inv, S, M_inv
else
S_inv = M^(-0.5)
S_inv
end
end
function square_root_inv(M::Diagonal{T, Vector{T}}; sq_rt::Bool=false, inverse::Bool=false,
full::Bool=false) where T <: Real
if sq_rt
S = sqrt(M)
S_inv = inv(S)
S_inv, S
elseif inverse
M_inv = inv(M)
S_inv = sqrt(M_inv)
S_inv, M_inv
elseif full
S = sqrt(M)
S_inv = inv(S)
M_inv = inv(M)
S_inv, S, M
else
S_inv = M.^(-0.5)
S_inv
end
end
function square_root_inv(M::Symmetric{T, Matrix{T}}; sq_rt::Bool=false, inverse::Bool=false,
full::Bool=false) where T <: Real
# stable square root inverse for close-to-singular inverse calculations
F = svd(M)
if sq_rt
# take advantage of the SVD calculation to produce both the square root inverse
# and square root simultaneously
S_inv = Symmetric(F.U * Diagonal(1.0 ./ sqrt.(F.S)) * F.Vt)
S = Symmetric(F.U * Diagonal(sqrt.(F.S)) * F.Vt)
S_inv, S
elseif inverse
# take advantage of the SVD calculation to produce the square root inverse
# and inverse calculations all at once
S_inv = Symmetric(F.U * Diagonal(1.0 ./ sqrt.(F.S)) * F.Vt)
M_inv = Symmetric(F.U * Diagonal(1.0 ./ F.S) * F.Vt)
S_inv, M_inv
elseif full
# take advantage of the SVD calculation to produce the square root inverse,
# square root and inverse calculations all at once
S_inv = Symmetric(F.U * Diagonal(1.0 ./ sqrt.(F.S)) * F.Vt)
S = Symmetric(F.U * Diagonal(sqrt.(F.S)) * F.Vt)
M_inv = Symmetric(F.U * Diagonal(1.0 ./ F.S) * F.Vt)
S_inv, S, M_inv
else
# only return the square root inverse, if other calculations are not necessary
S_inv = Symmetric(F.U * Diagonal(1.0 ./ sqrt.(F.S)) * F.Vt)
S_inv
end
end
##############################################################################################
"""
transform_R(analysis::String, ens::ArView(T), obs::VecA(T), H_obs::Function,
obs_cov::CovM(T), kwargs::StepKwargs; conditioning::ConM=1000.0I,
m_err::ArView(T)=(1.0 ./ zeros(1,1)),
tol::Float64 = 0.0001,
j_max::Int64=40,
Q::CovM(T)=1.0I) where T <: Float64
Computes the ensemble transform and related values for various flavors of ensemble
Kalman schemes. The output type is a tuple containing a right transform of the ensemble
anomalies, the weights for the mean update and a random orthogonal transformation
for stabilization:
```
return trans, w, U
```
where the tuple is of type [`TransM`](@ref).
`m_err`, `tol`, `j_max`, `Q` are optional arguments depending on the `analysis`, with
default values provided.
Serves as an auxilliary function for EnKF, ETKF(-N), EnKS, ETKS(-N), where
"analysis" is a string which determines the type of transform update. The observation
error covariance `obs_cov` is of type [`CovM`](@ref), the conditioning matrix `conditioning`
is of type [`ConM`](@ref), the keyword arguments dictionary `kwargs` is of type
[`StepKwargs`](@ref) and the model error covariance matrix `Q` is of type [`CovM`](@ref).
Currently validated `analysis` options:
* `analysis=="etkf" || analysis=="etks"` computes the deterministic ensemble transform
as in the ETKF described in [Grudzien et al.
2022](https://gmd.copernicus.org/articles/15/7641/2022/gmd-15-7641-2022.html).
* `analysis[1:7]=="mlef-ls" || analysis[1:7]=="mles-ls"` computes the maximum likelihood
ensemble filter transform described in [Grudzien et al.
2022](https://gmd.copernicus.org/articles/15/7641/2022/gmd-15-7641-2022.html),
optimizing the nonlinear
cost function with Newton-based
[line searches](https://julianlsolvers.github.io/LineSearches.jl/stable/).
* `analysis[1:4]=="mlef" || analysis[1:4]=="mles"` computes the maximum likelihood
ensemble filter transform described in
[Grudzien et al. 2022](https://gmd.copernicus.org/articles/15/7641/2022/gmd-15-7641-2022.html),
optimizing the nonlinear
cost function with simple Newton-based scheme.
* `analysis=="enkf-n-dual" || analysis=="enks-n-dual"`
computes the dual form of the EnKF-N transform as in [Bocquet et al.
2015](https://npg.copernicus.org/articles/22/645/2015/)
Note: this cannot be used with the nonlinear observation operator.
This uses the Brent method for the argmin problem as this
has been more reliable at finding a global minimum than Newton optimization.
* `analysis=="enkf-n-primal" || analysis=="enks-n-primal"`
computes the primal form of the EnKF-N transform as in [Bocquet et al.
2015](https://npg.copernicus.org/articles/22/645/2015/),
[Grudzien et al. 2022](https://gmd.copernicus.org/articles/15/7641/2022/gmd-15-7641-2022.html).
This differs from the MLEF/S-N in that there is no approximate linearization of
the observation operator in the EnKF-N, this only handles the approximation error
with respect to the adaptive inflation. This uses a simple Newton-based
minimization of the cost function for the adaptive inflation.
* `analysis=="enkf-n-primal-ls" || analysis=="enks-n-primal-ls"`
computes the primal form of the EnKF-N transform as in [Bocquet et al.
2015](https://npg.copernicus.org/articles/22/645/2015/),
[Grudzien et al. 2022](https://gmd.copernicus.org/articles/15/7641/2022/gmd-15-7641-2022.html).
This differs from the MLEF/S-N in that there is no approximate linearization of
the observation operator in the EnKF-N, this only handles the approximation error
with respect to the adaptive inflation. This uses a Newton-based
minimization of the cost function for the adaptive inflation with
[line searches](https://julianlsolvers.github.io/LineSearches.jl/stable/).
"""
function transform_R(analysis::String, ens::ArView(T), obs::VecA(T), H_obs::Function,
obs_cov::CovM(T), kwargs::StepKwargs; conditioning::ConM(T)=1000.0I,
m_err::ArView(T)=(1.0 ./ zeros(1,1)),
tol::Float64 = 0.0001,
j_max::Int64=40,
Q::CovM(T)=1.0I) where T <: Float64
if analysis=="etkf" || analysis=="etks"
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: compute the ensemble in observation space
Y = H_obs(ens, obs_dim, kwargs)
# step 2: compute the ensemble mean in observation space
y_mean = mean(Y, dims=2)
# step 3: compute the sensitivity matrix in observation space
obs_sqrt_inv = square_root_inv(obs_cov)
S = obs_sqrt_inv * (Y .- y_mean )
# step 4: compute the weighted innovation
δ = obs_sqrt_inv * ( obs - y_mean )
# step 5: compute the approximate hessian
hessian = Symmetric((N_ens - 1.0)*I + transpose(S) * S)
# step 6: compute the transform matrix, transform matrix inverse and
# hessian inverse simultaneously via the SVD for stability
trans, hessian_inv = square_root_inv(hessian, inverse=true)
# step 7: compute the analysis weights
w = hessian_inv * transpose(S) * δ
# step 8: generate mean preserving random orthogonal matrix as in sakov oke 08
U = rand_orth(N_ens)
elseif analysis[1:7]=="mlef-ls" || analysis[1:7]=="mles-ls"
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: set up inputs for the optimization
# step 1a: inial choice is no change to the mean state
ens_mean_0 = mean(ens, dims=2)
anom_0 = ens .- ens_mean_0
w = zeros(N_ens)
# step 1b: pre-compute the observation error covariance square root
obs_sqrt_inv = square_root_inv(obs_cov)
# step 1c: define the conditioning and parameters for finite size formalism if needed
if analysis[end-5:end] == "bundle"
trans = inv(conditioning)
trans_inv = conditioning
elseif analysis[end-8:end] == "transform"
trans = Symmetric(Matrix(1.0*I, N_ens, N_ens))
trans_inv = Symmetric(Matrix(1.0*I, N_ens, N_ens))
end
if analysis[8:9] == "-n"
# define the epsilon scaling and the effective ensemble size if finite size form
ϵ_N = 1.0 + (1.0 / N_ens)
N_effective = N_ens + 1.0
end
# step 1d: define the storage of the gradient and Hessian as global to the functions
grad_w = Vector{Float64}(undef, N_ens)
hess_w = Array{Float64}(undef, N_ens, N_ens)
cost_w = 0.0
# step 2: define the cost / gradient / hessian function to avoid repeated computations
function fgh!(G, H, C, trans::ConM(T1), trans_inv::ConM(T1),
w::Vector{T1}) where T1 <: Float64
# step 2a: define the linearization of the observation operator
ens_mean_iter = ens_mean_0 + anom_0 * w
ens = ens_mean_iter .+ anom_0 * trans
Y = H_obs(ens, obs_dim, kwargs)
y_mean = mean(Y, dims=2)
# step 2b: compute the weighted anomalies in observation space, conditioned
# with trans inverse
S = obs_sqrt_inv * (Y .- y_mean) * trans_inv
# step 2c: compute the weighted innovation
δ = obs_sqrt_inv * (obs - y_mean)
# step 2d: gradient, hessian and cost function definitions
if G != nothing
if analysis[8:9] == "-n"
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
G[:] = N_effective * ζ * w - transpose(S) * δ
else
G[:] = (N_ens - 1.0) * w - transpose(S) * δ
end
end
if H != nothing
if analysis[8:9] == "-n"
H .= Symmetric((N_effective - 1.0)*I + transpose(S) * S)
else
H .= Symmetric((N_ens - 1.0)*I + transpose(S) * S)
end
end
if C != nothing
if analysis[8:9] == "-n"
y_mean_iter = H_obs(ens_mean_iter, obs_dim, kwargs)
δ = obs_sqrt_inv * (obs - y_mean_iter)
return N_effective * log(ϵ_N + sum(w.^2.0)) + sum(δ.^2.0)
else
y_mean_iter = H_obs(ens_mean_iter, obs_dim, kwargs)
δ = obs_sqrt_inv * (obs - y_mean_iter)
return (N_ens - 1.0) * sum(w.^2.0) + sum(δ.^2.0)
end
end
nothing
end
function newton_ls!(grad_w, hess_w, trans::ConM(T1), trans_inv::ConM(T1),
w::Vector{T1}, linesearch) where T1 <: Float64
# step 2e: find the Newton direction and the transform update if needed
fx = fgh!(grad_w, hess_w, cost_w, trans, trans_inv, w)
p = -hess_w \ grad_w
if analysis[end-8:end] == "transform"
trans_tmp, trans_inv_tmp = square_root_inv(Symmetric(hess_w), sq_rt=true)
trans .= trans_tmp
trans_inv .= trans_inv_tmp
end
# step 2f: univariate line search functions
ϕ(α) = fgh!(nothing, nothing, cost_w, trans, trans_inv, w .+ α.*p)
function dϕ(α)
fgh!(grad_w, nothing, nothing, trans, trans_inv, w .+ α.*p)
return dot(grad_w, p)
end
function ϕdϕ(α)
phi = fgh!(grad_w, nothing, cost_w, trans, trans_inv, w .+ α.*p)
dphi = dot(grad_w, p)
return (phi, dphi)
end
# step 2g: define the linesearch
dϕ_0 = dot(p, grad_w)
α, fx = linesearch(ϕ, dϕ, ϕdϕ, 1.0, fx, dϕ_0)
Δw = α * p
w .= w + Δw
return Δw
end
# step 3: optimize
# step 3a: perform the optimization by Newton with linesearch
# we use StrongWolfe for RMSE performance as the default linesearch
#ln_search = HagerZhang()
ln_search = StrongWolfe()
j = 0
Δw = ones(N_ens)
while j < j_max && norm(Δw) > tol
Δw = newton_ls!(grad_w, hess_w, trans, trans_inv, w, ln_search)
end
if analysis[8:9] == "-n"
# peform a final inflation with the finite size cost function
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
hess_w = ζ * I - 2.0 * ζ^2.0 * w * transpose(w)
hess_w = Symmetric(transpose(S) * S + (N_ens + 1.0) * hess_w)
trans = square_root_inv(hess_w)
elseif analysis[end-5:end] == "bundle"
trans = square_root_inv(hess_w)
end
# step 3b: generate mean preserving random orthogonal matrix as in sakov oke 08
U = rand_orth(N_ens)
elseif analysis[1:4]=="mlef" || analysis[1:4]=="mles"
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: set up the optimization, inial choice is no change to the mean state
ens_mean_0 = mean(ens, dims=2)
anom_0 = ens .- ens_mean_0
w = zeros(N_ens)
# pre-compute the observation error covariance square root
obs_sqrt_inv = square_root_inv(obs_cov)
# define these variables as global compared to the while loop
grad_w = Vector{Float64}(undef, N_ens)
hess_w = Array{Float64}(undef, N_ens, N_ens)
S = Array{Float64}(undef, obs_dim, N_ens)
ens_mean_iter = copy(ens_mean_0)
# define the conditioning
if analysis[end-5:end] == "bundle"
trans = inv(conditioning)
trans_inv = conditioning
elseif analysis[end-8:end] == "transform"
trans = 1.0*I
trans_inv = 1.0*I
end
# step 2: perform the optimization by simple Newton
j = 0
if analysis[5:6] == "-n"
# define the epsilon scaling and the effective ensemble size if finite size form
ϵ_N = 1.0 + (1.0 / N_ens)
N_effective = N_ens + 1.0
end
while j < j_max
# step 2a: compute the observed ensemble and ensemble mean
ens_mean_iter = ens_mean_0 + anom_0 * w
ens = ens_mean_iter .+ anom_0 * trans
Y = H_obs(ens, obs_dim, kwargs)
y_mean = mean(Y, dims=2)
# step 2b: compute the weighted anomalies in observation space, conditioned
# with trans inverse
S = obs_sqrt_inv * (Y .- y_mean) * trans_inv
# step 2c: compute the weighted innovation
δ = obs_sqrt_inv * (obs - y_mean)
# step 2d: compute the gradient and hessian
if analysis[5:6] == "-n"
# for finite formalism, we follow the IEnKS-N convention where
# the gradient is computed with the finite-size cost function but we use the
# usual hessian, with the effective ensemble size
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
grad_w = N_effective * ζ * w - transpose(S) * δ
hess_w = Symmetric((N_effective - 1.0)*I + transpose(S) * S)
else
grad_w = (N_ens - 1.0) * w - transpose(S) * δ
hess_w = Symmetric((N_ens - 1.0)*I + transpose(S) * S)
end
# step 2e: perform Newton approximation, simultaneously computing
# the update transform trans with the SVD based inverse at once
if analysis[end-8:end] == "transform"
trans, trans_inv, hessian_inv = square_root_inv(hess_w, full=true)
Δw = hessian_inv * grad_w
else
Δw = hess_w \ grad_w
end
# 2f: update the weights
w -= Δw
if norm(Δw) < tol
break
else
# step 2g: update the iterative mean state
j+=1
end
end
if analysis[5:6] == "-n"
# peform a final inflation with the finite size cost function
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
hess_w = ζ * I - 2.0 * ζ^2.0 * w * transpose(w)
hess_w = Symmetric(transpose(S) * S + (N_ens + 1.0) * hess_w)
trans = square_root_inv(hess_w)
elseif analysis[end-5:end] == "bundle"
trans = square_root_inv(hess_w)
end
# step 7: generate mean preserving random orthogonal matrix as in sakov oke 08
U = rand_orth(N_ens)
elseif analysis=="etkf-sqrt-core" || analysis=="etks-sqrt-core"
### NOTE: STILL DEVELOPMENT CODE, NOT DEBUGGED
# needs to be revised for the calculation with unweighted anomalies
# Uses the contribution of the model error covariance matrix Q
# in the square root as in Raanes et al. 2015
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: compute the ensemble mean
x_mean = mean(ens, dims=2)
# step 2a: compute the normalized anomalies
A = (ens .- x_mean) / sqrt(N_ens - 1.0)
# step 2b: compute the SVD for the two-sided projected model error covariance
F = svd(A)
Σ_inv = Diagonal([1.0 ./ F.S[1:N_ens-1]; 0.0])
p_inv = F.V * Σ_inv * transpose(F.U)
## NOTE: want to
G = Symmetric(1.0I + (N_ens - 1.0) * p_inv * Q * transpose(p_inv))
# step 2c: compute the model error adjusted anomalies
A = A * square_root(G)
# step 3: compute the ensemble in observation space
Y = H_obs(ens, obs_dim, kwargs)
# step 4: compute the ensemble mean in observation space
y_mean = mean(Y, dims=2)
# step 5: compute the weighted anomalies in observation space
# first we find the observation error covariance inverse
obs_sqrt_inv = square_root_inv(obs_cov)
# then compute the weighted anomalies
S = (Y .- y_mean) / sqrt(N_ens - 1.0)
S = obs_sqrt_inv * S
# step 6: compute the weighted innovation
δ = obs_sqrt_inv * ( obs - y_mean )
# step 7: compute the transform matrix
trans = inv(Symmetric(1.0I + transpose(S) * S))
# step 8: compute the analysis weights
w = trans * transpose(S) * δ
# step 9: compute the square root of the transform
trans = sqrt(trans)
# step 10: generate mean preserving random orthogonal matrix as in sakov oke 08
U = rand_orth(N_ens)
elseif analysis=="enkf-n-dual" || analysis=="enks-n-dual"
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: compute the observed ensemble and ensemble mean
Y = H_obs(ens, obs_dim, kwargs)
y_mean = mean(Y, dims=2)
# step 2: compute the weighted anomalies in observation space
# first we find the observation error covariance inverse
obs_sqrt_inv = square_root_inv(obs_cov)
# then compute the sensitivity matrix in observation space
S = obs_sqrt_inv * (Y .- y_mean)
# step 5: compute the weighted innovation
δ = obs_sqrt_inv * (obs - y_mean)
# step 6: compute the SVD for the simplified cost function, gauge weights and range
F = svd(S)
ϵ_N = 1.0 + (1.0 / N_ens)
ζ_l = 0.000001
ζ_u = (N_ens + 1.0) / ϵ_N
# step 7: define the dual cost function derived in singular value form
function D(ζ::Float64)
cost = I - (F.U * Diagonal( F.S.^2.0 ./ (ζ .+ F.S.^2.0) ) * transpose(F.U) )
cost = transpose(δ) * cost * δ .+ ϵ_N * ζ .+
(N_ens + 1.0) * log((N_ens + 1.0) / ζ) .- (N_ens + 1.0)
cost[1]
end
# The below is defined for possible Hessian-based minimization
# NOTE: standard Brent's method appears to be more reliable at finding a
# global minimizer with some basic tests, may be tested further
#
#function D_v(ζ::Vector{Float64})
# ζ = ζ[1]
# cost = I - (F.U * Diagonal( F.S.^2.0 ./ (ζ .+ F.S.^2.0) ) * transpose(F.U) )
# cost = transpose(δ) * cost * δ .+ ϵ_N * ζ .+
# (N_ens + 1.0) * log((N_ens + 1.0) / ζ) .- (N_ens + 1.0)
# cost[1]
#end
#function D_prime!(storage::Vector{Float64}, ζ::Vector{Float64})
# ζ = ζ[1]
# grad = transpose(δ) * F.U * Diagonal( - F.S.^2.0 .* (ζ .+ F.S.^2.0).^(-2.0) ) *
# transpose(F.U) * δ
# storage[:, :] = grad .+ ϵ_N .- (N_ens + 1.0) / ζ
#end
#function D_hess!(storage::Array{Float64}, ζ::Vector{Float64})
# ζ = ζ[1]
# hess = transpose(δ) * F.U *
# Diagonal( 2.0 * F.S.^2.0 .* (ζ .+ F.S.^2.0).^(-3.0) ) * transpose(F.U) * δ
# storage[:, :] = hess .+ (N_ens + 1.0) * ζ^(-2.0)
#end
#lx = [ζ_l]
#ux = [ζ_u]
#ζ_0 = [(ζ_u + ζ_l)/2.0]
#df = TwiceDifferentiable(D_v, D_prime!, D_hess!, ζ_0)
#dfc = TwiceDifferentiableConstraints(lx, ux)
#ζ_b = optimize(D_v, D_prime!, D_hess!, ζ_0)
# step 8: find the argmin
ζ_a = optimize(D, ζ_l, ζ_u)
diag_vals = ζ_a.minimizer .+ F.S.^2.0
# step 9: compute the update weights
w = F.V * Diagonal( F.S ./ diag_vals ) * transpose(F.U) * δ
# step 10: compute the update transform
trans = Symmetric(Diagonal( F.S ./ diag_vals) * transpose(F.U) * δ *
transpose(δ) * F.U * Diagonal( F.S ./ diag_vals))
trans = Symmetric(Diagonal(diag_vals) -
( (2.0 * ζ_a.minimizer^2.0) / (N_ens + 1.0) ) * trans)
trans = Symmetric(F.V * square_root_inv(trans) * F.Vt)
# step 11: generate mean preserving random orthogonal matrix as in sakov oke 08
U = rand_orth(N_ens)
elseif analysis=="enkf-n-primal" || analysis=="enks-n-primal"
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: compute the observed ensemble and ensemble mean
Y = H_obs(ens, obs_dim, kwargs)
y_mean = mean(Y, dims=2)
# step 2: compute the weighted anomalies in observation space
# first we find the observation error covariance inverse
obs_sqrt_inv = square_root_inv(obs_cov)
# then compute the sensitivity matrix in observation space
S = obs_sqrt_inv * (Y .- y_mean)
# step 3: compute the weighted innovation
δ = obs_sqrt_inv * (obs - y_mean)
# step 4: define the epsilon scaling and the effective ensemble size
ϵ_N = 1.0 + (1.0 / N_ens)
N_effective = N_ens + 1.0
# step 5: set up the optimization
# step 5:a the inial choice is no change to the mean state
w = zeros(N_ens)
# step 5b: define the primal cost function
function P(w::Vector{Float64})
cost = (δ - S * w)
cost = sum(cost.^2.0) + N_effective * log(ϵ_N + sum(w.^2.0))
0.5 * cost
end
# step 5c: define the primal gradient
function ∇P!(grad::Vector{Float64}, w::Vector{Float64})
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
grad[:] = N_effective * ζ * w - transpose(S) * (δ - S * w)
end
# step 5d: define the primal hessian
function H_P!(hess::ArView(T1), w::Vector{Float64}) where T1 <: Float64
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
hess .= ζ * I - 2.0 * ζ^2.0 * w * transpose(w)
hess .= transpose(S) * S + N_effective * hess
end
# step 6: perform the optimization by simple Newton
j = 0
trans = Array{Float64}(undef, N_ens, N_ens)
grad_w = Array{Float64}(undef, N_ens)
hess_w = Array{Float64}(undef, N_ens, N_ens)
while j < j_max
# compute the gradient and hessian
∇P!(grad_w, w)
H_P!(hess_w, w)
# perform Newton approximation, simultaneously computing
# the update transform trans with the SVD based inverse at once
trans, hessian_inv = square_root_inv(Symmetric(hess_w), inverse=true)
Δw = hessian_inv * grad_w
w -= Δw
if norm(Δw) < tol
break
else
j+=1
end
end
# step 7: generate mean preserving random orthogonal matrix as in sakov oke 08
U = rand_orth(N_ens)
elseif analysis=="enkf-n-primal-ls" || analysis=="enks-n-primal-ls"
# step 0: infer the system, observation and ensemble dimensions
sys_dim, N_ens = size(ens)
obs_dim = length(obs)
# step 1: compute the observed ensemble and ensemble mean
Y = H_obs(ens, obs_dim, kwargs)
y_mean = mean(Y, dims=2)
# step 2: compute the weighted anomalies in observation space
# first we find the observation error covariance inverse
obs_sqrt_inv = square_root_inv(obs_cov)
# then compute the sensitivity matrix in observation space
S = obs_sqrt_inv * (Y .- y_mean)
# step 3: compute the weighted innovation
δ = obs_sqrt_inv * (obs - y_mean)
# step 4: define the epsilon scaling and the effective ensemble size
ϵ_N = 1.0 + (1.0 / N_ens)
N_effective = N_ens + 1.0
# step 5: set up the optimization
# step 5:a the inial choice is no change to the mean state
w = zeros(N_ens)
# step 5b: define the primal cost function
function J(w::Vector{Float64})
cost = (δ - S * w)
cost = sum(cost.^2.0) + N_effective * log(ϵ_N + sum(w.^2.0))
0.5 * cost
end
# step 5c: define the primal gradient
function ∇J!(grad::Vector{Float64}, w::Vector{Float64})
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
grad[:] = N_effective * ζ * w - transpose(S) * (δ - S * w)
end
# step 5d: define the primal hessian
function H_J!(hess::ArView(T1), w::Vector{Float64}) where T1 <: Float64
ζ = 1.0 / (ϵ_N + sum(w.^2.0))
hess .= ζ * I - 2.0 * ζ^2.0 * w * transpose(w)
hess .= transpose(S) * S + N_effective * hess
end
# step 6: find the argmin for the update weights
# step 6a: define the line search algorithm with Newton
# we use StrongWolfe for RMSE performance as the default linesearch
# method, see the LineSearches docs, alternative choice is commented below
# ln_search = HagerZhang()
ln_search = StrongWolfe()
opt_alg = Newton(linesearch = ln_search)
# step 6b: perform the optimization
w = Optim.optimize(J, ∇J!, H_J!, w, method=opt_alg, x_tol=tol).minimizer
# step 7: compute the update transform
trans = Symmetric(H_J!(Array{Float64}(undef, N_ens, N_ens), w))
trans = square_root_inv(trans)
# step 8: generate mean preserving random orthogonal matrix as in Sakov & Oke 08
U = rand_orth(N_ens)
end
return trans, w, U
end
##############################################################################################
"""
ens_gauss_newton(analysis::String, ens::ArView(T), obs::VecA(T),
H_obs::Function, obs_cov::CovM(T), kwargs::StepKwargs;
conditioning::ConM(T)=1000.0I,
m_err::ArView(T)=(1.0 ./ zeros(1,1)),
tol::Float64 = 0.0001,
j_max::Int64=40,
Q::CovM(T)=1.0I) where T <: Float64
Computes the ensemble estimated gradient and Hessian terms for nonlinear least-squares
```
return ∇_J, Hess_J
```
`m_err`, `tol`, `j_max`, `Q` are optional arguments depending on the `analysis`, with
default values provided.
Serves as an auxilliary function for IEnKS(-N), where "analysis" is a string which
determines the method of transform update ensemble Gauss-Newton calculation. The observation
error covariance `obs_cov` is of type [`CovM`](@ref), the conditioning matrix `conditioning`
is of type [`ConM`](@ref), the keyword arguments dictionary `kwargs` is of type
[`StepKwargs`](@ref) and the model error covariance matrix `Q` is of type [`CovM`](@ref).
Currently validated `analysis` options:
* `analysis == "ienks-bundle" || "ienks-n-bundle" || "ienks-transform" || "ienks-n-transform"`
computes the weighted observed anomalies as per the
bundle or transform version of the IEnKS, described in [Bocquet et al.
2013](https://rmets.onlinelibrary.wiley.com/doi/abs/10.1002/qj.2236),
[Grudzien et al. 2022](https://gmd.copernicus.org/articles/15/7641/2022/gmd-15-7641-2022.html).
Bundle versus tranform versions of the scheme are specified by the trailing
`analysis` string as `-bundle` or `-transform`. The bundle version uses a small uniform
scalar `ϵ`, whereas the transform version uses a matrix square root inverse as the
conditioning operator. This form of analysis differs from other schemes by returning a
sequential-in-time value for the cost function gradient and Hessian, which will is
utilized within the iterative smoother optimization. A finite-size inflation scheme,
based on the EnKF-N above, can be utilized by appending additionally a `-n` to the
`-bundle` or `-transform` version of the IEnKS scheme specified in `analysis`.
"""
function ens_gauss_newton(analysis::String, ens::ArView(T), obs::VecA(T),
H_obs::Function, obs_cov::CovM(T), kwargs::StepKwargs;
conditioning::ConM(T)=1000.0I,
m_err::ArView(T)=(1.0 ./ zeros(1,1)),
tol::Float64 = 0.0001,
j_max::Int64=40,
Q::CovM(T)=1.0I) where T <: Float64
if analysis[1:5]=="ienks"
# step 0: infer observation dimension
obs_dim = length(obs)
# step 1: compute the observed ensemble and ensemble mean
Y = H_obs(ens, obs_dim, kwargs)
y_mean = mean(Y, dims=2)
# step 2: compute the observed anomalies, proportional to the conditioning matrix
# here conditioning should be supplied as trans^(-1)
S = (Y .- y_mean) * conditioning
# step 3: compute the cost function gradient term
inv_obs_cov = inv(obs_cov)
∇J = transpose(S) * inv_obs_cov * (obs - y_mean)