-
Notifications
You must be signed in to change notification settings - Fork 135
/
Copy pathicepack_therm_itd.F90
2289 lines (1860 loc) · 88.9 KB
/
icepack_therm_itd.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
!=======================================================================
!
! Thermo calculations after call to coupler, related to ITD:
! ice thickness redistribution, lateral growth and melting.
!
! NOTE: The thermodynamic calculation is split in two for load balancing.
! First icepack_therm_vertical computes vertical growth rates and coupler
! fluxes. Then icepack_therm_itd does thermodynamic calculations not
! needed for coupling.
!
! authors William H. Lipscomb, LANL
! C. M. Bitz, UW
! Elizabeth C. Hunke, LANL
!
! 2003: Vectorized by Clifford Chen (Fujitsu) and William Lipscomb
! 2004: Block structure added by William Lipscomb
! 2006: Streamlined for efficiency by Elizabeth Hunke
! 2014: Column package created by Elizabeth Hunke
!
module icepack_therm_itd
use icepack_kinds
use icepack_parameters, only: c0, c1, c2, c3, c4, c6, c10
use icepack_parameters, only: p001, p1, p333, p5, p666, puny, bignum
use icepack_parameters, only: rhos, rhoi, Lfresh, ice_ref_salinity
use icepack_parameters, only: phi_init, dsin0_frazil, hs_ssl, salt_loss
use icepack_parameters, only: rhosi, conserv_check, rhosmin
use icepack_parameters, only: kitd, ktherm, heat_capacity
use icepack_parameters, only: z_tracers, solve_zsal, hfrazilmin
use icepack_tracers, only: ntrcr, nbtrcr
use icepack_tracers, only: nt_qice, nt_qsno, nt_fbri, nt_sice
use icepack_tracers, only: nt_apnd, nt_hpnd, nt_aero, nt_isosno, nt_isoice
use icepack_tracers, only: nt_Tsfc, nt_iage, nt_FY, nt_fsd, nt_rhos
use icepack_tracers, only: nt_alvl, nt_vlvl
use icepack_tracers, only: tr_pond_cesm, tr_pond_lvl, tr_pond_topo, tr_snow
use icepack_tracers, only: tr_iage, tr_FY, tr_lvl, tr_aero, tr_iso, tr_brine, tr_fsd
use icepack_tracers, only: n_aero, n_iso
use icepack_tracers, only: bio_index
use icepack_warnings, only: warnstr, icepack_warnings_add
use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
use icepack_fsd, only: fsd_weld_thermo, icepack_cleanup_fsd, get_subdt_fsd
use icepack_itd, only: reduce_area, cleanup_itd
use icepack_itd, only: aggregate_area, shift_ice
use icepack_itd, only: column_sum, column_conservation_check
use icepack_isotope, only: isoice_alpha, isotope_frac_method
use icepack_mushy_physics, only: liquidus_temperature_mush, enthalpy_mush
use icepack_therm_shared, only: hi_min
use icepack_zbgc, only: add_new_ice_bgc
use icepack_zbgc, only: lateral_melt_bgc
implicit none
private
public :: linear_itd, &
add_new_ice, &
lateral_melt, &
icepack_step_therm2
!=======================================================================
contains
!=======================================================================
!
! ITD scheme that shifts ice among categories
!
! See Lipscomb, W. H. Remapping the thickness distribution in sea
! ice models. 2001, J. Geophys. Res., Vol 106, 13989--14000.
!
! Using the thermodynamic "velocities", interpolate to find the
! velocities in thickness space at the category boundaries, and
! compute the new locations of the boundaries. Then for each
! category, compute the thickness distribution function, g(h),
! between hL and hR, the left and right boundaries of the category.
! Assume g(h) is a linear polynomial that satisfies two conditions:
!
! (1) The ice area implied by g(h) equals aicen(n).
! (2) The ice volume implied by g(h) equals aicen(n)*hicen(n).
!
! Given g(h), at each boundary compute the ice area and volume lying
! between the original and new boundary locations. Transfer area
! and volume across each boundary in the appropriate direction, thus
! restoring the original boundaries.
!
! authors: William H. Lipscomb, LANL
! Elizabeth C. Hunke, LANL
subroutine linear_itd (ncat, hin_max, &
nilyr, nslyr, &
ntrcr, trcr_depend, &
trcr_base, n_trcr_strata,&
nt_strata, &
aicen_init, vicen_init, &
aicen, trcrn, &
vicen, vsnon, &
aice, aice0, &
fpond )
integer (kind=int_kind), intent(in) :: &
ncat , & ! number of thickness categories
nilyr , & ! number of ice layers
nslyr , & ! number of snow layers
ntrcr ! number of tracers in use
real (kind=dbl_kind), dimension(0:ncat), intent(inout) :: &
hin_max ! category boundaries (m)
integer (kind=int_kind), dimension (:), intent(in) :: &
trcr_depend, & ! = 0 for aicen tracers, 1 for vicen, 2 for vsnon
n_trcr_strata ! number of underlying tracer layers
real (kind=dbl_kind), dimension (:,:), intent(in) :: &
trcr_base ! = 0 or 1 depending on tracer dependency
! argument 2: (1) aice, (2) vice, (3) vsno
integer (kind=int_kind), dimension (:,:), intent(in) :: &
nt_strata ! indices of underlying tracer layers
real (kind=dbl_kind), dimension(:), intent(in) :: &
aicen_init, & ! initial ice concentration (before vertical thermo)
vicen_init ! initial ice volume (m)
real (kind=dbl_kind), dimension (:), intent(inout) :: &
aicen , & ! ice concentration
vicen , & ! volume per unit area of ice (m)
vsnon ! volume per unit area of snow (m)
real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
trcrn ! ice tracers
real (kind=dbl_kind), intent(inout) :: &
aice , & ! concentration of ice
aice0 , & ! concentration of open water
fpond ! fresh water flux to ponds (kg/m^2/s)
! local variables
integer (kind=int_kind) :: &
n, nd , & ! category indices
k ! ice layer index
real (kind=dbl_kind) :: &
slope , & ! rate of change of dhice with hice
dh0 , & ! change in ice thickness at h = 0
da0 , & ! area melting from category 1
damax , & ! max allowed reduction in category 1 area
etamin, etamax,& ! left and right limits of integration
x1 , & ! etamax - etamin
x2 , & ! (etamax^2 - etamin^2) / 2
x3 , & ! (etamax^3 - etamin^3) / 3
wk1, wk2 ! temporary variables
real (kind=dbl_kind), dimension(0:ncat) :: &
hbnew ! new boundary locations
real (kind=dbl_kind), dimension(ncat) :: &
g0 , & ! constant coefficient in g(h)
g1 , & ! linear coefficient in g(h)
hL , & ! left end of range over which g(h) > 0
hR ! right end of range over which g(h) > 0
real (kind=dbl_kind), dimension(ncat) :: &
hicen , & ! ice thickness for each cat (m)
hicen_init , & ! initial ice thickness for each cat (pre-thermo)
dhicen , & ! thickness change for remapping (m)
daice , & ! ice area transferred across boundary
dvice ! ice volume transferred across boundary
real (kind=dbl_kind), dimension (ncat) :: &
eicen, & ! energy of melting for each ice layer (J/m^2)
esnon, & ! energy of melting for each snow layer (J/m^2)
vbrin, & ! ice volume defined by brine height (m)
sicen ! Bulk salt in h ice (ppt*m)
real (kind=dbl_kind) :: &
vice_init, vice_final, & ! ice volume summed over categories
vsno_init, vsno_final, & ! snow volume summed over categories
eice_init, eice_final, & ! ice energy summed over categories
esno_init, esno_final, & ! snow energy summed over categories
sice_init, sice_final, & ! ice bulk salinity summed over categories
vbri_init, vbri_final ! briny ice volume summed over categories
! NOTE: Third index of donor, daice, dvice should be ncat-1,
! except that compilers would have trouble when ncat = 1
integer (kind=int_kind), dimension(ncat) :: &
donor ! donor category index
logical (kind=log_kind) :: &
remap_flag ! remap ITD if remap_flag is true
character (len=char_len) :: &
fieldid ! field identifier
logical (kind=log_kind), parameter :: &
print_diags = .false. ! if true, prints when remap_flag=F
character(len=*),parameter :: subname='(linear_itd)'
!-----------------------------------------------------------------
! Initialize
!-----------------------------------------------------------------
hin_max(ncat) = 999.9_dbl_kind ! arbitrary big number
do n = 1, ncat
donor(n) = 0
daice(n) = c0
dvice(n) = c0
enddo
!-----------------------------------------------------------------
! Compute volume and energy sums that linear remapping should
! conserve.
!-----------------------------------------------------------------
if (conserv_check) then
do n = 1, ncat
eicen(n) = c0
esnon(n) = c0
vbrin(n) = c0
sicen(n) = c0
do k = 1, nilyr
eicen(n) = eicen(n) + trcrn(nt_qice+k-1,n) &
* vicen(n)/real(nilyr,kind=dbl_kind)
enddo
do k = 1, nslyr
esnon(n) = esnon(n) + trcrn(nt_qsno+k-1,n) &
* vsnon(n)/real(nslyr,kind=dbl_kind)
enddo
if (tr_brine) then
vbrin(n) = vbrin(n) + trcrn(nt_fbri,n) &
* vicen(n)/real(nilyr,kind=dbl_kind)
endif
do k = 1, nilyr
sicen(n) = sicen(n) + trcrn(nt_sice+k-1,n) &
* vicen(n)/real(nilyr,kind=dbl_kind)
enddo
enddo ! n
call column_sum (ncat, vicen, vice_init)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, vsnon, vsno_init)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, eicen, eice_init)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, esnon, esno_init)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, sicen, sice_init)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, vbrin, vbri_init)
if (icepack_warnings_aborted(subname)) return
endif ! conserv_check
!-----------------------------------------------------------------
! Initialize remapping flag.
! Remapping is done wherever remap_flag = .true.
! In rare cases the category boundaries may shift too far for the
! remapping algorithm to work, and remap_flag is set to .false.
! In these cases the simpler 'rebin' subroutine will shift ice
! between categories if needed.
!-----------------------------------------------------------------
remap_flag = .true.
!-----------------------------------------------------------------
! Compute thickness change in each category.
!-----------------------------------------------------------------
do n = 1, ncat
if (aicen_init(n) > puny) then
hicen_init(n) = vicen_init(n) / aicen_init(n)
else
hicen_init(n) = c0
endif ! aicen_init > puny
if (aicen (n) > puny) then
hicen (n) = vicen(n) / aicen(n)
dhicen(n) = hicen(n) - hicen_init(n)
else
hicen (n) = c0
dhicen(n) = c0
endif ! aicen > puny
enddo ! n
!-----------------------------------------------------------------
! Compute new category boundaries, hbnew, based on changes in
! ice thickness from vertical thermodynamics.
!-----------------------------------------------------------------
hbnew(0) = hin_max(0)
do n = 1, ncat-1
if (hicen_init(n) > puny .and. &
hicen_init(n+1) > puny) then
if ((hicen_init(n+1) - hicen_init(n))>0) then
! interpolate between adjacent category growth rates
slope = (dhicen(n+1) - dhicen(n)) / &
(hicen_init(n+1) - hicen_init(n))
hbnew(n) = hin_max(n) + dhicen(n) &
+ slope * (hin_max(n) - hicen_init(n))
else
write(warnstr,*) subname, &
'ITD Thermodynamics: hicen_init(n+1) <= hicen_init(n)'
call icepack_warnings_setabort(.true.)
call icepack_warnings_add(warnstr)
endif
elseif (hicen_init(n) > puny) then ! hicen_init(n+1)=0
hbnew(n) = hin_max(n) + dhicen(n)
elseif (hicen_init(n+1) > puny) then ! hicen_init(n)=0
hbnew(n) = hin_max(n) + dhicen(n+1)
else
hbnew(n) = hin_max(n)
endif
!-----------------------------------------------------------------
! Check that each boundary lies between adjacent values of hicen.
! If not, set remap_flag = .false.
! Write diagnosis outputs if remap_flag was changed to false
!-----------------------------------------------------------------
if (aicen(n) > puny .and. hicen(n) >= hbnew(n)) then
remap_flag = .false.
if (print_diags) then
write(warnstr,*) subname, 'ITD: hicen(n) > hbnew(n)'
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'cat ',n
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hicen(n) =', hicen(n)
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hbnew(n) =', hbnew(n)
call icepack_warnings_add(warnstr)
endif
elseif (aicen(n+1) > puny .and. hicen(n+1) <= hbnew(n)) then
remap_flag = .false.
if (print_diags) then
write(warnstr,*) subname, 'ITD: hicen(n+1) < hbnew(n)'
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'cat ',n
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hicen(n+1) =', hicen(n+1)
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hbnew(n) =', hbnew(n)
call icepack_warnings_add(warnstr)
endif
endif
!-----------------------------------------------------------------
! Check that hbnew(n) lies between hin_max(n-1) and hin_max(n+1).
! If not, set remap_flag = .false.
! (In principle we could allow this, but it would make the code
! more complicated.)
! Write diagnosis outputs if remap_flag was changed to false
!-----------------------------------------------------------------
if (hbnew(n) > hin_max(n+1)) then
remap_flag = .false.
if (print_diags) then
write(warnstr,*) subname, 'ITD hbnew(n) > hin_max(n+1)'
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'cat ',n
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hbnew(n) =', hbnew(n)
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hin_max(n+1) =', hin_max(n+1)
call icepack_warnings_add(warnstr)
endif
endif
if (hbnew(n) < hin_max(n-1)) then
remap_flag = .false.
if (print_diags) then
write(warnstr,*) subname, 'ITD: hbnew(n) < hin_max(n-1)'
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'cat ',n
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hbnew(n) =', hbnew(n)
call icepack_warnings_add(warnstr)
write(warnstr,*) subname, 'hin_max(n-1) =', hin_max(n-1)
call icepack_warnings_add(warnstr)
endif
endif
enddo ! boundaries, 1 to ncat-1
!-----------------------------------------------------------------
! Compute hbnew(ncat)
!-----------------------------------------------------------------
if (aicen(ncat) > puny) then
hbnew(ncat) = c3*hicen(ncat) - c2*hbnew(ncat-1)
else
hbnew(ncat) = hin_max(ncat)
endif
hbnew(ncat) = max(hbnew(ncat),hin_max(ncat-1))
!-----------------------------------------------------------------
! Identify cells where the ITD is to be remapped
!-----------------------------------------------------------------
if (remap_flag) then
!-----------------------------------------------------------------
! Compute g(h) for category 1 at start of time step
! (hicen = hicen_init)
!-----------------------------------------------------------------
call fit_line(aicen(1), hicen_init(1), &
hbnew(0), hin_max (1), &
g0 (1), g1 (1), &
hL (1), hR (1))
if (icepack_warnings_aborted(subname)) return
!-----------------------------------------------------------------
! Find area lost due to melting of thin (category 1) ice
!-----------------------------------------------------------------
if (aicen(1) > puny) then
dh0 = dhicen(1)
if (dh0 < c0) then ! remove area from category 1
dh0 = min(-dh0,hin_max(1)) ! dh0 --> |dh0|
!-----------------------------------------------------------------
! Integrate g(1) from 0 to dh0 to estimate area melted
!-----------------------------------------------------------------
! right integration limit (left limit = 0)
etamax = min(dh0,hR(1)) - hL(1)
if (etamax > c0) then
x1 = etamax
x2 = p5 * etamax*etamax
da0 = g1(1)*x2 + g0(1)*x1 ! ice area removed
! constrain new thickness <= hicen_init
damax = aicen(1) * (c1-hicen(1)/hicen_init(1)) ! damax > 0
da0 = min (da0, damax)
! remove area, conserving volume
hicen(1) = hicen(1) * aicen(1) / (aicen(1)-da0)
aicen(1) = aicen(1) - da0
if (tr_pond_topo) &
fpond = fpond - (da0 * trcrn(nt_apnd,1) &
* trcrn(nt_hpnd,1))
endif ! etamax > 0
else ! dh0 >= 0
hbnew(0) = min(dh0,hin_max(1)) ! shift hbnew(0) to right
endif
endif ! aicen(1) > puny
!-----------------------------------------------------------------
! Compute g(h) for each ice thickness category.
!-----------------------------------------------------------------
do n = 1, ncat
call fit_line(aicen(n), hicen(n), &
hbnew(n-1), hbnew(n), &
g0 (n), g1 (n), &
hL (n), hR (n))
if (icepack_warnings_aborted(subname)) return
!-----------------------------------------------------------------
! Compute area and volume to be shifted across each boundary.
!-----------------------------------------------------------------
donor(n) = 0
daice(n) = c0
dvice(n) = c0
enddo
do n = 1, ncat-1
if (hbnew(n) > hin_max(n)) then ! transfer from n to n+1
! left and right integration limits in eta space
etamin = max(hin_max(n), hL(n)) - hL(n)
etamax = min(hbnew(n), hR(n)) - hL(n)
donor(n) = n
else ! hbnew(n) <= hin_max(n); transfer from n+1 to n
! left and right integration limits in eta space
etamin = c0
etamax = min(hin_max(n), hR(n+1)) - hL(n+1)
donor(n) = n+1
endif ! hbnew(n) > hin_max(n)
if (etamax > etamin) then
x1 = etamax - etamin
wk1 = etamin*etamin
wk2 = etamax*etamax
x2 = p5 * (wk2 - wk1)
wk1 = wk1*etamin
wk2 = wk2*etamax
x3 = p333 * (wk2 - wk1)
nd = donor(n)
daice(n) = g1(nd)*x2 + g0(nd)*x1
dvice(n) = g1(nd)*x3 + g0(nd)*x2 + daice(n)*hL(nd)
endif ! etamax > etamin
! If daice or dvice is very small, shift no ice.
nd = donor(n)
if (daice(n) < aicen(nd)*puny) then
daice(n) = c0
dvice(n) = c0
donor(n) = 0
endif
if (dvice(n) < vicen(nd)*puny) then
daice(n) = c0
dvice(n) = c0
donor(n) = 0
endif
! If daice is close to aicen or dvice is close to vicen,
! shift entire category
if (daice(n) > aicen(nd)*(c1-puny)) then
daice(n) = aicen(nd)
dvice(n) = vicen(nd)
endif
if (dvice(n) > vicen(nd)*(c1-puny)) then
daice(n) = aicen(nd)
dvice(n) = vicen(nd)
endif
enddo ! boundaries, 1 to ncat-1
!-----------------------------------------------------------------
! Shift ice between categories as necessary
!-----------------------------------------------------------------
! maintain qsno negative definiteness
do n = 1, ncat
do k = nt_qsno, nt_qsno+nslyr-1
trcrn(k,n) = trcrn(k,n) + rhos*Lfresh
enddo
enddo
! maintain rhos_cmp positive definiteness
if (tr_snow) then
do n = 1, ncat
do k = nt_rhos, nt_rhos+nslyr-1
trcrn(k,n) = max(trcrn(k,n)-rhosmin, c0)
! trcrn(k,n) = trcrn(k,n) - rhosmin
enddo
enddo
endif
call shift_ice (ntrcr, ncat, &
trcr_depend, &
trcr_base, &
n_trcr_strata, &
nt_strata, &
aicen, trcrn, &
vicen, vsnon, &
hicen, donor, &
daice, dvice )
if (icepack_warnings_aborted(subname)) return
! maintain qsno negative definiteness
do n = 1, ncat
do k = nt_qsno, nt_qsno+nslyr-1
trcrn(k,n) = trcrn(k,n) - rhos*Lfresh
enddo
enddo
! maintain rhos_cmp positive definiteness
if (tr_snow) then
do n = 1, ncat
do k = nt_rhos, nt_rhos+nslyr-1
trcrn(k,n) = trcrn(k,n) + rhosmin
enddo
enddo
endif
!-----------------------------------------------------------------
! Make sure hice(1) >= minimum ice thickness hi_min.
!-----------------------------------------------------------------
if (hi_min > c0 .and. aicen(1) > puny .and. hicen(1) < hi_min) then
da0 = aicen(1) * (c1 - hicen(1)/hi_min)
aicen(1) = aicen(1) - da0
hicen(1) = hi_min
if (tr_pond_topo) &
fpond = fpond - (da0 * trcrn(nt_apnd,1) &
* trcrn(nt_hpnd,1))
endif
endif ! remap_flag
!-----------------------------------------------------------------
! Update fractional ice area in each grid cell.
!-----------------------------------------------------------------
call aggregate_area (ncat, aicen, aice, aice0)
if (icepack_warnings_aborted(subname)) return
!-----------------------------------------------------------------
! Check volume and energy conservation.
!-----------------------------------------------------------------
if (conserv_check) then
do n = 1, ncat
eicen(n) = c0
esnon(n) = c0
vbrin(n) = c0
sicen(n) = c0
do k = 1, nilyr
eicen(n) = eicen(n) + trcrn(nt_qice+k-1,n) &
* vicen(n)/real(nilyr,kind=dbl_kind)
enddo
do k = 1, nslyr
esnon(n) = esnon(n) + trcrn(nt_qsno+k-1,n) &
* vsnon(n)/real(nslyr,kind=dbl_kind)
enddo
if (tr_brine) then
vbrin(n) = vbrin(n) + trcrn(nt_fbri,n) &
* vicen(n)/real(nilyr,kind=dbl_kind)
endif
do k = 1, nilyr
sicen(n) = sicen(n) + trcrn(nt_sice+k-1,n) &
* vicen(n)/real(nilyr,kind=dbl_kind)
enddo
enddo ! n
call column_sum (ncat, vicen, vice_final)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, vsnon, vsno_final)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, eicen, eice_final)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, esnon, esno_final)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, sicen, sice_final)
if (icepack_warnings_aborted(subname)) return
call column_sum (ncat, vbrin, vbri_final)
if (icepack_warnings_aborted(subname)) return
fieldid = subname//':vice'
call column_conservation_check (fieldid, &
vice_init, vice_final, &
puny)
if (icepack_warnings_aborted(subname)) return
fieldid = subname//':vsno'
call column_conservation_check (fieldid, &
vsno_init, vsno_final, &
puny)
if (icepack_warnings_aborted(subname)) return
fieldid = subname//':eice'
call column_conservation_check (fieldid, &
eice_init, eice_final, &
puny*Lfresh*rhoi)
if (icepack_warnings_aborted(subname)) return
fieldid = subname//':esno'
call column_conservation_check (fieldid, &
esno_init, esno_final, &
puny*Lfresh*rhos)
if (icepack_warnings_aborted(subname)) return
fieldid = subname//':sicen'
call column_conservation_check (fieldid, &
sice_init, sice_final, &
puny)
if (icepack_warnings_aborted(subname)) return
fieldid = subname//':vbrin'
call column_conservation_check (fieldid, &
vbri_init, vbri_final, &
puny*c10)
if (icepack_warnings_aborted(subname)) return
endif ! conservation check
end subroutine linear_itd
!=======================================================================
!
! Fit g(h) with a line, satisfying area and volume constraints.
! To reduce roundoff errors caused by large values of g0 and g1,
! we actually compute g(eta), where eta = h - hL, and hL is the
! left boundary.
!
! authors: William H. Lipscomb, LANL
! Elizabeth C. Hunke, LANL
subroutine fit_line (aicen, hice, &
hbL, hbR, &
g0, g1, &
hL, hR)
real (kind=dbl_kind), intent(in) :: &
aicen ! concentration of ice
real (kind=dbl_kind), intent(in) :: &
hbL, hbR , & ! left and right category boundaries
hice ! ice thickness
real (kind=dbl_kind), intent(out):: &
g0, g1 , & ! coefficients in linear equation for g(eta)
hL , & ! min value of range over which g(h) > 0
hR ! max value of range over which g(h) > 0
! local variables
real (kind=dbl_kind) :: &
h13 , & ! hbL + 1/3 * (hbR - hbL)
h23 , & ! hbL + 2/3 * (hbR - hbL)
dhr , & ! 1 / (hR - hL)
wk1, wk2 ! temporary variables
character(len=*),parameter :: subname='(fit_line)'
!-----------------------------------------------------------------
! Compute g0, g1, hL, and hR for each category to be remapped.
!-----------------------------------------------------------------
if (aicen > puny .and. hbR - hbL > puny) then
! Initialize hL and hR
hL = hbL
hR = hbR
! Change hL or hR if hicen(n) falls outside central third of range
h13 = p333 * (c2*hL + hR)
h23 = p333 * (hL + c2*hR)
if (hice < h13) then
hR = c3*hice - c2*hL
elseif (hice > h23) then
hL = c3*hice - c2*hR
endif
! Compute coefficients of g(eta) = g0 + g1*eta
dhr = c1 / (hR - hL)
wk1 = c6 * aicen * dhr
wk2 = (hice - hL) * dhr
g0 = wk1 * (p666 - wk2)
g1 = c2*dhr * wk1 * (wk2 - p5)
else
g0 = c0
g1 = c0
hL = c0
hR = c0
endif ! aicen > puny
end subroutine fit_line
!=======================================================================
!
! Given some added new ice to the base of the existing ice, recalculate
! vertical tracer so that new grid cells are all the same size.
!
! author: A. K. Turner, LANL
!
subroutine update_vertical_tracers(nilyr, trc, h1, h2, trc0)
integer (kind=int_kind), intent(in) :: &
nilyr ! number of ice layers
real (kind=dbl_kind), dimension(:), intent(inout) :: &
trc ! vertical tracer
real (kind=dbl_kind), intent(in) :: &
h1, & ! old thickness
h2, & ! new thickness
trc0 ! tracer value of added ice on ice bottom
! local variables
real(kind=dbl_kind), dimension(nilyr) :: trc2 ! updated tracer temporary
! vertical indices for old and new grid
integer :: k1, k2
real (kind=dbl_kind) :: &
z1a, z1b, & ! upper, lower boundary of old cell/added new ice at bottom
z2a, z2b, & ! upper, lower boundary of new cell
overlap , & ! overlap between old and new cell
rnilyr
character(len=*),parameter :: subname='(update_vertical_tracers)'
rnilyr = real(nilyr,dbl_kind)
! loop over new grid cells
do k2 = 1, nilyr
! initialize new tracer
trc2(k2) = c0
! calculate upper and lower boundary of new cell
z2a = ((k2 - 1) * h2) / rnilyr
z2b = (k2 * h2) / rnilyr
! loop over old grid cells
do k1 = 1, nilyr
! calculate upper and lower boundary of old cell
z1a = ((k1 - 1) * h1) / rnilyr
z1b = (k1 * h1) / rnilyr
! calculate overlap between old and new cell
overlap = max(min(z1b, z2b) - max(z1a, z2a), c0)
! aggregate old grid cell contribution to new cell
trc2(k2) = trc2(k2) + overlap * trc(k1)
enddo ! k1
! calculate upper and lower boundary of added new ice at bottom
z1a = h1
z1b = h2
! calculate overlap between added ice and new cell
overlap = max(min(z1b, z2b) - max(z1a, z2a), c0)
! aggregate added ice contribution to new cell
trc2(k2) = trc2(k2) + overlap * trc0
! renormalize new grid cell
trc2(k2) = (rnilyr * trc2(k2)) / h2
enddo ! k2
! update vertical tracer array with the adjusted tracer
trc = trc2
end subroutine update_vertical_tracers
!=======================================================================
!
! Given the fraction of ice melting laterally in each grid cell
! (computed in subroutine frzmlt_bottom_lateral), melt ice.
!
! author: C. M. Bitz, UW
! 2003: Modified by William H. Lipscomb and Elizabeth C. Hunke, LANL
! 2016 Lettie Roach, NIWA/VUW, added floe size dependence
!
subroutine lateral_melt (dt, ncat, &
nilyr, nslyr, &
n_aero, &
fpond, &
fresh, fsalt, &
fhocn, faero_ocn, &
fiso_ocn, &
rside, meltl, &
fside, sss, &
aicen, vicen, &
vsnon, trcrn, &
fzsal, flux_bio, &
nbtrcr, nblyr, &
nfsd, d_afsd_latm,&
floe_rad_c, floe_binwidth)
real (kind=dbl_kind), intent(in) :: &
dt ! time step (s)
integer (kind=int_kind), intent(in) :: &
ncat , & ! number of thickness categories
nfsd , & ! number of floe size categories
nilyr , & ! number of ice layers
nblyr , & ! number of bio layers
nslyr , & ! number of snow layers
n_aero , & ! number of aerosol tracers
nbtrcr ! number of bio tracers
real (kind=dbl_kind), dimension (:), intent(inout) :: &
aicen , & ! concentration of ice
vicen , & ! volume per unit area of ice (m)
vsnon ! volume per unit area of snow (m)
real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
trcrn ! tracer array
real (kind=dbl_kind), intent(in) :: &
rside , & ! fraction of ice that melts laterally
fside ! lateral heat flux (W/m^2)
real (kind=dbl_kind), intent(inout) :: &
fpond , & ! fresh water flux to ponds (kg/m^2/s)
fresh , & ! fresh water flux to ocean (kg/m^2/s)
fsalt , & ! salt flux to ocean (kg/m^2/s)
fhocn , & ! net heat flux to ocean (W/m^2)
meltl , & ! lateral ice melt (m/step-->cm/day)
fzsal ! salt flux from zsalinity (kg/m2/s)
real (kind=dbl_kind), dimension (:), intent(in) :: &
floe_rad_c , & ! fsd size bin centre in m (radius)
floe_binwidth ! fsd size bin width in m (radius)
real (kind=dbl_kind), dimension (:), intent(out) :: &
d_afsd_latm ! change in fsd due to lateral melt (m)
real (kind=dbl_kind), dimension(nbtrcr), &
intent(inout) :: &
flux_bio ! biology tracer flux from layer bgc (mmol/m^2/s)
real (kind=dbl_kind), dimension(:), intent(inout) :: &
faero_ocn ! aerosol flux to ocean (kg/m^2/s)
real (kind=dbl_kind), dimension(:), intent(inout) :: &
fiso_ocn ! isotope flux to ocean (kg/m^2/s)
! local variables
integer (kind=int_kind) :: &
n , & ! thickness category index
k , & ! layer index
nsubt ! sub timesteps for FSD tendency
real (kind=dbl_kind) :: &
dfhocn , & ! change in fhocn
dfpond , & ! change in fpond
dfresh , & ! change in fresh
dfsalt , & ! change in fsalt
dvssl , & ! snow surface layer volume
dvint , & ! snow interior layer
cat1_arealoss, tmp !
logical (kind=log_kind) :: &
flag ! .true. if there could be lateral melting
real (kind=dbl_kind), dimension (ncat) :: &
aicen_init, & ! initial area fraction
vicen_init, & ! volume per unit area of ice (m)
G_radialn , & ! rate of lateral melt (m/s)
delta_an , & ! change in the ITD
qin , & ! enthalpy
rsiden ! delta_an/aicen
real (kind=dbl_kind), dimension (nfsd,ncat) :: &
afsdn , & ! floe size distribution tracer
afsdn_init ! initial value
real (kind=dbl_kind), dimension (nfsd) :: &
df_flx, & ! finite difference for FSD
afsd_tmp, d_afsd_tmp
real (kind=dbl_kind), dimension(nfsd+1) :: &
f_flx !
!echmod - for average qin
real (kind=dbl_kind), intent(in) :: &
sss
real (kind=dbl_kind) :: &
Ti, Si0, qi0, &
elapsed_t, & ! FSD subcycling
subdt ! FSD timestep (s)
character(len=*), parameter :: subname='(lateral_melt)'
flag = .false.
dfhocn = c0
dfpond = c0
dfresh = c0
dfsalt = c0
dvssl = c0
dvint = c0
cat1_arealoss = c0