-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathyelmox.f90
1244 lines (916 loc) · 50.2 KB
/
yelmox.f90
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
program yelmox
use nml
use ncio
use timer
use timeout
use yelmo
use yelmo_tools, only : smooth_gauss_2D
use ice_optimization
use basal_dragging, only : calc_lambda_bed_lin, calc_lambda_bed_exp
! External libraries
use fastisostasy ! also reexports barysealevel
use snapclim
use marine_shelf
use smbpal
use sediments
use geothermal
implicit none
type(yelmo_class) :: yelmo1
type(bsl_class) :: bsl
type(snapclim_class) :: snp1
type(marshelf_class) :: mshlf1
type(smbpal_class) :: smbpal1
type(sediments_class) :: sed1
type(geothermal_class) :: gthrm1
type(isos_class) :: isos1
character(len=256) :: outfldr, file1D, file2D, file2D_small, domain
character(len=256) :: file_isos, file_bsl
character(len=512) :: path_par
character(len=512) :: path_lgm
real(wp) :: time, time_bp, time_elapsed
real(wp) :: dT_now
integer :: n
type(timeout_class) :: tm_1D, tm_2D, tm_2Dsm
! Model timing
type(timer_class) :: tmr
type(timer_class) :: tmrs
character(len=512) :: tmr_file
type reg_def_class
character(len=56) :: name
character(len=512) :: fnm
logical, allocatable :: mask(:,:)
logical :: write
end type
type(reg_def_class) :: reg1
type(reg_def_class) :: reg2
type(reg_def_class) :: reg3
character(len=512) :: regions_mask_fnm
real(wp), allocatable :: regions_mask(:,:)
type ctrl_params
character(len=56) :: run_step
real(wp) :: time_init
real(wp) :: time_end
real(wp) :: time_equil ! Only for spinup
real(wp) :: time_const ! Only for spinup
real(wp) :: time_lgm_step
real(wp) :: time_pd_step
real(wp) :: dtt
real(wp) :: dt_restart
real(wp) :: dt_clim
logical :: transient_clim
logical :: use_lgm_step
logical :: use_pd_step
logical :: with_ice_sheet
character(len=56) :: equil_method
real(wp) :: isos_tau_1
real(wp) :: isos_tau_2
real(wp) :: isos_sigma
end type
type negis_params
logical :: use_negis_par
real(wp) :: cf_0
real(wp) :: cf_1
real(wp) :: cf_centre
real(wp) :: cf_north
real(wp) :: cf_south
real(wp) :: cf_x
end type
type stats_class
logical :: defined
real(wp) :: pd_rmse_H
real(wp) :: pd_rmse_H_flt
real(wp) :: V
real(wp) :: A
real(wp) :: A_grnd
real(wp) :: A_flt
end type
type(ctrl_params) :: ctl
type(ice_opt_params) :: opt
type(stats_class) :: stats_pd_1, stats_lgm, stats_pd_2
type(negis_params) :: ngs
! Internal parameters
logical :: running_laurentide
logical :: laurentide_init_const_H
real(wp) :: laurentide_time_equil
logical :: running_greenland
logical :: greenland_init_marine_H
logical :: scale_glacial_smb
real(wp), parameter :: time_lgm = -19050.0_wp ! [yr CE] == 21 kyr ago
real(wp), parameter :: time_pd = 1950.0_wp ! [yr CE] == 0 kyr ago
! Determine the parameter file from the command line
call yelmo_load_command_line_args(path_par)
! Timing and other parameters
call nml_read(path_par,"ctrl","time_init", ctl%time_init) ! [yr] Starting time
call nml_read(path_par,"ctrl","time_end", ctl%time_end) ! [yr] Ending time
call nml_read(path_par,"ctrl","time_equil", ctl%time_equil) ! [yr] Years to equilibrate first
call nml_read(path_par,"ctrl","time_const", ctl%time_const)
call nml_read(path_par,"ctrl","time_lgm_step", ctl%time_lgm_step)
call nml_read(path_par,"ctrl","time_pd_step", ctl%time_pd_step)
call nml_read(path_par,"ctrl","dtt", ctl%dtt) ! [yr] Main loop time step
call nml_read(path_par,"ctrl","dt_restart", ctl%dt_restart)
call nml_read(path_par,"ctrl","transient_clim", ctl%transient_clim) ! Calculate transient climate?
call nml_read(path_par,"ctrl","use_lgm_step", ctl%use_lgm_step) ! Use lgm_step?
call nml_read(path_par,"ctrl","use_pd_step", ctl%use_pd_step) ! Use pd_step?
call nml_read(path_par,"ctrl","with_ice_sheet", ctl%with_ice_sheet) ! Include an active ice sheet
call nml_read(path_par,"ctrl","equil_method", ctl%equil_method) ! What method should be used for spin-up?
! Get output times
call timeout_init(tm_1D, path_par,"tm_1D", "small", ctl%time_init,ctl%time_end)
call timeout_init(tm_2D, path_par,"tm_2D", "heavy", ctl%time_init,ctl%time_end)
call timeout_init(tm_2Dsm,path_par,"tm_2Dsm","medium", ctl%time_init,ctl%time_end)
! Hard-coded for now:
ctl%dt_clim = 10.0 ! [yrs] Frequency to update snapclim snapshot
! Consistency checks ===
! transient_clim overrides use_lgm_step
if (ctl%transient_clim) ctl%use_lgm_step = .FALSE.
! lgm step should only come after time_equil is finished...
if (ctl%time_lgm_step .lt. ctl%time_equil) then
write(5,*) ""
write(5,*) "yelmox:: time_lgm_step must be greater than time_equil."
write(5,*) "Try again."
stop "Program stopped."
end if
! pd step should only come after lgm step
if (ctl%time_pd_step .lt. ctl%time_lgm_step) then
write(5,*) ""
write(5,*) "yelmox:: time_pd_step must be greater than time_lgm_step."
write(5,*) "Try again."
stop "Program stopped."
end if
if (trim(ctl%equil_method) .eq. "opt") then
! Load optimization parameters
! Initially set to zero
opt%tf_basins = 0
call optimize_par_load(opt,path_par,"opt")
end if
! Set initial time
time = ctl%time_init
time_bp = time - 1950.0_wp
! Start timing
call timer_step(tmr,comp=-1)
! Assume program is running from the output folder
outfldr = "./"
! Define input and output locations
file1D = trim(outfldr)//"yelmo1D.nc"
file2D = trim(outfldr)//"yelmo2D.nc"
file2D_small = trim(outfldr)//"yelmo2Dsm.nc"
file_isos = trim(outfldr)//"fastisostasy.nc"
file_bsl = trim(outfldr)//"bsl.nc"
tmr_file = trim(outfldr)//"timer_table.txt"
! Print summary of run settings
write(*,*)
write(*,*) "transient_clim: ", ctl%transient_clim
write(*,*) "use_lgm_step: ", ctl%use_lgm_step
write(*,*) "use_pd_step: ", ctl%use_pd_step
write(*,*) "with_ice_sheet: ", ctl%with_ice_sheet
write(*,*) "equil_method: ", trim(ctl%equil_method)
write(*,*)
write(*,*) "time_init: ", ctl%time_init
write(*,*) "time_end: ", ctl%time_end
write(*,*) "dtt: ", ctl%dtt
write(*,*) "dt_restart: ", ctl%dt_restart
write(*,*)
if (.not. ctl%transient_clim) then
write(*,*) "time_equil: ", ctl%time_equil
write(*,*) "time_const: ", ctl%time_const
if (ctl%use_lgm_step) then
write(*,*) "time_lgm_step: ", ctl%time_lgm_step
end if
if (ctl%use_pd_step) then
write(*,*) "time_pd_step: ", ctl%time_pd_step
end if
! Set time before present == to constant climate time
time_bp = ctl%time_const - 1950.0_wp
end if
write(*,*) "time = ", time
write(*,*) "time_bp = ", time_bp
write(*,*)
! === Initialize ice sheet model =====
! Initialize data objects and load initial topography
call yelmo_init(yelmo1,filename=path_par,grid_def="file",time=time)
! Store domain name as a shortcut
domain = yelmo1%par%domain
! Ensure optimization fields are allocated and preassigned
allocate(opt%cf_min(yelmo1%grd%nx,yelmo1%grd%ny))
allocate(opt%cf_max(yelmo1%grd%nx,yelmo1%grd%ny))
opt%cf_min = opt%cf_min_par
opt%cf_max = yelmo1%dyn%par%till_cf_ref
! Define specific regions of interest =====================
! Shortcut switches for later use
running_laurentide = .FALSE.
running_greenland = .FALSE.
select case(trim(domain))
case("Antarctica")
! Define base regions for whole domain first
regions_mask_fnm = "ice_data/Antarctica/"//trim(yelmo1%par%grid_name)//&
"/"//trim(yelmo1%par%grid_name)//"_BASINS-nasa.nc"
allocate(regions_mask(yelmo1%grd%nx,yelmo1%grd%ny))
! Load mask from file
call nc_read(regions_mask_fnm,"mask_regions",regions_mask)
! APIS region (region=3.0 in regions map)
reg1%write = .TRUE.
reg1%name = "APIS"
reg1%fnm = trim(outfldr)//"yelmo1D_"//trim(reg1%name)//".nc"
allocate(reg1%mask(yelmo1%grd%nx,yelmo1%grd%ny))
reg1%mask = .FALSE.
where(abs(regions_mask - 3.0) .lt. 1e-3) reg1%mask = .TRUE.
! WAIS region (region=4.0 in regions map)
reg2%write = .TRUE.
reg2%name = "WAIS"
reg2%fnm = trim(outfldr)//"yelmo1D_"//trim(reg2%name)//".nc"
allocate(reg2%mask(yelmo1%grd%nx,yelmo1%grd%ny))
reg2%mask = .FALSE.
where(abs(regions_mask - 4.0) .lt. 1e-3) reg2%mask = .TRUE.
! EAIS region (region=5.0 in regions map)
reg3%write = .TRUE.
reg3%name = "EAIS"
reg3%fnm = trim(outfldr)//"yelmo1D_"//trim(reg3%name)//".nc"
allocate(reg3%mask(yelmo1%grd%nx,yelmo1%grd%ny))
reg3%mask = .FALSE.
where(abs(regions_mask - 5.0) .lt. 1e-3) reg3%mask = .TRUE.
case("Laurentide")
running_laurentide = .TRUE.
ctl%use_lgm_step = .FALSE.
ctl%use_pd_step = .FALSE.
laurentide_init_const_H = .FALSE.
if (laurentide_init_const_H) then
! Make sure relaxation spinup is short, but transient spinup
! with modified positive smb over North America is reasonably long.
ctl%time_equil = 10.0
laurentide_time_equil = 5e3
else
! When starting from ice-6g, positive smb spinup is not necessary.
! however ctl%time_equil should not be changed, as it may be useful
! for spinning up thermodynamics.
laurentide_time_equil = 0.0
end if
! Make sure to set ice_allowed to prevent ice from growing in
! Greenland (and on grid borders)
where(abs(yelmo1%bnd%regions - 1.30) .lt. 1e-3) yelmo1%bnd%ice_allowed = .FALSE.
yelmo1%bnd%ice_allowed(1,:) = .FALSE.
yelmo1%bnd%ice_allowed(yelmo1%grd%nx,:) = .FALSE.
yelmo1%bnd%ice_allowed(:,1) = .FALSE.
yelmo1%bnd%ice_allowed(:,yelmo1%grd%ny) = .FALSE.
! Hudson region (region=1.12 in regions map)
reg1%write = .TRUE.
reg1%name = "Hudson"
reg1%fnm = trim(outfldr)//"yelmo1D_"//trim(reg1%name)//".nc"
allocate(reg1%mask(yelmo1%grd%nx,yelmo1%grd%ny))
reg1%mask = .FALSE.
where(abs(yelmo1%bnd%regions - 1.12) .lt. 1e-3) reg1%mask = .TRUE.
reg2%write = .FALSE.
reg3%write = .FALSE.
case("Greenland")
running_greenland = .TRUE.
! Should extra ice be imposed over continental shelf to mimic LGM state to start
greenland_init_marine_H = .TRUE.
! Should glacial smb be modified to reduce negative smb values
scale_glacial_smb = .FALSE.
! Should NEGIS parameter modifications be used
ngs%use_negis_par = .TRUE.
if (.FALSE.) then
! ajr, 2025-01-15, missing these parameters in param files - ask Ilaria!!
! Load NEGIS parameters from file, if used
if (ngs%use_negis_par) then
call nml_read(path_par,"negis","cf_0", ngs%cf_0)
call nml_read(path_par,"negis","cf_1", ngs%cf_1)
call nml_read(path_par,"negis","cf_centre", ngs%cf_centre)
call nml_read(path_par,"negis","cf_north", ngs%cf_north)
call nml_read(path_par,"negis","cf_south", ngs%cf_south)
end if
end if
! Make sure to set ice_allowed to prevent ice from growing in
! Iceland and Svaalbard (on grid borders)
where(abs(yelmo1%bnd%regions - 1.20) .lt. 1e-3) yelmo1%bnd%ice_allowed = .FALSE.
where(abs(yelmo1%bnd%regions - 1.23) .lt. 1e-3) yelmo1%bnd%ice_allowed = .FALSE.
where(abs(yelmo1%bnd%regions - 1.31) .lt. 1e-3) yelmo1%bnd%ice_allowed = .FALSE.
if (yelmo1%dyn%par%till_method .eq. -1) then
! Initialize cb_ref to constant value to start
! (will be updated in the time loop)
yelmo1%dyn%now%cb_ref = yelmo1%dyn%par%till_cf_ref
end if
case DEFAULT
reg1%write = .FALSE.
reg2%write = .FALSE.
reg3%write = .FALSE.
end select
! === Initialize external models (forcing for ice sheet) ======
! Initialize barysealevel model
call bsl_init(bsl, path_par, time_bp)
! Initialize fastisosaty
call isos_init(isos1, path_par, "isos", yelmo1%grd%nx, yelmo1%grd%ny, &
yelmo1%grd%dx, yelmo1%grd%dy)
! Initialize "climate" model (climate and ocean forcing)
call snapclim_init(snp1,path_par,domain,yelmo1%par%grid_name,yelmo1%grd%nx,yelmo1%grd%ny,yelmo1%bnd%basins)
! Initialize surface mass balance model (bnd%smb, bnd%T_srf)
call smbpal_init(smbpal1,path_par,x=yelmo1%grd%xc,y=yelmo1%grd%yc,lats=yelmo1%grd%lat)
! Initialize marine melt model (bnd%bmb_shlf)
call marshelf_init(mshlf1,path_par,"marine_shelf",yelmo1%grd%nx,yelmo1%grd%ny,domain,yelmo1%par%grid_name,yelmo1%bnd%regions,yelmo1%bnd%basins)
! Load other constant boundary variables (bnd%H_sed, bnd%Q_geo)
call sediments_init(sed1,path_par,yelmo1%grd%nx,yelmo1%grd%ny,domain,yelmo1%par%grid_name)
yelmo1%bnd%H_sed = sed1%now%H
call geothermal_init(gthrm1,path_par,yelmo1%grd%nx,yelmo1%grd%ny,domain,yelmo1%par%grid_name)
yelmo1%bnd%Q_geo = gthrm1%now%ghf
! === Update initial boundary conditions for current time and yelmo state =====
! ybound: z_bed, z_sl, H_sed, smb, T_srf, bmb_shlf , Q_geo
! Barystatic sea level
call bsl_update(bsl, year_bp=time_bp)
call bsl_write_init(bsl, file_bsl, time)
! Initialize the isostasy reference state using reference topography fields
call isos_init_ref(isos1, yelmo1%bnd%z_bed_ref, yelmo1%bnd%H_ice_ref)
call isos_init_state(isos1, yelmo1%bnd%z_bed, yelmo1%tpo%now%H_ice, time, bsl)
call isos_write_init_extended(isos1, file_isos, time)
yelmo1%bnd%z_bed = isos1%out%z_bed
yelmo1%bnd%z_sl = isos1%out%z_ss
! Update snapclim
call snapclim_update(snp1,z_srf=yelmo1%tpo%now%z_srf,time=time_bp,domain=domain,dx=yelmo1%grd%dx,basins=yelmo1%bnd%basins)
! Equilibrate snowpack for itm
if (trim(smbpal1%par%abl_method) .eq. "itm") then
call smbpal_update_monthly_equil(smbpal1,snp1%now%tas,snp1%now%pr, &
yelmo1%tpo%now%z_srf,yelmo1%tpo%now%H_ice,time_bp,time_equil=100.0)
end if
! Testing related to present-day surface mass balance
! snp1%now%tas = snp1%now%tas + 0.0
! snp1%now%pr = snp1%now%pr*exp(0.05*(0.0))
! call smbpal_update_monthly(smbpal1,snp1%now%tas,snp1%now%pr, &
! yelmo1%tpo%now%z_srf,yelmo1%tpo%now%H_ice,time, &
! file_out=trim(outfldr)//"smbpal.nc",write_now=.TRUE.,write_init=.TRUE.)
! stop
call smbpal_update_monthly(smbpal1,snp1%now%tas,snp1%now%pr, &
yelmo1%tpo%now%z_srf,yelmo1%tpo%now%H_ice,time_bp)
yelmo1%bnd%smb = smbpal1%ann%smb*yelmo1%bnd%c%conv_we_ie*1e-3 ! [mm we/a] => [m ie/a]
yelmo1%bnd%T_srf = smbpal1%ann%tsrf
if (trim(yelmo1%par%domain) .eq. "Greenland" .and. scale_glacial_smb) then
! Modify glacial smb
call calc_glacial_smb(yelmo1%bnd%smb,yelmo1%grd%lat,snp1%now%ta_ann,snp1%clim0%ta_ann)
end if
call marshelf_update_shelf(mshlf1,yelmo1%tpo%now%H_ice,yelmo1%bnd%z_bed,yelmo1%tpo%now%f_grnd, &
yelmo1%bnd%basins,yelmo1%bnd%z_sl,yelmo1%grd%dx,snp1%now%depth, &
snp1%now%to_ann,snp1%now%so_ann,dto_ann=snp1%now%to_ann-snp1%clim0%to_ann)
call marshelf_update(mshlf1,yelmo1%tpo%now%H_ice,yelmo1%bnd%z_bed,yelmo1%tpo%now%f_grnd, &
yelmo1%bnd%regions,yelmo1%bnd%basins,yelmo1%bnd%z_sl,dx=yelmo1%grd%dx)
yelmo1%bnd%bmb_shlf = mshlf1%now%bmb_shlf
yelmo1%bnd%T_shlf = mshlf1%now%T_shlf
call yelmo_print_bound(yelmo1%bnd)
! Initialize state variables (dyn,therm,mat)
! (initialize temps with robin method with a cold base)
call yelmo_init_state(yelmo1,time=time,thrm_method="robin-cold")
! ===== basal friction optimization ======
if (trim(ctl%equil_method) .eq. "opt") then
! Ensure that cb_ref will be optimized (till_method == set externally)
yelmo1%dyn%par%till_method = -1
! If not using restart, prescribe cb_ref to initial guess
if (.not. yelmo1%par%use_restart) then
yelmo1%dyn%now%cb_ref = opt%cf_init
end if
end if
! ========================================
if (.not. yelmo1%par%use_restart) then
! No restart file used, perform various initialization steps
if (running_laurentide) then
! Start with some ice thickness for testing
! Load LGM reconstruction into reference ice thickness
path_lgm = "ice_data/Laurentide/"//trim(yelmo1%par%grid_name)//&
"/"//trim(yelmo1%par%grid_name)//"_TOPO-ICE-6G_C.nc"
call nc_read(path_lgm,"dz",yelmo1%bnd%H_ice_ref,start=[1,1,1], &
count=[yelmo1%tpo%par%nx,yelmo1%tpo%par%ny,1])
! Start with some ice cover to speed up initialization
if (laurentide_init_const_H) then
yelmo1%tpo%now%H_ice = 0.0
where (yelmo1%bnd%regions .eq. 1.1 .and. yelmo1%bnd%z_bed .gt. 0.0) yelmo1%tpo%now%H_ice = 1000.0
where (yelmo1%bnd%regions .eq. 1.12) yelmo1%tpo%now%H_ice = 1000.0
else
! Set LGM reconstruction as initial ice thickness over North America
where ( yelmo1%bnd%z_bed .gt. -500.0 .and. &
( yelmo1%bnd%regions .eq. 1.1 .or. &
yelmo1%bnd%regions .eq. 1.11 .or. &
yelmo1%bnd%regions .eq. 1.12) )
yelmo1%tpo%now%H_ice = yelmo1%bnd%H_ice_ref
end where
! Apply Gaussian smoothing to keep things stable
call smooth_gauss_2D(yelmo1%tpo%now%H_ice,dx=yelmo1%grd%dx,f_sigma=2.0)
end if
! Load sediment mask
path_lgm = "ice_data/Laurentide/"//trim(yelmo1%par%grid_name)//&
"/"//trim(yelmo1%par%grid_name)//"_SED-L97.nc"
call nc_read(path_lgm,"z_sed",yelmo1%bnd%H_sed)
if (ctl%with_ice_sheet) then
! Run Yelmo for briefly to update surface topography
call yelmo_update_equil(yelmo1,time,time_tot=1.0_prec,dt=1.0,topo_fixed=.TRUE.)
! Addtional cleanup - remove floating ice
where( yelmo1%tpo%now%mask_bed .eq. 5) yelmo1%tpo%now%H_ice = 0.0
call yelmo_update_equil(yelmo1,time,time_tot=1.0_prec,dt=1.0,topo_fixed=.TRUE.)
end if
! Update snapclim to reflect new topography
call snapclim_update(snp1,z_srf=yelmo1%tpo%now%z_srf,time=time,domain=domain,dx=yelmo1%grd%dx,basins=yelmo1%bnd%basins)
! Update smbpal
call smbpal_update_monthly(smbpal1,snp1%now%tas,snp1%now%pr, &
yelmo1%tpo%now%z_srf,yelmo1%tpo%now%H_ice,time_bp)
yelmo1%bnd%smb = smbpal1%ann%smb*yelmo1%bnd%c%conv_we_ie*1e-3 ! [mm we/a] => [m ie/a]
yelmo1%bnd%T_srf = smbpal1%ann%tsrf
if (laurentide_init_const_H) then
! Additionally ensure smb is postive for land above 50degN in Laurentide region
! to make sure ice grows everywhere needed (Coridilleran ice sheet mainly)
where (yelmo1%bnd%regions .eq. 1.1 .and. yelmo1%grd%lat .gt. 50.0 .and. &
yelmo1%bnd%z_bed .gt. 0.0 .and. yelmo1%bnd%smb .lt. 0.0 ) yelmo1%bnd%smb = 0.5
if (ctl%with_ice_sheet) then
! Run yelmo for several years to ensure stable central ice dome
call yelmo_update_equil(yelmo1,time,time_tot=5e3,dt=5.0,topo_fixed=.FALSE.)
end if
else
if (ctl%with_ice_sheet) then
! Run yelmo for several years with constant boundary conditions to stabilize fields
call yelmo_update_equil(yelmo1,time,time_tot=1e2,dt=5.0,topo_fixed=.FALSE.)
end if
end if
else if (running_greenland) then
! Special start-up steps for Greenland
if (ngs%use_negis_par) then
! Ensure till method is correct, since we are updating cb_ref externally
yelmo1%dyn%par%till_method = -1
! Update cb_ref using negis parameters
call negis_update_cb_ref(yelmo1,ngs,time)
end if
if (greenland_init_marine_H) then
! Add extra ice-thickness over continental shelf to start with
! an LGM-like state
where(yelmo1%bnd%ice_allowed .and. yelmo1%tpo%now%H_ice .lt. 600.0 &
.and. yelmo1%bnd%z_bed .gt. -500.0)
yelmo1%tpo%now%H_ice = 800.0
end where
if (ctl%with_ice_sheet) then
! Run yelmo for a few years with constant boundary conditions
! to synchronize all model fields a bit
call yelmo_update_equil(yelmo1,time,time_tot=10.0_prec, dt=1.0_prec,topo_fixed=.FALSE.)
end if
end if
else
! Run simple startup equilibration step
if (ctl%with_ice_sheet) then
! Run yelmo for a few years with constant boundary conditions
! to synchronize all model fields a bit
call yelmo_update_equil(yelmo1,time,time_tot=10.0_prec, dt=1.0_prec,topo_fixed=.FALSE.)
end if
end if
end if
! ===== Initialize output files =====
call yelmo_write_init(yelmo1,file2D,time_init=time,units="years")
call yelmo_write_reg_init(yelmo1,file1D,time_init=time,units="years",mask=yelmo1%bnd%ice_allowed)
call yelmo_write_init(yelmo1,file2D_small,time_init=time,units="years")
if (reg1%write) then
call yelmo_write_reg_init(yelmo1,reg1%fnm,time_init=time,units="years",mask=reg1%mask)
end if
if (reg2%write) then
call yelmo_write_reg_init(yelmo1,reg2%fnm,time_init=time,units="years",mask=reg2%mask)
end if
if (reg3%write) then
call yelmo_write_reg_init(yelmo1,reg3%fnm,time_init=time,units="years",mask=reg3%mask)
end if
! Set stats
stats_pd_1%defined = .FALSE.
stats_lgm%defined = .FALSE.
stats_pd_2%defined = .FALSE.
call timer_step(tmr,comp=1,label="initialization")
call timer_step(tmrs,comp=-1)
! ==== Begin main time loop =====
! Advance timesteps
do n = 0, ceiling((ctl%time_end-ctl%time_init)/ctl%dtt)
! Get current time
time = ctl%time_init + n*ctl%dtt
if (ctl%transient_clim) then
time_bp = time - 1950.0_wp
else
if (ctl%use_lgm_step .and. time .ge. ctl%time_lgm_step) then
ctl%time_const = time_lgm
end if
if (ctl%use_pd_step .and. time .ge. ctl%time_pd_step) then
ctl%time_const = time_pd
end if
time_bp = ctl%time_const - 1950.0_wp
end if
time_elapsed = time - ctl%time_init
! Spin-up procedure - only relevant for time-time_init <= time_equil
select case(trim(ctl%equil_method))
case("opt")
if (time_elapsed .le. opt%rel_time2) then
! Apply relaxation to the model
! Update model relaxation time scale and error scaling (in [m])
call optimize_set_transient_param(opt%rel_tau,time_elapsed,time1=opt%rel_time1,time2=opt%rel_time2, &
p1=opt%rel_tau1,p2=opt%rel_tau2,m=opt%rel_m)
! Set model tau, and set yelmo relaxation switch (4: gl line and grounding zone relaxing; 0: no relaxation)
yelmo1%tpo%par%topo_rel_tau = opt%rel_tau
yelmo1%tpo%par%topo_rel = 4
else
! Turn-off relaxation now
yelmo1%tpo%par%topo_rel = 0
end if
! === Optimization update step =========
if (opt%opt_cf .and. &
(time_elapsed .ge. opt%cf_time_init .and. time_elapsed .le. opt%cf_time_end) ) then
! Perform cf_ref optimization
! Update cb_ref based on error metric(s)
call optimize_cb_ref(yelmo1%dyn%now%cb_ref,yelmo1%tpo%now%H_ice, &
yelmo1%tpo%now%dHidt,yelmo1%bnd%z_bed,yelmo1%bnd%z_sl,yelmo1%dyn%now%ux_s,yelmo1%dyn%now%uy_s, &
yelmo1%dta%pd%H_ice,yelmo1%dta%pd%uxy_s,yelmo1%dta%pd%H_grnd, &
opt%cf_min,opt%cf_max,yelmo1%tpo%par%dx,opt%sigma_err,opt%sigma_vel,opt%tau_c,opt%H0, &
dt=ctl%dtt,fill_method=opt%fill_method,fill_dist=opt%sigma_err, &
cb_tgt=yelmo1%dyn%now%cb_tgt)
end if
if (opt%opt_tf .and. &
(time_elapsed .ge. opt%tf_time_init .and. time_elapsed .le. opt%tf_time_end) ) then
! Perform tf_corr optimization
call optimize_tf_corr(mshlf1%now%tf_corr,yelmo1%tpo%now%H_ice,yelmo1%tpo%now%H_grnd,yelmo1%tpo%now%dHidt, &
yelmo1%dta%pd%H_ice,yelmo1%dta%pd%H_grnd,opt%H_grnd_lim,opt%tau_m,opt%m_temp, &
opt%tf_min,opt%tf_max,yelmo1%tpo%par%dx,sigma=opt%tf_sigma,dt=ctl%dtt)
! call optimize_tf_corr(mshlf1%now%tf_corr,yelmo1%tpo%now%H_ice,yelmo1%tpo%now%H_grnd,yelmo1%tpo%now%dHidt, &
! yelmo1%dta%pd%H_ice,yelmo1%bnd%basins,opt%H_grnd_lim, &
! opt%tau_m,opt%m_temp,opt%tf_min,opt%tf_max,opt%tf_basins,dt=ctl%dtt)
end if
case("relax")
! ===== relaxation spinup ==================
if (time_elapsed .lt. ctl%time_equil) then
! Turn on relaxation for now, to let thermodynamics equilibrate
! without changing the topography too much. Important when
! effective pressure = f(thermodynamics).
yelmo1%tpo%par%topo_rel = 3
yelmo1%tpo%par%topo_rel_tau = 50.0
write(*,*) "timelog, tau = ", yelmo1%tpo%par%topo_rel_tau
else if (time_elapsed .eq. ctl%time_equil) then
! Disable relaxation now...
yelmo1%tpo%par%topo_rel = 0
write(*,*) "timelog, relation off..."
end if
case DEFAULT ! == "none", etc
! Pass - do nothing
end select
call timer_step(tmrs,comp=0)
! == ISOSTASY and SEA LEVEL ======================================================
call bsl_update(bsl, time_bp)
call isos_update(isos1, yelmo1%tpo%now%H_ice, time, bsl, dwdt_corr=yelmo1%bnd%dzbdt_corr)
yelmo1%bnd%z_bed = isos1%out%z_bed
yelmo1%bnd%z_sl = isos1%out%z_ss
call timer_step(tmrs,comp=1,time_mod=[time-ctl%dtt,time]*1e-3,label="isostasy")
! == ICE SHEET ===================================================
if (running_greenland .and. ngs%use_negis_par) then
! Update cb_ref using negis parameters
call negis_update_cb_ref(yelmo1,ngs,time)
end if
! Update Yelmo
if (ctl%with_ice_sheet .and. (.not. (n .eq. 0 .and. yelmo1%par%use_restart)) ) then
call yelmo_update(yelmo1,time)
end if
call timer_step(tmrs,comp=2,time_mod=[time-ctl%dtt,time]*1e-3,label="yelmo")
! == CLIMATE (ATMOSPHERE AND OCEAN) ====================================
if (mod(nint(time*100),nint(ctl%dt_clim*100))==0) then
! Update snapclim
call snapclim_update(snp1,z_srf=yelmo1%tpo%now%z_srf,time=time_bp,domain=domain,dx=yelmo1%grd%dx,basins=yelmo1%bnd%basins)
end if
! == SURFACE MASS BALANCE ==============================================
! ajr: just testing...
dT_now = 0.0
!if (time .gt. 7000.0) dT_now = 4.0
call smbpal_update_monthly(smbpal1,snp1%now%tas,snp1%now%pr, &
yelmo1%tpo%now%z_srf,yelmo1%tpo%now%H_ice,time_bp)
yelmo1%bnd%smb = smbpal1%ann%smb*yelmo1%bnd%c%conv_we_ie*1e-3 ! [mm we/a] => [m ie/a]
yelmo1%bnd%T_srf = smbpal1%ann%tsrf
if (trim(yelmo1%par%domain) .eq. "Greenland" .and. scale_glacial_smb) then
! Modify glacial smb
call calc_glacial_smb(yelmo1%bnd%smb,yelmo1%grd%lat,snp1%now%ta_ann,snp1%clim0%ta_ann)
end if
! yelmo1%bnd%smb = yelmo1%dta%pd%smb
! yelmo1%bnd%T_srf = yelmo1%dta%pd%t2m
if (running_laurentide .and. laurentide_init_const_H &
.and. (time-ctl%time_init) .lt. laurentide_time_equil ) then
! Additionally ensure smb is postive for land above 50degN in Laurentide region
! to make sure ice grows everywhere needed (Coridilleran ice sheet mainly)
where (yelmo1%bnd%regions .eq. 1.1 .and. yelmo1%grd%lat .gt. 50.0 .and. &
yelmo1%bnd%z_bed .gt. 0.0 .and. yelmo1%bnd%smb .lt. 0.0 ) yelmo1%bnd%smb = 0.5
end if
! == MARINE AND TOTAL BASAL MASS BALANCE ===============================
call marshelf_update_shelf(mshlf1,yelmo1%tpo%now%H_ice,yelmo1%bnd%z_bed,yelmo1%tpo%now%f_grnd, &
yelmo1%bnd%basins,yelmo1%bnd%z_sl,yelmo1%grd%dx,snp1%now%depth, &
snp1%now%to_ann,snp1%now%so_ann,dto_ann=snp1%now%to_ann-snp1%clim0%to_ann)
call marshelf_update(mshlf1,yelmo1%tpo%now%H_ice,yelmo1%bnd%z_bed,yelmo1%tpo%now%f_grnd, &
yelmo1%bnd%regions,yelmo1%bnd%basins,yelmo1%bnd%z_sl,dx=yelmo1%grd%dx)
yelmo1%bnd%bmb_shlf = mshlf1%now%bmb_shlf
yelmo1%bnd%T_shlf = mshlf1%now%T_shlf
call timer_step(tmrs,comp=3,time_mod=[time-ctl%dtt,time]*1e-3,label="climate")
! == MODEL OUTPUT =======================================================
if (timeout_check(tm_2D,time)) then
call yelmox_write_step(yelmo1,snp1,mshlf1,smbpal1,file2D,time=time)
end if
if (timeout_check(tm_2Dsm,time)) then
call yelmo_write_step(yelmo1,file2D_small,time,compare_pd=.FALSE.)
end if
if (timeout_check(tm_1D,time)) then
call yelmo_write_reg_step(yelmo1,file1D,time=time)
if (reg1%write) then
call yelmo_write_reg_step(yelmo1,reg1%fnm,time=time,mask=reg1%mask)
end if
if (reg2%write) then
call yelmo_write_reg_step(yelmo1,reg2%fnm,time=time,mask=reg2%mask)
end if
if (reg3%write) then
call yelmo_write_reg_step(yelmo1,reg3%fnm,time=time,mask=reg3%mask)
end if
end if
if (mod(nint(time*100),nint(ctl%dt_restart*100))==0) then
call yelmox_restart_write(bsl,isos1,yelmo1,time)
end if
call timer_step(tmrs,comp=4,time_mod=[time-ctl%dtt,time]*1e-3,label="io")
if (mod(time_elapsed,10.0)==0) then
! Print timestep timing info and write log table
call timer_write_table(tmrs,[time,ctl%dtt]*1e-3,"m",tmr_file,init=time_elapsed .eq. 0.0)
end if
if (mod(time,10.0)==0 .and. (.not. yelmo_log)) then
write(*,"(a,f14.4)") "yelmo:: time = ", time
end if
end do
! Stop timing
call timer_step(tmr,comp=2,time_mod=[ctl%time_init,time]*1e-3,label="timeloop")
! Write the restart snapshot for the end of the simulation
call yelmox_restart_write(bsl,isos1,yelmo1,time)
! Finalize program
call yelmo_end(yelmo1,time=time)
! Print timing summary
call timer_print_summary(tmr,units="m",units_mod="kyr",time_mod=time*1e-3)
contains
subroutine calc_glacial_smb(smb,lat2D,ta_ann,ta_ann_pd)
implicit none
real(wp), intent(INOUT) :: smb(:,:)
real(wp), intent(IN) :: lat2D(:,:)
real(wp), intent(IN) :: ta_ann(:,:)
real(wp), intent(IN) :: ta_ann_pd(:,:)
! Local variables
integer :: i, j, nx, ny
real(wp) :: t0, tnow
real(wp) :: at
real(wp) :: fac
real(wp), parameter :: dt_lgm = -8.0
real(wp), parameter :: lat_lim = 55.0
real(wp), parameter :: fac_lim = 0.90
nx = size(smb,1)
ny = size(smb,2)
! Determine a quasi glacial-interglacial index
! 0: interglacial
! 1: glacial
tnow = sum(ta_ann) / real(nx*ny,wp)
t0 = sum(ta_ann_pd) / real(nx*ny,wp)
at = (tnow-t0)/dt_lgm
if (at .lt. 0.0) at = 0.0
if (at .gt. 1.0) at = 1.0
! Now determine smb scaling as a function of glacial index and latitude
! fac==0: no change to smb
! 0<fac<=1: scale smb up to value of fac_lim
do j = 1, ny
do i = 1, nx
if (smb(i,j) .lt. 0.0 .and. lat2D(i,j) .gt. lat_lim) then
! Calculate scaling here
smb(i,j) = smb(i,j) - smb(i,j) * at * fac_lim
end if
end do
end do
return
end subroutine calc_glacial_smb
! ===== NEGIS routines ============================
subroutine negis_update_cb_ref(ylmo,ngs,time)
implicit none
type(yelmo_class), intent(INOUT) :: ylmo
type(negis_params), intent(INOUT) :: ngs
real(wp), intent(IN) :: time
! Local variables
integer :: i, j, nx, ny
nx = ylmo%grd%nx
ny = ylmo%grd%ny
if (time .lt. -11e3) then
ngs%cf_x = ngs%cf_0
else
! Linear interpolation from cf_0 to cf_1
ngs%cf_x = ngs%cf_0 + (time - (-11e3))/ ((0.0) - (-11e3)) * (ngs%cf_1 - ngs%cf_0)
end if
if (time .lt. -4e3) then
ngs%cf_south = 1.0
else
ngs%cf_north = 1.0
end if
! Update bed roughness coefficients cb_ref and c_bed (which are independent of velocity)
! like normal, using the default function defined in Yelmo:
call calc_cb_ref(ylmo%dyn%now%cb_ref,ylmo%bnd%z_bed,ylmo%bnd%z_bed_sd,ylmo%bnd%z_sl, &
ylmo%bnd%H_sed,ylmo%dyn%par%till_f_sed,ylmo%dyn%par%till_sed_min, &
ylmo%dyn%par%till_sed_max,ylmo%dyn%par%till_cf_ref,ylmo%dyn%par%till_cf_min, &
ylmo%dyn%par%till_z0,ylmo%dyn%par%till_z1,ylmo%dyn%par%till_n_sd, &
ylmo%dyn%par%till_scale,ylmo%dyn%par%till_method)
! === Finally, apply NEGIS scaling =============================
do j = 1, ny
do i = 1, nx
if(ylmo%bnd%basins(i,j) .eq. 9.1) yelmo1%dyn%now%cb_ref(i,j) = yelmo1%dyn%now%cb_ref(i,j) * ngs%cf_centre
if(ylmo%bnd%basins(i,j) .eq. 9.2) yelmo1%dyn%now%cb_ref(i,j) = yelmo1%dyn%now%cb_ref(i,j) * ngs%cf_south
if(ylmo%bnd%basins(i,j) .eq. 9.3) yelmo1%dyn%now%cb_ref(i,j) = yelmo1%dyn%now%cb_ref(i,j) * ngs%cf_north
end do
end do
return
end subroutine negis_update_cb_ref
subroutine yelmox_write_step(ylmo,snp,mshlf,srf,filename,time)
implicit none
type(yelmo_class), intent(IN) :: ylmo
type(snapclim_class), intent(IN) :: snp
type(marshelf_class), intent(IN) :: mshlf
type(smbpal_class), intent(IN) :: srf
!type(sediments_class), intent(IN) :: sed
!type(geothermal_class), intent(IN) :: gthrm
!type(isos_class), intent(IN) :: isos
character(len=*), intent(IN) :: filename
real(wp), intent(IN) :: time