/
MadsSensitivityAnalysis.jl
1592 lines (1487 loc) · 60.8 KB
/
MadsSensitivityAnalysis.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
import MetaProgTools
import ProgressMeter
import Distributions
import OrderedCollections
import DataFrames
import StatsBase
import JSON
import JLD2
import FileIO
import Random
import Distributed
import LinearAlgebra
"""
Make gradient function needed for local sensitivity analysis
$(DocumentFunction.documentfunction(makelocalsafunction;
argtext=Dict("madsdata"=>"MADS problem dictionary"),
keytext=Dict("multiplycenterbyweights"=>"multiply center by observation weights [default=`true`]")))
Returns:
- gradient function
"""
function makelocalsafunction(madsdata::AbstractDict; multiplycenterbyweights::Bool=true)
f = makemadscommandfunction(madsdata)
restartdir = getrestartdir(madsdata)
obskeys = Mads.getobskeys(madsdata)
weights = Mads.getobsweight(madsdata, obskeys)
nO = length(obskeys)
optparamkeys = Mads.getoptparamkeys(madsdata)
lineardx = getparamsstep(madsdata, optparamkeys)
nP = length(optparamkeys)
initparams = Mads.getparamdict(madsdata)
function f_sa(arrayparameters::AbstractVector)
parameters = copy(initparams)
for i = 1:length(arrayparameters)
parameters[optparamkeys[i]] = arrayparameters[i]
end
resultdict = f(parameters)
results = Array{Float64}(undef, 0)
for obskey in obskeys
push!(results, resultdict[obskey]) # preserve the expected order
end
return results .* weights
end
function inner_grad(arrayparameters_dx_center_tuple::Tuple)
arrayparameters = arrayparameters_dx_center_tuple[1]
dx = arrayparameters_dx_center_tuple[2]
center = arrayparameters_dx_center_tuple[3]
if sizeof(center) > 0 && multiplycenterbyweights
center = center .* weights
end
if sizeof(dx) == 0
dx = lineardx
end
p = Vector{Float64}[]
for i in 1:nP
a = copy(arrayparameters)
a[i] += dx[i]
push!(p, a)
end
if sizeof(center) == 0
filename = ReusableFunctions.gethashfilename(restartdir, arrayparameters)
center = ReusableFunctions.loadresultfile(filename)
center_computed = (center !== nothing) && length(center) == nO
if !center_computed
push!(p, arrayparameters)
end
else
center_computed = true
end
local fevals
try
fevals = RobustPmap.rpmap(f_sa, p)
catch errmsg
printerrormsg(errmsg)
Mads.madswarn("RobustPmap executions for localsa fails!")
return nothing
end
if !center_computed
center = fevals[nP+1]
if restartdir != ""
ReusableFunctions.saveresultfile(restartdir, center, arrayparameters)
end
end
jacobian = Array{Float64}(undef, nO, nP)
for j in 1:nO
for i in 1:nP
jacobian[j, i] = (fevals[i][j] - center[j]) / dx[i]
end
end
return jacobian
end
reusable_inner_grad = makemadsreusablefunction(madsdata, inner_grad, "g_lm"; usedict=false)
"""
Gradient function for the forward model used for local sensitivity analysis
"""
function g_sa(arrayparameters::AbstractVector{Float64}; dx::Array{Float64,1}=Array{Float64}(undef, 0), center::Array{Float64,1}=Array{Float64}(undef, 0))
return reusable_inner_grad(tuple(arrayparameters, dx, center))
end
return f_sa, g_sa
end
"""
Local sensitivity analysis based on eigen analysis of the parameter covariance matrix
$(DocumentFunction.documentfunction(localsa;
argtext=Dict("madsdata"=>"MADS problem dictionary"),
keytext=Dict("sinspace"=>"apply sin transformation [default=`true`]",
"keyword"=>"keyword to be added in the filename root",
"filename"=>"output file name",
"format"=>"output plot format (`png`, `pdf`, etc.)",
"datafiles"=>"flag to write data files [default=`true`]",
"imagefiles"=>"flag to create image files [default=`Mads.graphoutput`]",
"par"=>"parameter set",
"obs"=>"observations for the parameter set",
"J"=>"Jacobian matrix")))
Dumps:
- `filename` : output plot file
"""
function localsa(madsdata::AbstractDict; sinspace::Bool=true, keyword::AbstractString="", filename::AbstractString="", format::AbstractString="", datafiles::Bool=true, imagefiles::Bool=Mads.graphoutput, par::AbstractVector=Array{Float64}(undef, 0), obs::AbstractVector=Array{Float64}(undef, 0), J::AbstractMatrix=Array{Float64}(undef, 0, 0))
f_sa, g_sa = Mads.makelocalsafunction(madsdata)
if haskey(ENV, "MADS_NO_PLOT") || haskey(ENV, "MADS_NO_GADFLY") || !isdefined(Mads, :Gadfly)
imagefiles = false
end
if filename == ""
rootname = Mads.getmadsrootname(madsdata)
ext = ""
else
rootname = Mads.getrootname(filename)
ext = "." * Mads.getextension(filename)
end
if keyword != ""
rootname = string(rootname, "-", keyword)
end
paramkeys = getoptparamkeys(madsdata)
nPall = length(getparamkeys(madsdata))
obskeys = getobskeys(madsdata)
plotlabels = getparamlabels(madsdata, paramkeys)
nP = length(paramkeys)
nPi = length(par)
if nPi == 0
param = getparamsinit(madsdata, paramkeys)
elseif nPi == nP
param = par
elseif nPi == nPall
param = Mads.getoptparams(madsdata, par)
else
param = getoptparams(madsdata, par, paramkeys)
end
nO = length(obskeys)
if sizeof(J) == 0
if sinspace
lowerbounds = Mads.getparamsmin(madsdata, paramkeys)
upperbounds = Mads.getparamsmax(madsdata, paramkeys)
logtransformed = Mads.getparamslog(madsdata, paramkeys)
indexlogtransformed = findall(logtransformed)
lowerbounds[indexlogtransformed] = log10.(lowerbounds[indexlogtransformed])
upperbounds[indexlogtransformed] = log10.(upperbounds[indexlogtransformed])
sinparam = asinetransform(param, lowerbounds, upperbounds, indexlogtransformed)
sindx = Mads.getsindx(madsdata)
g_sa_sin = Mads.sinetransformgradient(g_sa, lowerbounds, upperbounds, indexlogtransformed, sindx=sindx)
J = g_sa_sin(sinparam, center=obs)
else
J = g_sa(param, center=obs)
end
end
if J === nothing
Mads.madswarn("Jacobian computation failed")
return
end
if any(isnan, J)
Mads.madswarn("Local sensitivity analysis cannot be performed; provided Jacobian matrix contains NaN's")
Base.display(J)
Mads.madscritical("Mads quits!")
end
if length(obskeys) != size(J, 1) && length(paramkeys) != size(J, 2)
Mads.madscritical("Jacobian matrix size does not match the problem: J $(size(J))")
end
f = Mads.forward(madsdata, param)
ofval = Mads.of(madsdata, f)
datafiles && DelimitedFiles.writedlm("$(rootname)-jacobian.dat", [transposevector(["Obs"; paramkeys]); obskeys J])
mscale = max(abs(minimumnan(J)), abs(maximumnan(J)))
if imagefiles
jacmat = Gadfly.spy(J, Gadfly.Scale.x_discrete(labels = i->plotlabels[i]), Gadfly.Scale.y_discrete(labels = i->obskeys[i]),
Gadfly.Guide.YLabel("Observations"), Gadfly.Guide.XLabel("Parameters"),
Gadfly.Theme(point_size=20Gadfly.pt, major_label_font_size=14Gadfly.pt, minor_label_font_size=12Gadfly.pt, key_title_font_size=16Gadfly.pt, key_label_font_size=12Gadfly.pt),
Gadfly.Scale.ContinuousColorScale(Gadfly.Scale.lab_gradient(Base.parse(Colors.Colorant, "green"), Base.parse(Colors.Colorant, "yellow"), Base.parse(Colors.Colorant, "red")), minvalue = -mscale, maxvalue = mscale))
filename = "$(rootname)-jacobian" * ext
filename, format = setplotfileformat(filename, format)
try
Gadfly.draw(Gadfly.eval(Symbol(format))(filename, 3Gadfly.inch+0.25Gadfly.inch*nP, 3Gadfly.inch+0.25Gadfly.inch*nO), jacmat)
catch
madswarn("Gadfly could not plot!")
end
Mads.madsinfo("Jacobian matrix plot saved in $filename")
end
JpJ = J' * J
local covar
try
u, s, v = LinearAlgebra.svd(JpJ)
covar = v * inv(LinearAlgebra.Diagonal(s)) * u'
catch errmsg1
try
covar = inv(JpJ)
catch errmsg2
printerrormsg(errmsg1)
printerrormsg(errmsg2)
Mads.madswarn("JpJ inversion fails")
return nothing
end
end
stddev = sqrt.(abs.(LinearAlgebra.diag(covar)))
if datafiles
DelimitedFiles.writedlm("$(rootname)-covariance.dat", covar)
f = open("$(rootname)-stddev.dat", "w")
for i in 1:nP
write(f, "$(paramkeys[i]) $(param[i]) $(stddev[i])\n")
end
close(f)
end
correl = covar ./ LinearAlgebra.diag(covar)
datafiles && DelimitedFiles.writedlm("$(rootname)-correlation.dat", correl)
z = LinearAlgebra.eigen(covar)
eigenv = z.values
eigenm = z.vectors
eigenv = abs.(eigenv)
index = sortperm(eigenv)
sortedeigenv = eigenv[index]
sortedeigenm = real(eigenm[:,index])
datafiles && DelimitedFiles.writedlm("$(rootname)-eigenmatrix.dat", [paramkeys sortedeigenm])
datafiles && DelimitedFiles.writedlm("$(rootname)-eigenvalues.dat", sortedeigenv)
if imagefiles
eigenmat = Gadfly.spy(sortedeigenm, Gadfly.Scale.y_discrete(labels = i->plotlabels[i]), Gadfly.Scale.x_discrete,
Gadfly.Guide.YLabel("Parameters"), Gadfly.Guide.XLabel("Eigenvectors"),
Gadfly.Theme(point_size=20Gadfly.pt, major_label_font_size=14Gadfly.pt, minor_label_font_size=12Gadfly.pt, key_title_font_size=16Gadfly.pt, key_label_font_size=12Gadfly.pt),
Gadfly.Scale.ContinuousColorScale(Gadfly.Scale.lab_gradient(Base.parse(Colors.Colorant, "green"), Base.parse(Colors.Colorant, "yellow"), Base.parse(Colors.Colorant, "red"))))
# eigenval = plot(x=1:length(sortedeigenv), y=sortedeigenv, Scale.x_discrete, Scale.y_log10, Geom.bar, Guide.YLabel("Eigenvalues"), Guide.XLabel("Eigenvectors"))
filename = "$(rootname)-eigenmatrix" * ext
filename, format = setplotfileformat(filename, format)
Gadfly.draw(Gadfly.eval(Symbol(format))(filename, 4Gadfly.inch+0.25Gadfly.inch*nP, 4Gadfly.inch+0.25Gadfly.inch*nP), eigenmat)
Mads.madsinfo("Eigen matrix plot saved in $filename")
eigenval = Gadfly.plot(x=1:length(sortedeigenv), y=sortedeigenv, Gadfly.Scale.x_discrete, Gadfly.Scale.y_log10,
Gadfly.Geom.line(),
Gadfly.Theme(line_width=4Gadfly.pt, major_label_font_size=14Gadfly.pt, minor_label_font_size=12Gadfly.pt, key_title_font_size=16Gadfly.pt, key_label_font_size=12Gadfly.pt),
Gadfly.Guide.YLabel("Eigenvalues"), Gadfly.Guide.XLabel("Eigenvectors"))
filename = "$(rootname)-eigenvalues" * ext
filename, format = setplotfileformat(filename, format)
Gadfly.draw(Gadfly.eval(Symbol(format))(filename, 4Gadfly.inch+0.25Gadfly.inch*nP, 4Gadfly.inch), eigenval)
Mads.madsinfo("Eigen values plot saved in $filename")
end
Dict("of"=>ofval, "jacobian"=>J, "covar"=>covar, "stddev"=>stddev, "eigenmatrix"=>sortedeigenm, "eigenvalues"=>sortedeigenv)
end
"""
$(DocumentFunction.documentfunction(sampling;
argtext=Dict("param"=>"Parameter vector",
"J"=>"Jacobian matrix",
"numsamples"=>"Number of samples"),
keytext=Dict("seed"=>"random esee [default=`0`]",
"scale"=>"data scaling [default=`1`]")))
Returns:
- generated samples (vector or array)
- vector of log-likelihoods
"""
function sampling(param::AbstractVector, J::Array, numsamples::Number; seed::Integer=-1, scale::Number=1)
u, d, v = LinearAlgebra.svd(J' * J)
done = false
vo = copy(v)
local gooddirections
local dist
numdirections = length(d)
numgooddirections = numdirections
while !done
try
covmat = (v * LinearAlgebra.Diagonal(1 ./ d) * permutedims(u)) .* scale
dist = Distributions.MvNormal(zeros(numgooddirections), covmat)
done = true
catch errmsg
# printerrormsg(errmsg)
numgooddirections -= 1
if numgooddirections <= 0
madscritical("Reduction in sampling directions failed!")
end
gooddirections = vo[:, 1:numgooddirections]
newJ = J * gooddirections
u, d, v = LinearAlgebra.svd(newJ' * newJ)
end
end
madsinfo("Reduction in sampling directions ... (from $(numdirections) to $(numgooddirections))")
setseed(seed)
gooddsamples = Distributions.rand(dist, numsamples)
llhoods = map(i->Distributions.loglikelihood(dist, gooddsamples[:, i:i]), 1:numsamples)
if numdirections > numgooddirections
samples = gooddirections * gooddsamples
else
samples = gooddsamples
end
for i = 1:size(samples, 2)
samples[:, i] += param
end
return samples, llhoods
end
"""
Reweigh samples using importance sampling -- returns a vector of log-likelihoods after reweighing
$(DocumentFunction.documentfunction(reweighsamples;
argtext=Dict("madsdata"=>"MADS problem dictionary",
"predictions"=>"the model predictions for each of the samples",
"oldllhoods"=>"the log likelihoods of the parameters in the old distribution")))
Returns:
- vector of log-likelihoods after reweighing
"""
function reweighsamples(madsdata::AbstractDict, predictions::Array, oldllhoods::AbstractVector)
obskeys = getobskeys(madsdata)
weights = getobsweight(madsdata)
targets = getobstarget(madsdata)
newllhoods = -oldllhoods
j = 1
for okey in obskeys
if haskey(madsdata["Observations"][okey], "weight")
newllhoods -= .5 * (weights[j] * (predictions[j, :] .- targets[j])) .^ 2
end
j += 1
end
return newllhoods .- maximumnan(newllhoods) # normalize likelihoods so the most likely thing has likelihood 1
end
"""
Get important samples
$(DocumentFunction.documentfunction(getimportantsamples;
argtext=Dict("samples"=>"array of samples",
"llhoods"=>"vector of log-likelihoods")))
Returns:
- array of important samples
"""
function getimportantsamples(samples::Array, llhoods::AbstractVector)
sortedlhoods = sort(exp.(llhoods), rev=true)
sortedprobs = sortedlhoods / sum(sortedlhoods)
cumprob = 0.
i = 1
while cumprob < .95
cumprob += sortedprobs[i]
i += 1
end
thresholdllhood = log(sortedlhoods[i - 1])
imp_samples = Array{Float64}(undef, size(samples, 2), 0)
for i = 1:length(llhoods)
if llhoods[i] > thresholdllhood
imp_samples = hcat(imp_samples, vec(samples[i, :]))
end
end
return imp_samples
end
"""
Get weighted mean and variance samples
$(DocumentFunction.documentfunction(weightedstats;
argtext=Dict("samples"=>"array of samples",
"llhoods"=>"vector of log-likelihoods")))
Returns:
- vector of sample means
- vector of sample variances
"""
function weightedstats(samples::Array, llhoods::AbstractVector)
wv = StatsBase.Weights(exp.(llhoods))
return Statistics.mean(samples, wv, 1), Statistics.var(samples, wv, 1; corrected=false)
end
function getparamrandom(madsdata::AbstractDict, numsamples::Integer=1, parameterkey::AbstractString=""; init_dist::Bool=false)
if parameterkey != ""
return getparamrandom(madsdata, parameterkey; numsamples=numsamples, init_dist=init_dist)
else
sample = OrderedCollections.OrderedDict()
paramkeys = getoptparamkeys(madsdata)
paramdist = getparamdistributions(madsdata; init_dist=init_dist)
for k in paramkeys
if numsamples == 1 # if only one sample, don't write a 1-element array to each dictionary key. write a scalar.
sample[k] = getparamrandom(madsdata, k; numsamples=numsamples, paramdist=paramdist)[1]
else
sample[k] = getparamrandom(madsdata, k; numsamples=numsamples, paramdist=paramdist)
end
end
return sample
end
end
function getparamrandom(madsdata::AbstractDict, parameterkey::AbstractString; numsamples::Integer=1, paramdist::AbstractDict=Dict(), init_dist::Bool=false)
if haskey(madsdata, "Parameters") && haskey(madsdata["Parameters"], parameterkey)
if length(paramdist) == 0
paramdist = getparamdistributions(madsdata; init_dist=init_dist)
end
if Mads.islog(madsdata, parameterkey)
dist = paramdist[parameterkey]
if typeof(dist) <: Distributions.Uniform
a = log10(dist.a)
b = log10(dist.b)
return 10. .^(a .+ (b .- a) .* Distributions.rand(numsamples))
elseif typeof(dist) <: Distributions.Normal
μ = log10(dist.μ)
return 10. .^(μ .+ dist.σ .* Distributions.randn(numsamples))
end
end
return Distributions.rand(paramdist[parameterkey], numsamples)
end
return nothing
end
@doc """
Get independent sampling of model parameters defined in the MADS problem dictionary
$(DocumentFunction.documentfunction(getparamrandom;
argtext=Dict("madsdata"=>"MADS problem dictionary",
"parameterkey"=>"model parameter key",
"numsamples"=>"number of samples, [default=`1`]"),
keytext=Dict("init_dist"=>"if `true` use the distribution set for initialization in the MADS problem dictionary (defined using `init_dist` parameter field); if `false` (default) use the regular distribution set in the MADS problem dictionary (defined using `dist` parameter field)",
"numsamples"=>"number of samples",
"paramdist"=>"dictionary of parameter distributions")))
Returns:
- generated sample
""" getparamrandom
"""
Saltelli sensitivity analysis (brute force)
$(DocumentFunction.documentfunction(saltellibrute;
argtext=Dict("madsdata"=>"MADS problem dictionary"),
keytext=Dict("N"=>"number of samples [default=`1000`]",
"seed"=>"random seed [default=`0`]",
"restartdir"=>"directory where files will be stored containing model results for fast simulation restarts")))
"""
function saltellibrute(madsdata::AbstractDict; N::Integer=1000, seed::Integer=-1, restartdir::AbstractString="") # TODO Saltelli (brute force) does not seem to work; not sure
setseed(seed)
numsamples = round(Int,sqrt(N))
numoneparamsamples = numsamples
nummanyparamsamples = numsamples
# convert the distribution strings into actual distributions
paramkeys = getoptparamkeys(madsdata)
# find the mean and variance
f = makemadscommandfunction(madsdata)
distributions = getparamdistributions(madsdata)
results = Array{OrderedCollections.OrderedDict}(undef, numsamples)
paramdict = Mads.getparamdict(madsdata)
for i = 1:numsamples
for j in 1:length(paramkeys)
paramdict[paramkeys[j]] = Distributions.rand(distributions[paramkeys[j]]) # TODO use parametersample
end
results[i] = f(paramdict) # this got to be slow to process
end
obskeys = getobskeys(madsdata)
sum = OrderedCollections.OrderedDict{String,Float64}()
for i = 1:length(obskeys)
sum[obskeys[i]] = 0.
end
for j = 1:numsamples
for i = 1:length(obskeys)
sum[obskeys[i]] += results[j][obskeys[i]]
end
end
mean = OrderedCollections.OrderedDict{String,Float64}()
for i = 1:length(obskeys)
mean[obskeys[i]] = sum[obskeys[i]] / numsamples
end
for i = 1:length(paramkeys)
sum[paramkeys[i]] = 0.
end
for j = 1:numsamples
for i = 1:length(obskeys)
sum[obskeys[i]] += (results[j][obskeys[i]] - mean[obskeys[i]]) ^ 2
end
end
variance = OrderedCollections.OrderedDict{String,Float64}()
for i = 1:length(obskeys)
variance[obskeys[i]] = sum[obskeys[i]] / (numsamples - 1)
end
madsinfo("Compute the main effect (first order) sensitivities (indices)")
mes = OrderedCollections.OrderedDict{String,OrderedCollections.OrderedDict}()
for k = 1:length(obskeys)
mes[obskeys[k]] = OrderedCollections.OrderedDict{String,Float64}()
end
for i = 1:length(paramkeys)
madsinfo("Parameter : $(paramkeys[i])")
cond_means = Array{OrderedCollections.OrderedDict}(undef, numoneparamsamples)
@ProgressMeter.showprogress 1 "Computing ... " for j = 1:numoneparamsamples
cond_means[j] = OrderedCollections.OrderedDict()
for k = 1:length(obskeys)
cond_means[j][obskeys[k]] = 0.
end
paramdict[paramkeys[i]] = Distributions.rand(distributions[paramkeys[i]]) # TODO use parametersample
for k = 1:nummanyparamsamples
for m = 1:length(paramkeys)
if m != i
paramdict[paramkeys[m]] = Distributions.rand(distributions[paramkeys[m]]) # TODO use parametersample
end
end
results = f(paramdict)
for k = 1:length(obskeys)
cond_means[j][obskeys[k]] += results[obskeys[k]]
end
end
for k = 1:length(obskeys)
cond_means[j][obskeys[k]] /= nummanyparamsamples
end
end
v = Array{Float64}(undef, numoneparamsamples)
for k = 1:length(obskeys)
for j = 1:numoneparamsamples
v[j] = cond_means[j][obskeys[k]]
end
mes[obskeys[k]][paramkeys[i]] = Statistics.std(v) ^ 2 / variance[obskeys[k]]
end
end
madsinfo("Compute the total effect sensitivities (indices)") # TODO we should use the same samples for total and main effect
tes = OrderedCollections.OrderedDict{String,OrderedCollections.OrderedDict}()
var = OrderedCollections.OrderedDict{String,OrderedCollections.OrderedDict}()
for k = 1:length(obskeys)
tes[obskeys[k]] = OrderedCollections.OrderedDict{String,Float64}()
var[obskeys[k]] = OrderedCollections.OrderedDict{String,Float64}()
end
for i = 1:length(paramkeys)
madsinfo("Parameter : $(paramkeys[i])")
cond_vars = Array{OrderedCollections.OrderedDict}(undef, nummanyparamsamples)
cond_means = Array{OrderedCollections.OrderedDict}(undef, nummanyparamsamples)
@ProgressMeter.showprogress 1 "Computing ... " for j = 1:nummanyparamsamples
cond_vars[j] = OrderedCollections.OrderedDict{String,Float64}()
cond_means[j] = OrderedCollections.OrderedDict{String,Float64}()
for m = 1:length(obskeys)
cond_means[j][obskeys[m]] = 0.
cond_vars[j][obskeys[m]] = 0.
end
for m = 1:length(paramkeys)
if m != i
paramdict[paramkeys[m]] = Distributions.rand(distributions[paramkeys[m]]) # TODO use parametersample
end
end
results = Array{OrderedCollections.OrderedDict}(undef, numoneparamsamples)
for k = 1:numoneparamsamples
paramdict[paramkeys[i]] = Distributions.rand(distributions[paramkeys[i]]) # TODO use parametersample
results[k] = f(paramdict)
for m = 1:length(obskeys)
cond_means[j][obskeys[m]] += results[k][obskeys[m]]
end
end
for m = 1:length(obskeys)
cond_means[j][obskeys[m]] /= numoneparamsamples
end
for k = 1:numoneparamsamples
for m = 1:length(obskeys)
cond_vars[j][obskeys[m]] += (results[k][obskeys[m]] - cond_means[j][obskeys[m]]) ^ 2
end
end
for m = 1:length(obskeys)
cond_vars[j][obskeys[m]] /= numoneparamsamples - 1
end
end
for j = 1:length(obskeys)
runningsum = 0.
for m = 1:nummanyparamsamples
runningsum += cond_vars[m][obskeys[j]]
end
tes[obskeys[j]][paramkeys[i]] = runningsum / nummanyparamsamples / variance[obskeys[j]]
var[obskeys[j]][paramkeys[i]] = runningsum / nummanyparamsamples
end
end
Dict("mes" => mes, "tes" => tes, "var" => var, "samplesize" => N, "seed" => seed, "method" => "saltellibrute")
end
"""
Load Saltelli sensitivity analysis results for fast simulation restarts
$(DocumentFunction.documentfunction(loadsaltellirestart!;
argtext=Dict("evalmat"=>"loaded array",
"matname"=>"matrix (array) name (defines the name of the loaded file)",
"restartdir"=>"directory where files will be stored containing model results for fast simulation restarts")))
Returns:
- `true` when successfully loaded, `false` when it is not
"""
function loadsaltellirestart!(evalmat::Array, matname::AbstractString, restartdir::AbstractString)
if restartdir == ""
return false
end
filename = joinpath(restartdir, string(matname, "_", Distributed.myid(), ".jld"))
if !isfile(filename)
return false
end
mat = FileIO.load(filename, "mat")
copy!(evalmat, mat)
return true
end
"""
Save Saltelli sensitivity analysis results for fast simulation restarts
$(DocumentFunction.documentfunction(savesaltellirestart;
argtext=Dict("evalmat"=>"saved array",
"matname"=>"matrix (array) name (defines the name of the loaded file)",
"restartdir"=>"directory where files will be stored containing model results for fast simulation restarts")))
"""
function savesaltellirestart(evalmat::Array, matname::AbstractString, restartdir::AbstractString)
if restartdir != ""
Mads.mkdir(restartdir)
FileIO.save(joinpath(restartdir, string(matname, "_", Distributed.myid(), ".jld2")), "mat", evalmat)
end
return nothing
end
"""
Saltelli sensitivity analysis
$(DocumentFunction.documentfunction(saltelli;
argtext=Dict("madsdata"=>"MADS problem dictionary"),
keytext=Dict("N"=>"number of samples [default=`100`]",
"seed"=>"random seed [default=`0`]",
"restartdir"=>"directory where files will be stored containing model results for fast simulation restarts",
"parallel"=>"set to true if the model runs should be performed in parallel [default=`false`]",
"checkpointfrequency"=>"check point frequency [default=`N`]")))
"""
function saltelli(madsdata::AbstractDict; N::Integer=100, seed::Integer=-1, restartdir::AbstractString="", parallel::Bool=false, checkpointfrequency::Integer=N)
Mads.setseed(seed)
Mads.madsoutput("Number of samples: $N\n");
paramallkeys = Mads.getparamkeys(madsdata)
paramalldict = OrderedCollections.OrderedDict{String,Float64}(zip(paramallkeys, Mads.getparamsinit(madsdata)))
paramoptkeys = Mads.getoptparamkeys(madsdata)
nP = length(paramoptkeys)
Mads.madsoutput("Number of model paramters to be analyzed: $(nP) \n");
Mads.madsoutput("Number of model evaluations to be perforemed: $(N * 2 + N * nP) \n");
obskeys = Mads.getobskeys(madsdata)
nO = length(obskeys)
distributions = Mads.getparamdistributions(madsdata)
f = Mads.makemadscommandfunction(madsdata)
A = Array{Float64}(undef, N, 0)
B = Array{Float64}(undef, N, 0)
C = Array{Float64}(undef, N, nP)
variance = OrderedCollections.OrderedDict{String, OrderedCollections.OrderedDict{String, Float64}}() # variance
mes = OrderedCollections.OrderedDict{String, OrderedCollections.OrderedDict{String, Float64}}() # main effect (first order) sensitivities
tes = OrderedCollections.OrderedDict{String, OrderedCollections.OrderedDict{String, Float64}}() # total effect sensitivities
for i = 1:nO
variance[obskeys[i]] = OrderedCollections.OrderedDict{String, Float64}()
mes[obskeys[i]] = OrderedCollections.OrderedDict{String, Float64}()
tes[obskeys[i]] = OrderedCollections.OrderedDict{String, Float64}()
end
for key in paramoptkeys
delete!(paramalldict, key)
end
for j = 1:nP
s1 = Mads.getparamrandom(madsdata, N, paramoptkeys[j])
s2 = Mads.getparamrandom(madsdata, N, paramoptkeys[j])
A = [A s1]
B = [B s2]
end
Mads.madsoutput( "Computing model outputs to calculate total output mean and variance ... Sample A ...\n" );
function farray(Ai)
feval = f(merge(paramalldict, OrderedCollections.OrderedDict{String,Float64}(zip(paramoptkeys, Ai))))
result = Array{Float64}(undef, length(obskeys))
for i = 1:length(obskeys)
result[i] = feval[obskeys[i]]
end
return result
end
yA = Array{Float64}(undef, N, length(obskeys))
if restartdir == ""
restartdir = getrestartdir(madsdata)
end
flagrestart = restartdir != ""
if parallel
Avecs = Array{Array{Float64, 1}}(undef, size(A, 1))
for i = 1:N
Avecs[i] = vec(A[i, :])
end
if flagrestart
pmapresult = RobustPmap.crpmap(farray, checkpointfrequency, joinpath(restartdir, "yA"), Avecs; t=Array{Float64, 1})
else
pmapresult = RobustPmap.rpmap(farray, Avecs; t=Array{Float64, 1})
end
for i = 1:N
for j = 1:length(obskeys)
yA[i, j] = pmapresult[i][j]
end
end
else
if !loadsaltellirestart!(yA, "yA", restartdir)
for i = 1:N
feval = f(merge(paramalldict, OrderedCollections.OrderedDict{String, Float64}(zip(paramoptkeys, A[i, :]))))
for j = 1:length(obskeys)
yA[i, j] = feval[obskeys[j]]
end
end
savesaltellirestart(yA, "yA", restartdir)
end
end
Mads.madsoutput( "Computing model outputs to calculate total output mean and variance ... Sample B ...\n" );
yB = Array{Float64}(undef, N, length(obskeys))
if parallel
Bvecs = Array{Array{Float64, 1}}(undef, size(B, 1))
for i = 1:N
Bvecs[i] = vec(B[i, :])
end
if flagrestart
pmapresult = RobustPmap.crpmap(farray, checkpointfrequency, joinpath(restartdir, "yB"), Bvecs; t=Array{Float64, 1})
else
pmapresult = RobustPmap.rpmap(farray, Bvecs; t=Array{Float64, 1})
end
for i = 1:N
for j = 1:length(obskeys)
yB[i, j] = pmapresult[i][j]
end
end
else
if !loadsaltellirestart!(yB, "yB", restartdir)
for i = 1:N
feval = f(merge(paramalldict, OrderedCollections.OrderedDict{String,Float64}(zip(paramoptkeys, B[i, :]))))
for j = 1:length(obskeys)
yB[i, j] = feval[obskeys[j]]
end
end
savesaltellirestart(yB, "yB", restartdir)
end
end
for i = 1:nP
for j = 1:N
for k = 1:nP
if k == i
C[j, k] = B[j, k]
else
C[j, k] = A[j, k]
end
end
end
Mads.madsoutput( "Computing model outputs to calculate total output mean and variance ... Sample C ... Parameter $(paramoptkeys[i])\n" );
yC = Array{Float64}(undef, N, length(obskeys))
if parallel
Cvecs = Array{Array{Float64, 1}}(undef, size(C, 1))
for j = 1:N
Cvecs[j] = vec(C[j, :])
end
if flagrestart
pmapresult = RobustPmap.crpmap(farray, checkpointfrequency, joinpath(restartdir, "yC$i"), Cvecs; t=Array{Float64, 1})
else
pmapresult = RobustPmap.rpmap(farray, Cvecs; t=Array{Float64, 1})
end
for j = 1:N
for k = 1:length(obskeys)
yC[j, k] = pmapresult[j][k]
end
end
else
if !loadsaltellirestart!(yC, "yC$(i)", restartdir)
for j = 1:N
feval = f(merge(paramalldict, OrderedCollections.OrderedDict{String,Float64}(zip(paramoptkeys, C[j, :]))))
for k = 1:length(obskeys)
yC[j, k] = feval[obskeys[k]]
end
end
savesaltellirestart(yC, "yC$(i)", restartdir)
end
end
maxnnans = 0
for j = 1:nO
yAnonan = isnan.(yA[:,j])
yBnonan = isnan.(yB[:,j])
yCnonan = isnan.(yC[:,j])
nonan = ( yAnonan .+ yBnonan .+ yCnonan ) .== 0
yT = vcat( yA[nonan,j], yB[nonan,j] ) # this should not include C
nanindices = findall(map(~, nonan))
# println("$nanindices")
nnans = length(nanindices)
if nnans > maxnnans
maxnnans = nnans
end
nnonnans = N - nnans
# f0T = Statistics.mean(yT)
# f0A = Statistics.mean( yA[nonan,j] )
# f0B = Statistics.mean( yB[nonan,j] )
# f0C = Statistics.mean( yC[nonan,j] )
varT = Statistics.var(yT)
# varA = abs( ( LinearAlgebra.dot( yA[nonan,j], yA[nonan,j] ) - f0A^2 * nnonnans ) / ( nnonnans - 1 ) )
# varA = Statistics.var( yA[nonan,j] ) # this is faster
# varB = abs( ( LinearAlgebra.dot( yB[nonan,j], yB[nonan,j] ) - f0B^2 * nnonnans ) / ( nnonnans - 1 ) )
# varB = Statistics.var( yB[nonan,j] ) # this is faster
varC = Statistics.var(yC[nonan,j])
# varP = abs( ( LinearAlgebra.dot( yB[nonan, j], yC[nonan, j] ) / nnonnans - f0B * f0C ) ) # Orignial
# varP2 = abs( ( LinearAlgebra.dot( yB[nonan, j], yC[nonan, j] ) - f0B * f0C * nnonnans ) / ( nnonnans - 1 ) ) # Imporved
varP3 = abs( Statistics.mean( yB[nonan,j] .* ( yC[nonan, j] - yA[nonan,j] ) ) ) # Recommended
# varP4 = Statistics.mean( ( yB[nonan,j] - yC[nonan, j] ).^2 ) / 2 # Mads.c; very different from all the other estimates
# println("varP $varP varP2 $varP2 varP3 $varP3 varP4 $varP4")
# varPnot = abs( ( LinearAlgebra.dot( yA[nonan, j], yC[nonan, j] ) / nnonnans - f0A * f0C ) ) # Orignial
# varPnot2 = abs( ( LinearAlgebra.dot( yA[nonan, j], yC[nonan, j] ) - f0A * f0C * nnonnans ) / ( nnonnans - 1 ) ) # Imporved
# varPnot3 = Statistics.mean( ( yA[nonan,j] - yC[nonan, j] ).^2 ) / 2 # Recommended; also used in Mads.c
expPnot = Statistics.mean( ( yA[nonan,j] - yC[nonan, j] ).^2 ) / 2
# println("varPnot $varPnot varPnot2 $varPnot2 varPnot3 $varPnot3")
variance[obskeys[j]][paramoptkeys[i]] = varC
# if varA < eps(Float64) && varP < eps(Float64)
# tes[obskeys[j]][paramoptkeys[i]] = mes[obskeys[j]][paramoptkeys[i]] = NaN
# else
# mes[obskeys[j]][paramoptkeys[i]] = min( 1, max( 0, varP / varA ) ) # varT or varA? i think it should be varA
# tes[obskeys[j]][paramoptkeys[i]] = min( 1, max( 0, 1 - varPnot / varB) ) # varT or varA; i think it should be varA; i do not think should be varB?
# end
# mes[obskeys[j]][paramoptkeys[i]] = varP / varT
# tes[obskeys[j]][paramoptkeys[i]] = 1 - varPnot / varT
# varianceA[obskeys[j]][paramoptkeys[i]] = varT
# varianceB[obskeys[j]][paramoptkeys[i]] = 1 - varPnot / varB
# varianceB[obskeys[j]][paramoptkeys[i]] = varC
mes[obskeys[j]][paramoptkeys[i]] = varP3 / varT
# tes[obskeys[j]][paramoptkeys[i]] = 1 - varPnot / varB
tes[obskeys[j]][paramoptkeys[i]] = expPnot / varT
# println("N $N nnonnans $nnonnans f0A $f0A f0B $f0B varA $varA varB $varB varP $varP varPnot $varPnot mes $(varP / varA) tes $(1 - varPnot / varB)")
end
if maxnnans > 0
Mads.madswarn("There are $(maxnnans) NaN's")
end
end
Dict("mes" => mes, "tes" => tes, "var" => variance, "samplesize" => N, "seed" => seed, "method" => "saltelli")
end
"""
Compute sensitivities for each model parameter; averaging the sensitivity indices over the entire observation range
$(DocumentFunction.documentfunction(computeparametersensitities;
argtext=Dict("madsdata"=>"MADS problem dictionary",
"saresults"=>"dictionary with sensitivity analysis results")))
"""
function computeparametersensitities(madsdata::AbstractDict, saresults::AbstractDict)
paramkeys = getoptparamkeys(madsdata)
obskeys = getobskeys(madsdata)
mes = saresults["mes"]
tes = saresults["tes"]
var = saresults["var"]
pvar = OrderedCollections.OrderedDict{String, Float64}() # parameter variance
pmes = OrderedCollections.OrderedDict{String, Float64}() # parameter main effect (first order) sensitivities
ptes = OrderedCollections.OrderedDict{String, Float64}() # parameter total effect sensitivities
for i = 1:length(paramkeys)
pv = pm = pt = 0
for j = 1:length(obskeys)
m = typeof(saresults["mes"][obskeys[j]][paramkeys[i]]) == Nothing ? 0 : saresults["mes"][obskeys[j]][paramkeys[i]]
t = typeof(saresults["tes"][obskeys[j]][paramkeys[i]]) == Nothing ? 0 : saresults["tes"][obskeys[j]][paramkeys[i]]
v = typeof(saresults["var"][obskeys[j]][paramkeys[i]]) == Nothing ? 0 : saresults["var"][obskeys[j]][paramkeys[i]]
pv += isnan.(v) ? 0 : v
pm += isnan.(m) ? 0 : m
pt += isnan.(t) ? 0 : t
end
pvar[paramkeys[i]] = pv / length(obskeys)
pmes[paramkeys[i]] = pm / length(obskeys)
ptes[paramkeys[i]] = pt / length(obskeys)
end
Dict("var" => pvar, "mes" => pmes, "tes" => ptes)
end
# Parallelization of Saltelli functions
saltelli_functions = ["saltelli", "saltellibrute"]
global index = 0
for mi = 1:length(saltelli_functions)
global index = mi
q = quote
"""
Parallel version of $(saltelli_functions[index])
"""
function $(Symbol(string(saltelli_functions[index], "parallel")))(madsdata::AbstractDict, numsaltellis::Integer; N::Integer=100, seed::Integer=-1, restartdir::AbstractString="")
Mads.setseed(seed)
if numsaltellis < 1
madserror("Number of parallel sensitivity runs must be > 0 ($numsaltellis < 1)")
return
end
results = RobustPmap.rpmap(i->$(Symbol(saltelli_functions[index]))(madsdata; N=N, seed=seed+i, restartdir=restartdir), 1:numsaltellis)
mesall = results[1]["mes"]
tesall = results[1]["tes"]
varall = results[1]["var"]
for i = 2:numsaltellis
mes, tes, var = results[i]["mes"], results[i]["tes"], results[i]["var"]
for obskey in keys(mes)
for paramkey in keys(mes[obskey])
#meanall[obskey][paramkey] += mean[obskey][paramkey]
#varianceall[obskey][paramkey] += variance[obskey][paramkey]
mesall[obskey][paramkey] += mes[obskey][paramkey]
tesall[obskey][paramkey] += tes[obskey][paramkey]
varall[obskey][paramkey] += var[obskey][paramkey]
end
end
end
for obskey in keys(mesall)
for paramkey in keys(mesall[obskey])
#meanall[obskey][paramkey] /= numsaltellis
#varianceall[obskey][paramkey] /= numsaltellis
mesall[obskey][paramkey] /= numsaltellis
tesall[obskey][paramkey] /= numsaltellis
varall[obskey][paramkey] /= numsaltellis
end
end
Dict("mes" => mesall, "tes" => tesall, "var" => varall, "samplesize" => N * numsaltellis, "seed" => seed, "method" => $(saltelli_functions[index])*"_parallel")
end # end function
end # end quote
Core.eval(Mads, q)
end
"""
Print sensitivity analysis results
$(DocumentFunction.documentfunction(printSAresults;
argtext=Dict("madsdata"=>"MADS problem dictionary",
"results"=>"dictionary with sensitivity analysis results")))
"""
function printSAresults(madsdata::AbstractDict, results::AbstractDict)
mes = results["mes"]
tes = results["tes"]
N = results["samplesize"]
#=
madsoutput("Mean\n")
madsoutput("\t")
obskeys = getobskeys(madsdata)
paramkeys = getparamkeys(madsdata)
for paramkey in paramkeys
madsoutput("\t$(paramkey)")
end
madsoutput("\n")
for obskey in obskeys
madsoutput(obskey)
for paramkey in paramkeys
madsoutput("\t$(mean[obskey][paramkey])")
end
madsoutput("\n")
end
madsoutput("\nVariance\n")
madsoutput("\t")
obskeys = getobskeys(madsdata)
paramkeys = getparamkeys(madsdata)
for paramkey in paramkeys
madsoutput("\t$(paramkey)")
end
madsoutput("\n")
for obskey in obskeys
madsoutput(obskey)
for paramkey in paramkeys
madsoutput("\t$(variance[obskey][paramkey])")
end
madsoutput("\n")
end
=#
print("\nMain Effect Indices\n")
print("obs")
obskeys = getobskeys(madsdata)
paramkeys = getoptparamkeys(madsdata)
for paramkey in paramkeys
print("\t$(paramkey)")
end
print("\tSum")
print("\n")
for obskey in obskeys
print(obskey)
sum = 0
for paramkey in paramkeys
sum += mes[obskey][paramkey]
print("\t$(@Printf.sprintf("%f",mes[obskey][paramkey]))")
end
print("\t$(sum)")
print("\n")
end
print("\nTotal Effect Indices\n")
print("obs")
for paramkey in paramkeys
print("\t$(paramkey)")
end
print("\tSum")
print("\n")
for obskey in obskeys
print(obskey)
sum = 0
for paramkey in paramkeys
sum += tes[obskey][paramkey]
print("\t$(@Printf.sprintf("%f",tes[obskey][paramkey]))")
end
print("\t$(sum)")
print("\n")
end
print("\n")
end
"""
Print sensitivity analysis results (method 2)
$(DocumentFunction.documentfunction(printSAresults2;
argtext=Dict("madsdata"=>"MADS problem dictionary",
"results"=>"dictionary with sensitivity analysis results")))
"""
function printSAresults2(madsdata::AbstractDict, results::AbstractDict)