forked from NOAA-GFDL/SIS2
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathSIS_fast_thermo.F90
1178 lines (1020 loc) · 56.2 KB
/
SIS_fast_thermo.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
!***********************************************************************
!* GNU General Public License *
!* This file is a part of SIS2. *
!* *
!* SIS2 is free software; you can redistribute it and/or modify it and *
!* are expected to follow the terms of the GNU General Public License *
!* as published by the Free Software Foundation; either version 2 of *
!* the License, or (at your option) any later version. *
!* *
!* SIS2 is distributed in the hope that it will be useful, but WITHOUT *
!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
!* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
!* License for more details. *
!* *
!* For the full text of the GNU General Public License, *
!* write to: Free Software Foundation, Inc., *
!* 675 Mass Ave, Cambridge, MA 02139, USA. *
!* or see: http://www.gnu.org/licenses/gpl.html *
!***********************************************************************
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
! SIS2 is a SEA ICE MODEL for coupling through the GFDL exchange grid. SIS2 !
! is a revision of the original SIS with have extended capabilities, including !
! the option of using a B-grid or C-grid spatial discretization. The SIS2 !
! software has been extensively reformulated from SIS for greater consistency !
! with the Modular Ocean Model, version 6 (MOM6), and to permit might tighter !
! dynamical coupling between the ocean and sea-ice. !
! This module handles the rapid thermodynamic interactions between the ice !
! and the atmosphere, including heating and the accumulation of fluxes, but !
! not changes to the ice or snow mass. !
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
module SIS_fast_thermo
use SIS_diag_mediator, only : SIS_diag_ctrl
! ! use SIS_diag_mediator, only : enable_SIS_averaging, disable_SIS_averaging
! ! use SIS_diag_mediator, only : query_SIS_averaging_enabled, post_SIS_data
! ! use SIS_diag_mediator, only : register_diag_field=>register_SIS_diag_field
use SIS_debugging, only : hchksum
use MOM_error_handler, only : SIS_error=>MOM_error, FATAL, WARNING, SIS_mesg=>MOM_mesg
use MOM_error_handler, only : callTree_enter, callTree_leave, callTree_waypoint
use MOM_file_parser, only : get_param, log_param, log_version, param_file_type
! use MOM_hor_index, only : hor_index_type, hor_index_init
! use MOM_obsolete_params, only : obsolete_logical
! use MOM_string_functions, only : uppercase
use MOM_time_manager, only : time_type, time_type_to_real
use MOM_time_manager, only : set_date, set_time, operator(+), operator(-)
use MOM_time_manager, only : operator(>), operator(*), operator(/), operator(/=)
use coupler_types_mod, only : coupler_3d_bc_type
use SIS_types, only : ice_state_type, IST_chksum, IST_bounds_check
use SIS_types, only : fast_ice_avg_type, ice_rad_type, simple_OSS_type, total_sfc_flux_type
use SIS_types, only : VIS_DIR, VIS_DIF, NIR_DIR, NIR_DIF
use ice_boundary_types, only : atmos_ice_boundary_type ! , land_ice_boundary_type
use SIS_hor_grid, only : SIS_hor_grid_type
use ice_grid, only : ice_grid_type
use SIS2_ice_thm, only : SIS2_ice_thm_CS, SIS2_ice_thm_init, SIS2_ice_thm_end
use SIS2_ice_thm, only : ice_temp_SIS2
use SIS2_ice_thm, only : get_SIS2_thermo_coefs, enth_from_TS, Temp_from_En_S
implicit none ; private
public :: do_update_ice_model_fast, SIS_fast_thermo_init, SIS_fast_thermo_end
public :: accumulate_deposition_fluxes, convert_frost_to_snow
public :: fast_thermo_CS, avg_top_quantities, total_top_quantities, infill_array
public :: redo_update_ice_model_fast, rescale_shortwave, find_excess_fluxes
type fast_thermo_CS ; private
! These two arrarys are used with column_check when evaluating the enthalpy
! conservation with the fast thermodynamics code.
real, pointer, dimension(:,:,:) :: &
enth_prev, heat_in
logical :: debug ! If true, write verbose checksums for debugging purposes.
logical :: column_check ! If true, enable the heat check column by column.
real :: imb_tol ! The tolerance for imbalances to be flagged by
! column_check, nondim.
logical :: bounds_check ! If true, check for sensible values of thicknesses
! temperatures, fluxes, etc.
integer :: n_fast = 0 ! The number of times update_ice_model_fast
! has been called.
logical :: Reorder_0C_heatflux=.false. ! If true, rearrange the calculation
! of the heat fluxes projected back to 0C to work
! on each contribution separately, so that they
! can be indentically replicated if there is
! a single fast timestep per coupled timestep and
! REDO_FAST_ICE_UPDATE=True
! These are pointers to the control structures for subsidiary modules.
type(SIS2_ice_thm_CS), pointer :: ice_thm_CSp => NULL()
end type fast_thermo_CS
contains
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> sum_top_quantities does a running sum of fluxes for later use by the slow ice
!! physics and the ocean. Nothing here will be exposed to other modules until
!! after it has passed through avg_top_quantities.
subroutine sum_top_quantities (FIA, ABT, flux_u, flux_v, flux_sh, evap, &
flux_sw, flux_lw, lprec, fprec, flux_lh, t_skin, SST, &
sh_T0, evap_T0, lw_T0, dshdt, devapdt, dlwdt, G, IG)
type(fast_ice_avg_type), intent(inout) :: FIA
type(atmos_ice_boundary_type), intent(in) :: ABT
type(SIS_hor_grid_type), intent(in) :: G
type(ice_grid_type), intent(in) :: IG
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce), intent(in) :: &
flux_u, &
flux_v, &
flux_sh, & ! The upward sensible heat flux from the top of the ice into
! the atmosphere in W m-2.
evap, & ! The upward flux of water due to sublimation or evaporation
! from the top of the ice to the atmosphere, in kg m-2 s-1.
flux_lw, & ! The net longwave heat flux from the atmosphere into the
! ice or ocean, in W m-2.
lprec, fprec, & ! The liquid and frozen precipitation onto the ice
! in kg m-2 s-1.
flux_lh, & ! The upward latent heat flux associated with sublimation or
! evaporation, in W m-2.
sh_T0, & ! The upward sensible heat flux from the top of the ice into
! the atmosphere when the skin temperature is 0 C, in W m-2.
evap_T0, & ! The sublimation rate when the skin temperature is 0 C,
! in kg m-2 s-1.
lw_T0, & ! The downward longwave heat flux from the atmosphere into the
! ice or ocean when the skin temperature is 0 C, in W m-2.
dshdt, & ! The derivative of the upward sensible heat flux from the
! the top of the ice into the atmosphere with ice skin
! temperature in W m-2 K-1.
devapdt, & ! The derivative of the sublimation rate with the surface
! temperature, in kg m-2 s-1 K-1.
dlwdt ! The derivative of the longwave heat flux from the atmosphere
! into the ice or ocean with ice skin temperature, in W m-2 K-1.
real, dimension(G%isd:G%ied,G%jsd:G%jed,IG%CatIce), intent(in) :: &
t_skin
real, dimension(G%isd:G%ied,G%jsd:G%jed), intent(in) :: &
SST
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce,size(FIA%flux_sw_top,4)), &
intent(in) :: flux_sw ! The downward shortwave heat fluxes in W m-2. The 4th
! dimension is a combination of angular orientation & frequency.
real :: t_sfc
integer :: i, j, k, m, n, b, nb, i2, j2, k2, isc, iec, jsc, jec, i_off, j_off, ncat
integer :: ind
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
nb = size(FIA%flux_sw_top,4)
i_off = LBOUND(ABT%t_flux,1) - G%isc
j_off = LBOUND(ABT%t_flux,2) - G%jsc
if (FIA%num_tr_fluxes < 0) then
! Determine how many atmospheric boundary fluxes have been passed in, and
! set up both an indexing array for these and a space to take their average.
! This code is only exercised the first time that sum_top_quantities is called.
FIA%num_tr_fluxes = 0
if (ABT%fluxes%num_bcs > 0) then
do n=1,ABT%fluxes%num_bcs
FIA%num_tr_fluxes = FIA%num_tr_fluxes + ABT%fluxes%bc(n)%num_fields
enddo
if (FIA%num_tr_fluxes > 0) then
allocate(FIA%tr_flux_top(G%isd:G%ied, G%jsd:G%jed, 0:IG%CatIce, FIA%num_tr_fluxes))
FIA%tr_flux_top(:,:,:,:) = 0.0
endif
endif
endif
if (FIA%avg_count == 0) then
! zero_top_quantities - zero fluxes to begin summing in ice fast physics.
FIA%flux_u_top(:,:,:) = 0.0 ; FIA%flux_v_top(:,:,:) = 0.0
FIA%flux_sh_top(:,:,:) = 0.0 ; FIA%evap_top(:,:,:) = 0.0
FIA%flux_lw_top(:,:,:) = 0.0 ; FIA%flux_lh_top(:,:,:) = 0.0
FIA%flux_sw_top(:,:,:,:) = 0.0
FIA%lprec_top(:,:,:) = 0.0 ; FIA%fprec_top(:,:,:) = 0.0
FIA%flux_sw_dn(:,:) = 0.0 ; FIA%Tskin_avg(:,:) = 0.0
if (allocated(FIA%flux_sh0)) then
FIA%dshdt(:,:,:) = 0.0 ; FIA%devapdt(:,:,:) = 0.0 ; FIA%dlwdt(:,:,:) = 0.0
FIA%flux_sh0(:,:,:) = 0.0 ; FIA%evap0(:,:,:) = 0.0
FIA%flux_lw0(:,:,:) = 0.0 ; FIA%Tskin_cat(:,:,:) = 0.0
endif
if (FIA%num_tr_fluxes > 0) FIA%tr_flux_top(:,:,:,:) = 0.0
endif
!$OMP parallel do default(shared)
do j=jsc,jec ; do k=0,ncat ; do i=isc,iec
FIA%flux_u_top(i,j,k) = FIA%flux_u_top(i,j,k) + flux_u(i,j,k)
FIA%flux_v_top(i,j,k) = FIA%flux_v_top(i,j,k) + flux_v(i,j,k)
FIA%flux_sh_top(i,j,k) = FIA%flux_sh_top(i,j,k) + flux_sh(i,j,k)
FIA%evap_top(i,j,k) = FIA%evap_top(i,j,k) + evap(i,j,k)
do b=1,nb ; FIA%flux_sw_top(i,j,k,b) = FIA%flux_sw_top(i,j,k,b) + flux_sw(i,j,k,b) ; enddo
FIA%flux_lw_top(i,j,k) = FIA%flux_lw_top(i,j,k) + flux_lw(i,j,k)
FIA%lprec_top(i,j,k) = FIA%lprec_top(i,j,k) + lprec(i,j,k)
FIA%fprec_top(i,j,k) = FIA%fprec_top(i,j,k) + fprec(i,j,k)
FIA%flux_lh_top(i,j,k) = FIA%flux_lh_top(i,j,k) + flux_lh(i,j,k)
enddo ; enddo ; enddo
! FIA%flux_sw_dn is accumulated where the fast radiation diagnostics are output
! because it depends on arrays that are stored in the public ice_data_type.
ind = 0
do n=1,ABT%fluxes%num_bcs ; do m=1,ABT%fluxes%bc(n)%num_fields
ind = ind + 1
!Do not handle air_sea_deposition fluxes here, they need to be handled after atmos_down
if(ABT%fluxes%bc(n)%flux_type /= 'air_sea_deposition') then
do k=0,ncat ; do j=jsc,jec ; do i=isc,iec
i2 = i+i_off ; j2 = j+j_off ; k2 = k+1
FIA%tr_flux_top(i,j,k,ind) = FIA%tr_flux_top(i,j,k,ind) + &
ABT%fluxes%bc(n)%field(m)%values(i2,j2,k2)
enddo ; enddo ; enddo
endif
enddo ; enddo
if (allocated(FIA%flux_sh0)) then
!$OMP parallel do default(shared) private(t_sfc)
do j=jsc,jec ; do k=0,ncat ; do i=isc,iec
FIA%dshdt(i,j,k) = FIA%dshdt(i,j,k) + dshdt(i,j,k)
FIA%devapdt(i,j,k) = FIA%devapdt(i,j,k) + devapdt(i,j,k)
FIA%dlwdt(i,j,k) = FIA%dlwdt(i,j,k) + dlwdt(i,j,k)
FIA%flux_sh0(i,j,k) = FIA%flux_sh0(i,j,k) + sh_T0(i,j,k)
FIA%evap0(i,j,k) = FIA%evap0(i,j,k) + evap_T0(i,j,k)
FIA%flux_lw0(i,j,k) = FIA%flux_lw0(i,j,k) + lw_T0(i,j,k)
t_sfc = SST(i,j) ; if (k>0) t_sfc = t_skin(i,j,k)
FIA%Tskin_cat(i,j,k) = FIA%Tskin_cat(i,j,k) + t_sfc
enddo ; enddo ; enddo
endif
FIA%avg_count = FIA%avg_count + 1
end subroutine sum_top_quantities
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> avg_top_quantities determines time average fluxes for later use by the
!! slow ice physics and by the ocean.
subroutine avg_top_quantities(FIA, Rad, IST, G, IG)
type(fast_ice_avg_type), intent(inout) :: FIA
type(ice_rad_type), intent(in) :: Rad
type(ice_state_type), intent(in) :: IST
type(SIS_hor_grid_type), intent(inout) :: G
type(ice_grid_type), intent(in) :: IG
real :: u, v, divid, sign
real :: I_avc ! The inverse of the number of contributions.
real :: I_wts ! 1.0 / ice_cover or 0 if ice_cover is 0, nondim.
integer :: i, j, k, m, n, b, nb, isc, iec, jsc, jec, ncat
integer :: isd, ied, jsd, jed
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
nb = size(FIA%flux_sw_top,4)
!
! compute average fluxes
!
if (FIA%avg_count == 0) call SIS_error(FATAL,'avg_top_quantities: '//&
'no ocean model fluxes have been averaged')
! Rotate the stress from lat/lon to ocean coordinates and possibly change the
! sign to positive for downward fluxes of positive momentum.
sign = 1.0 ; if (FIA%atmos_winds) sign = -1.0
I_avc = 1.0/real(FIA%avg_count)
!$OMP parallel do default(shared) private(u,v)
do j=jsc,jec
do k=0,ncat ; do i=isc,iec
u = FIA%flux_u_top(i,j,k) * (sign*I_avc)
v = FIA%flux_v_top(i,j,k) * (sign*I_avc)
FIA%flux_u_top(i,j,k) = u*G%cos_rot(i,j)-v*G%sin_rot(i,j) ! rotate stress from lat/lon
FIA%flux_v_top(i,j,k) = v*G%cos_rot(i,j)+u*G%sin_rot(i,j) ! to ocean coordinates
FIA%flux_sh_top(i,j,k) = FIA%flux_sh_top(i,j,k) * I_avc
FIA%evap_top(i,j,k) = FIA%evap_top(i,j,k) * I_avc
do b=1,nb ; FIA%flux_sw_top(i,j,k,b) = FIA%flux_sw_top(i,j,k,b) * I_avc ; enddo
FIA%flux_lw_top(i,j,k) = FIA%flux_lw_top(i,j,k) * I_avc
FIA%fprec_top(i,j,k) = FIA%fprec_top(i,j,k) * I_avc
FIA%lprec_top(i,j,k) = FIA%lprec_top(i,j,k) * I_avc
FIA%flux_lh_top(i,j,k) = FIA%flux_lh_top(i,j,k) * I_avc
! Copy radiation fields from the fast to the slow states.
if (k>0) FIA%sw_abs_ocn(i,j,k) = Rad%sw_abs_ocn(i,j,k)
! Convert frost forming atop sea ice into frozen precip.
! if ((k>0) .and. (FIA%evap_top(i,j,k) < 0.0)) then
! FIA%fprec_top(i,j,k) = FIA%fprec_top(i,j,k) - FIA%evap_top(i,j,k)
! FIA%evap_top(i,j,k) = 0.0
! endif
do n=1,FIA%num_tr_fluxes
FIA%tr_flux_top(i,j,k,n) = FIA%tr_flux_top(i,j,k,n) * I_avc
enddo
enddo ; enddo
do i=isc,iec
FIA%flux_sw_dn(i,j) = FIA%flux_sw_dn(i,j)*I_avc
FIA%Tskin_avg(i,j) = FIA%Tskin_avg(i,j) * I_avc
enddo
enddo
! Determine the fractional ice coverage and the wind stresses averaged
! across all the ice thickness categories on an A-grid.
FIA%WindStr_x(:,:) = 0.0 ; FIA%WindStr_y(:,:) = 0.0 ; FIA%ice_cover(:,:) = 0.0
!$OMP parallel do default(none) shared(isc,iec,jsc,jec,ncat,FIA,Rad,IST) &
!$OMP private(I_wts)
do j=jsc,jec
do k=1,ncat ; do i=isc,iec
FIA%WindStr_x(i,j) = FIA%WindStr_x(i,j) + IST%part_size(i,j,k) * FIA%flux_u_top(i,j,k)
FIA%WindStr_y(i,j) = FIA%WindStr_y(i,j) + IST%part_size(i,j,k) * FIA%flux_v_top(i,j,k)
FIA%ice_cover(i,j) = FIA%ice_cover(i,j) + IST%part_size(i,j,k)
enddo ; enddo
do i=isc,iec
FIA%WindStr_ocn_x(i,j) = FIA%flux_u_top(i,j,0)
FIA%WindStr_ocn_y(i,j) = FIA%flux_v_top(i,j,0)
if (FIA%ice_cover(i,j) > 0.0) then
I_wts = 1.0 / FIA%ice_cover(i,j)
FIA%WindStr_x(i,j) = FIA%WindStr_x(i,j) * I_wts
FIA%WindStr_y(i,j) = FIA%WindStr_y(i,j) * I_wts
if (FIA%ice_cover(i,j) > 1.0) FIA%ice_cover(i,j) = 1.0
! The max with 0 in the following line is here for safety; the only known
! instance where it has been required is when reading a SIS-1-derived
! restart file with tiny negative concentrations. SIS2 should not need it.
FIA%ice_free(i,j) = max(IST%part_size(i,j,0), 0.0)
! Rescale to add up to 1?
! I_wts = 1.0 / (FIA%ice_free(i,j) + FIA%ice_cover(i,j))
! FIA%ice_free(i,j) = FIA%ice_free(i,j) * I_wts
! FIA%ice_cover(i,j) = FIA%ice_cover(i,j) * I_wts
else
FIA%ice_free(i,j) = 1.0 ; FIA%ice_cover(i,j) = 0.0
FIA%WindStr_x(i,j) = FIA%flux_u_top(i,j,0)
FIA%WindStr_y(i,j) = FIA%flux_v_top(i,j,0)
endif
enddo
enddo
if (allocated(FIA%flux_sh0)) then
!$OMP parallel do default(shared)
do j=jsc,jec ; do k=0,ncat ; do i=isc,iec
FIA%dshdt(i,j,k) = FIA%dshdt(i,j,k) * I_avc
FIA%devapdt(i,j,k) = FIA%devapdt(i,j,k) * I_avc
FIA%dlwdt(i,j,k) = FIA%dlwdt(i,j,k) * I_avc
FIA%flux_sh0(i,j,k) = FIA%flux_sh0(i,j,k) * I_avc
FIA%evap0(i,j,k) = FIA%evap0(i,j,k) * I_avc
FIA%flux_lw0(i,j,k) = FIA%flux_lw0(i,j,k) * I_avc
FIA%Tskin_cat(i,j,k) = FIA%Tskin_cat(i,j,k) * I_avc
enddo ; enddo ; enddo
! Fill in the information to reconstruct the fluxes for any area-less categories.
! The open-ocean category must always be calculated for this to work properly.
call infill_array(IST, FIA%dshdt(:,:,0), FIA%dshdt(:,:,1:), G, IG)
call infill_array(IST, FIA%devapdt(:,:,0), FIA%devapdt(:,:,1:), G, IG)
call infill_array(IST, FIA%dlwdt(:,:,0), FIA%dlwdt(:,:,1:), G, IG)
call infill_array(IST, FIA%flux_sh0(:,:,0), FIA%flux_sh0(:,:,1:), G, IG)
call infill_array(IST, FIA%evap0(:,:,0), FIA%evap0(:,:,1:), G, IG)
call infill_array(IST, FIA%flux_lw0(:,:,0), FIA%flux_lw0(:,:,1:), G, IG)
call infill_array(IST, FIA%Tskin_cat(:,:,0), FIA%Tskin_cat(:,:,1:), G, IG)
endif
! set count to zero and fluxes will be zeroed before the next sum
FIA%avg_count = 0
end subroutine avg_top_quantities
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> total_top_quantities determines the sum across partitions of various fluxes
!! for later use on a potentially different ice state on the slow side.
subroutine total_top_quantities(FIA, TSF, part_size, G, IG)
type(fast_ice_avg_type), intent(in) :: FIA
type(total_sfc_flux_type), intent(inout) :: TSF
type(SIS_hor_grid_type), intent(inout) :: G
type(ice_grid_type), intent(in) :: IG
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce), &
intent(in) :: part_size
integer :: i, j, k, m, n, b, nb, isc, iec, jsc, jec, ncat
integer :: isd, ied, jsd, jed
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
nb = size(FIA%flux_sw_top,4)
if (TSF%num_tr_fluxes < 0) then
! Allocate the arrays to hold the tracer fluxes. This code is only exercised
! the first time that total_top_quantities is called.
TSF%num_tr_fluxes = FIA%num_tr_fluxes
if (TSF%num_tr_fluxes > 0) then
allocate(TSF%tr_flux(G%isd:G%ied, G%jsd:G%jed, TSF%num_tr_fluxes))
endif
endif
TSF%flux_u(:,:) = 0.0 ; TSF%flux_v(:,:) = 0.0
TSF%flux_sh(:,:) = 0.0 ; TSF%flux_lw(:,:) = 0.0 ; TSF%flux_lh(:,:) = 0.0
TSF%flux_sw(:,:,:) = 0.0
TSF%evap(:,:) = 0.0 ; TSF%fprec(:,:) = 0.0 ; TSF%lprec(:,:) = 0.0
if (TSF%num_tr_fluxes > 0) TSF%tr_flux(:,:,:) = 0.0
do k=0,ncat ; do j=jsc,jec ; do i=isc,iec
TSF%flux_u(i,j) = TSF%flux_u(i,j) + part_size(i,j,k) * FIA%flux_u_top(i,j,k)
TSF%flux_v(i,j) = TSF%flux_v(i,j) + part_size(i,j,k) * FIA%flux_v_top(i,j,k)
TSF%flux_sh(i,j) = TSF%flux_sh(i,j) + part_size(i,j,k) * FIA%flux_sh_top(i,j,k)
TSF%flux_lw(i,j) = TSF%flux_lw(i,j) + part_size(i,j,k) * FIA%flux_lw_top(i,j,k)
TSF%flux_lh(i,j) = TSF%flux_lh(i,j) + part_size(i,j,k) * FIA%flux_lh_top(i,j,k)
do b=1,nb
TSF%flux_sw(i,j,b) = TSF%flux_sw(i,j,b) + &
part_size(i,j,k) * FIA%flux_sw_top(i,j,k,b)
enddo
TSF%evap(i,j) = TSF%evap(i,j) + part_size(i,j,k) * FIA%evap_top(i,j,k)
TSF%fprec(i,j) = TSF%fprec(i,j) + part_size(i,j,k) * FIA%fprec_top(i,j,k)
TSF%lprec(i,j) = TSF%lprec(i,j) + part_size(i,j,k) * FIA%lprec_top(i,j,k)
do n=1,TSF%num_tr_fluxes
TSF%tr_flux(i,j,n) = TSF%tr_flux(i,j,n) + part_size(i,j,k) * FIA%tr_flux_top(i,j,k,n)
enddo
enddo ; enddo ; enddo
! If the sum of part_size across all the ice and ocean categories is not
! exactly 1, rescaling might be advisable, but for now it is assumed that
! part_size is properly scaled.
end subroutine total_top_quantities
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> find_excess_fluxes determines the difference between the sum across
!! partitions of various fluxes amd the sum previously found by total_top_quantities.
subroutine find_excess_fluxes(FIA, TSF, XSF, part_size, G, IG)
type(fast_ice_avg_type), intent(in) :: FIA
type(total_sfc_flux_type), intent(in) :: TSF
type(total_sfc_flux_type), intent(inout) :: XSF
type(SIS_hor_grid_type), intent(inout) :: G
type(ice_grid_type), intent(in) :: IG
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce), &
intent(in) :: part_size
integer :: i, j, k, m, n, b, nb, isc, iec, jsc, jec, ncat
integer :: isd, ied, jsd, jed
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
nb = size(FIA%flux_sw_top,4)
if (XSF%num_tr_fluxes < 0) then
! Allocate the arrays to hold the tracer fluxes. This code is only exercised
! the first time that total_top_quantities is called.
XSF%num_tr_fluxes = FIA%num_tr_fluxes
if (XSF%num_tr_fluxes > 0) then
allocate(XSF%tr_flux(G%isd:G%ied, G%jsd:G%jed, XSF%num_tr_fluxes))
endif
endif
! Note that XSF%flux_u and XSF%flux_v are not necessary as the changing ice cover is
! already being taken into account.
! This could be done more succinctly by initializing XSF = -TSF, but the
! answers would then be more accutely impacted by roundoff.
XSF%flux_sh(:,:) = 0.0 ; XSF%flux_lw(:,:) = 0.0 ; XSF%flux_lh(:,:) = 0.0
XSF%flux_sw(:,:,:) = 0.0
XSF%evap(:,:) = 0.0 ; XSF%fprec(:,:) = 0.0 ; XSF%lprec(:,:) = 0.0
if (XSF%num_tr_fluxes > 0) XSF%tr_flux(:,:,:) = 0.0
do k=0,ncat ; do j=jsc,jec ; do i=isc,iec
XSF%flux_sh(i,j) = XSF%flux_sh(i,j) + part_size(i,j,k) * FIA%flux_sh_top(i,j,k)
XSF%flux_lw(i,j) = XSF%flux_lw(i,j) + part_size(i,j,k) * FIA%flux_lw_top(i,j,k)
XSF%flux_lh(i,j) = XSF%flux_lh(i,j) + part_size(i,j,k) * FIA%flux_lh_top(i,j,k)
do b=1,nb
XSF%flux_sw(i,j,b) = XSF%flux_sw(i,j,b) + &
part_size(i,j,k) * FIA%flux_sw_top(i,j,k,b)
enddo
XSF%evap(i,j) = XSF%evap(i,j) + part_size(i,j,k) * FIA%evap_top(i,j,k)
XSF%fprec(i,j) = XSF%fprec(i,j) + part_size(i,j,k) * FIA%fprec_top(i,j,k)
XSF%lprec(i,j) = XSF%lprec(i,j) + part_size(i,j,k) * FIA%lprec_top(i,j,k)
do n=1,XSF%num_tr_fluxes
XSF%tr_flux(i,j,n) = XSF%tr_flux(i,j,n) + part_size(i,j,k) * FIA%tr_flux_top(i,j,k,n)
enddo
enddo ; enddo ; enddo
do j=jsc,jec ; do i=isc,iec
XSF%flux_sh(i,j) = XSF%flux_sh(i,j) - TSF%flux_sh(i,j)
XSF%flux_lw(i,j) = XSF%flux_lw(i,j) - TSF%flux_lw(i,j)
XSF%flux_lh(i,j) = XSF%flux_lh(i,j) - TSF%flux_lh(i,j)
do b=1,nb ; XSF%flux_sw(i,j,b) = XSF%flux_sw(i,j,b) - TSF%flux_sw(i,j,b) ; enddo
XSF%evap(i,j) = XSF%evap(i,j) - TSF%evap(i,j)
XSF%fprec(i,j) = XSF%fprec(i,j) - TSF%fprec(i,j)
XSF%lprec(i,j) = XSF%lprec(i,j) - TSF%lprec(i,j)
do n=1,XSF%num_tr_fluxes
XSF%tr_flux(i,j,n) = XSF%tr_flux(i,j,n) - TSF%tr_flux(i,j,n)
enddo
enddo ; enddo
end subroutine find_excess_fluxes
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> rescale_shortwave determines whether any shortwave frequency bands exceed their
!! intensity during the atmospheric steps and if so scales them back for energy
!! conservation. If the shortwave heating in any bands have decreased, the
!! difference will later be applied to the ocean.
subroutine rescale_shortwave(FIA, TSF, part_size, G, IG)
type(fast_ice_avg_type), intent(inout) :: FIA
type(total_sfc_flux_type), intent(in) :: TSF
type(SIS_hor_grid_type), intent(inout) :: G
type(ice_grid_type), intent(in) :: IG
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce), &
intent(in) :: part_size
real, dimension(G%isd:G%ied,size(FIA%flux_sw_top,4)) :: &
sw_tot_band
real :: rescale ! A rescaling factor between 0 and 1.
integer :: i, j, k, m, n, b, nb, nfb, isc, iec, jsc, jec, ncat
integer :: isd, ied, jsd, jed
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
nb = size(FIA%flux_sw_top,4) ; nfb = nb/2
do j=jsc,jec
sw_tot_band(:,:) = 0.0
do k=0,ncat ; do b=1,nb ; do i=isc,iec
sw_tot_band(i,b) = sw_tot_band(i,b) + &
part_size(i,j,k) * FIA%flux_sw_top(i,j,k,b)
enddo ; enddo ; enddo
do i=isc,iec ; do b=1,nb-1,2
! Ice can scatter direct shortwave into diffuse without loss of energy
! conservation, so it only the total of the shortwave in each frequency
! band that needs to be considered. If there are more than 2 angular
! bands (direct and diffuse), this code will need to be modified.
if (abs(sw_tot_band(i,b) + sw_tot_band(i,b+1)) > &
abs(TSF%flux_sw(i,j,b) + TSF%flux_sw(i,j,b+1))) then
rescale = abs(TSF%flux_sw(i,j,b) + TSF%flux_sw(i,j,b+1)) / &
abs(sw_tot_band(i,b) + sw_tot_band(i,b+1))
do k=0,ncat
FIA%flux_sw_top(i,j,k,b) = rescale * FIA%flux_sw_top(i,j,k,b)
FIA%flux_sw_top(i,j,k,b+1) = rescale * FIA%flux_sw_top(i,j,k,b+1)
enddo
endif
enddo ; enddo
enddo
end subroutine rescale_shortwave
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> infill_array fills in an array with actual, interpolated or plausible values
!! from other ice categories.
subroutine infill_array(IST, ta_ocn, ta, G, IG)
type(ice_state_type), intent(in) :: IST !< The ice state type whose values of part_size and
!! mH_ice are being used to control the infilling.
type(SIS_hor_grid_type), intent(in) :: G !< The sea-ice lateral grid type.
type(ice_grid_type), intent(in) :: IG !< The sea-ice grid type.
real, dimension(G%isd:G%ied, G%jsd:G%jed), &
intent(in) :: ta_ocn !< The value of the array for ocean partitions.
real, dimension(G%isd:G%ied, G%jsd:G%jed, IG%CatIce), &
intent(inout) :: ta !< The array that is being infilled.
integer, dimension(G%isd:G%ied) :: k_thick, k_thin
real, dimension(G%isd:G%ied) :: &
ta_thick, ta_thin, mH_thin
real :: wt2, mH_cat_tgt
logical, dimension(G%isd:G%ied) :: infill
logical :: any_ice
integer :: i, j, k, k2, k3, isc, iec, jsc, jec, ncat, NkIce
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec
ncat = IG%CatIce ; NkIce = IG%NkIce
do j=jsc,jec
! Determine which categories exist.
do i=isc,iec ; k_thick(i) = -1 ; k_thin(i) = ncat+1 ; mH_thin(i) = 0.0 ; enddo
any_ice = .false.
do k=1,ncat ; do i=isc,iec ; if (IST%part_size(i,j,k) > 0.0) then
k_thick(i) = k
ta_thick(i) = ta(i,j,k)
any_ice = .true.
endif ; enddo ; enddo
do i=isc,iec
if (k_thick(i) < 1) ta_thick(i) = G%mask2dT(i,j) * ta_ocn(i,j)
enddo
if (.not.any_ice) then
! This entire row is ice free.
do k=1,ncat ; do i=isc,iec
ta(i,j,k) = G%mask2dT(i,j)*ta_ocn(i,j)
enddo ; enddo
else
do k=ncat,1,-1 ; do i=isc,iec ; if (IST%part_size(i,j,k) > 0.0) then
k_thin(i) = k
ta_thin(i) = ta(i,j,k)
mH_thin(i) = IST%mH_ice(i,j,k)
endif ; enddo ; enddo
! The logic for the k=1 (thinest) category is particularly simple.
do i=isc,iec ; if (IST%part_size(i,j,1) <= 0.0) then
ta(i,j,1) = G%mask2dT(i,j)*ta_ocn(i,j)
endif ; enddo
do k=2,ncat ; do i=isc,iec
if (G%mask2dT(i,j) > 0.0) then
if (k > k_thick(i)) then
! This is the most common case, since it applies to ice-free water.
ta(i,j,k) = ta_thick(i)
elseif (IST%part_size(i,j,k) > 0.0) then
! This is a valid value - no interpolation is needed.
cycle
elseif (k < k_thin(i)) then
! Linearly interpolate to the ocean's freezing temperature.
mH_cat_tgt = 0.5*(IG%mH_cat_bound(K) + IG%mH_cat_bound(K+1))
wt2 = 0.0 ! If in doubt, use ta_thin. This should not happen.
if (mH_thin(i) > mH_cat_tgt) wt2 = (mH_thin(i) - mH_cat_tgt) / mH_thin(i)
ta(i,j,k) = wt2*ta_ocn(i,j) + (1.0-wt2)*ta_thin(i)
else
! Linearly interpolate between bracketing neighboring
! categories with valid values.
do k2=k+1,ncat ; if (IST%part_size(i,j,k2) > 0.0) exit ; enddo
do k3=k-1,1,-1 ; if (IST%part_size(i,j,k3) > 0.0) exit ; enddo
mH_cat_tgt = 0.5*(IG%mH_cat_bound(K) + IG%mH_cat_bound(K+1))
wt2 = 0.5
if ((IST%mH_ice(i,j,k3) > mH_cat_tgt) .and. &
(IST%mH_ice(i,j,k2) < mH_cat_tgt)) &
wt2 = (IST%mH_ice(i,j,k3) - mH_cat_tgt) / &
(IST%mH_ice(i,j,k3) - IST%mH_ice(i,j,k2))
ta(i,j,k) = wt2*ta(i,j,k2) + (1.0-wt2) * ta(i,j,k3)
endif
else ! This is a land point.
ta(i,j,k) = 0.0
endif
enddo ; enddo
endif ! .not.any_ice
enddo
end subroutine infill_array
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> do_update_ice_model_fast applies the surface heat fluxes, shortwave radiation
!! diffusion of heat to the sea-ice to implicitly determine a new temperature
!! profile, subject to the constraint that ice and snow temperatures are never
!! above freezing. Melting and freezing occur elsewhere.
subroutine do_update_ice_model_fast(Atmos_boundary, IST, sOSS, Rad, FIA, &
Time_step, CS, G, IG )
type(atmos_ice_boundary_type), intent(in) :: Atmos_boundary
type(ice_state_type), intent(inout) :: IST
type(simple_OSS_type), intent(in) :: sOSS
type(ice_rad_type), intent(inout) :: Rad
type(fast_ice_avg_type), intent(inout) :: FIA
type(time_type), intent(in) :: Time_step ! The amount of time over which to advance the ice.
type(fast_thermo_CS), pointer :: CS
type(SIS_hor_grid_type), intent(inout) :: G
type(ice_grid_type), intent(in) :: IG
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce) :: &
flux_sh, & ! The upward sensible heat flux from the ice to the atmosphere
! at the surface of the ice, in W m-2.
evap, & ! The upward flux of water due to sublimation or evaporation
! from the top of the ice to the atmosphere, in kg m-2 s-1.
flux_lh, & ! The upward latent heat flux associated with sublimation or
! evaporation, in W m-2.
flux_lw, & ! The net downward longwave heat flux into the ice, in W m-2.
flux_u, flux_v, lprec, fprec, &
sh_T0, & ! The upward sensible heat flux from the top of the ice into
! the atmosphere when the skin temperature is 0 C, in W m-2.
evap_T0, & ! The sublimation rate when the skin temperature is 0 C,
! in kg m-2 s-1.
lw_T0, & ! The downward longwave heat flux from the atmosphere into the
! ice or ocean when the skin temperature is 0 C, in W m-2.
dshdt, & ! The derivative of the upward sensible heat flux with the surface
! temperature in W m-2 K-1.
devapdt, & ! The derivative of the sublimation rate with the surface
! temperature, in kg m-2 s-1 K-1.
dlwdt ! The derivative of the downward radiative heat flux with surface
! temperature (i.e. d(flux_lw)/d(surf_temp)) in W m-2 K-1.
real, dimension(G%isd:G%ied,G%jsd:G%jed,0:IG%CatIce,size(FIA%flux_sw_top,4)) :: &
flux_sw ! The downward shortwave heat fluxes in W m-2. The fourth
! dimension is a combination of angular orientation and frequency.
real, dimension(0:IG%NkIce) :: T_col ! The temperature of a column of ice and snow in degC.
real, dimension(IG%NkIce) :: S_col ! The thermodynamic salinity of a column of ice, in g/kg.
real, dimension(0:IG%NkIce) :: enth_col ! The enthalpy of a column of snow and ice, in enth_unit (J/kg?).
real, dimension(0:IG%NkIce) :: SW_abs_col
real :: dt_fast, ts_new, dts, hf, hfd, latent
real :: hf_0 ! The positive upward surface heat flux when T_sfc = 0 C, in W m-2.
real :: dhf_dt ! The deriviative of the upward surface heat flux with Ts, in W m-2 C-1.
real :: sw_tot ! sum over all shortwave (dir/dif and vis/nir) components
real :: LatHtFus ! The latent heat of fusion of ice in J/kg.
real :: LatHtVap ! The latent heat of vaporization of water at 0C in J/kg.
real :: H_to_m_ice ! The specific volumes of ice and snow times the
real :: H_to_m_snow ! conversion factor from thickness units, in m H-1.
logical :: slab_ice ! If true, use the very old slab ice thermodynamics,
! with effectively zero heat capacity of ice and snow.
type(time_type) :: Dt_ice
logical :: sent
integer :: i, j, k, m, i2, j2, k2, isc, iec, jsc, jec, ncat, i_off, j_off, NkIce
real :: tot_heat_in, enth_here, enth_imb, norm_enth_imb, SW_absorbed
real :: enth_liq_0 ! The value of enthalpy for liquid fresh water at 0 C, in
! enthalpy units (sometimes J kg-1).
real :: enth_units ! A conversion factor from Joules kg-1 to enthalpy units.
real :: I_enth_unit ! The inverse of enth_units, in J kg-1 enth_unit-1.
real :: I_Nk
real :: kg_H_Nk ! The conversion factor from units of H to kg/m2 over Nk.
if (.not.associated(CS)) call SIS_error(FATAL, &
"SIS_fast_thermo: Module must be initialized before it is used.")
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
i_off = LBOUND(Atmos_boundary%t_flux,1) - G%isc
j_off = LBOUND(Atmos_boundary%t_flux,2) - G%jsc
NkIce = IG%NkIce ; I_Nk = 1.0 / NkIce ; kg_H_Nk = IG%H_to_kg_m2 * I_Nk
CS%n_fast = CS%n_fast + 1
if (CS%column_check .and. (FIA%avg_count==0)) then
CS%heat_in(:,:,:) = 0.0
CS%enth_prev(:,:,:) = 0.0
do k=1,ncat ; do j=jsc,jec ; do i=isc,iec ; if (IST%mH_ice(i,j,k)>0.0) then
CS%enth_prev(i,j,k) = (IST%mH_snow(i,j,k)*IG%H_to_kg_m2) * IST%enth_snow(i,j,k,1)
do m=1,IG%NkIce
CS%enth_prev(i,j,k) = CS%enth_prev(i,j,k) + &
(IST%mH_ice(i,j,k)*kg_H_Nk) * IST%enth_ice(i,j,k,m)
enddo
endif ; enddo ; enddo ; enddo
endif
if (CS%debug) &
call IST_chksum("Start do_update_ice_model_fast", IST, G, IG)
!$OMP parallel do default(shared) private(i2,j2,k2)
do j=jsc,jec
! Set up local copies of fluxes. The Atmos_boundary arrays may have
! different index conventions than are used internally in this component.
do k=0,ncat ; do i=isc,iec
i2 = i+i_off ; j2 = j+j_off ; k2 = k+1
flux_u(i,j,k) = Atmos_boundary%u_flux(i2,j2,k2)
flux_v(i,j,k) = Atmos_boundary%v_flux(i2,j2,k2)
flux_sh(i,j,k) = Atmos_boundary%t_flux(i2,j2,k2)
evap(i,j,k) = Atmos_boundary%q_flux(i2,j2,k2)
flux_lw(i,j,k) = Atmos_boundary%lw_flux(i2,j2,k2)
flux_sw(i,j,k,nir_dir) = Atmos_boundary%sw_flux_nir_dir(i2,j2,k2)
flux_sw(i,j,k,nir_dif) = Atmos_boundary%sw_flux_nir_dif(i2,j2,k2)
flux_sw(i,j,k,vis_dir) = Atmos_boundary%sw_flux_vis_dir(i2,j2,k2)
flux_sw(i,j,k,vis_dif) = Atmos_boundary%sw_flux_vis_dif(i2,j2,k2)
lprec(i,j,k) = Atmos_boundary%lprec(i2,j2,k2)
fprec(i,j,k) = Atmos_boundary%fprec(i2,j2,k2)
dshdt(i,j,k) = Atmos_boundary%dhdt(i2,j2,k2)
devapdt(i,j,k) = Atmos_boundary%dedt(i2,j2,k2)
dlwdt(i,j,k) = -1.*Atmos_boundary%drdt(i2,j2,k2)
enddo ; enddo
do i=isc,iec
sh_T0(i,j,0) = flux_sh(i,j,0) - dshdt(i,j,0) * sOSS%SST_C(i,j)
evap_T0(i,j,0) = evap(i,j,0) - devapdt(i,j,0) * sOSS%SST_C(i,j)
lw_T0(i,j,0) = flux_lw(i,j,0) - dlwdt(i,j,0) * sOSS%SST_C(i,j)
enddo
do k=1,ncat ; do i=isc,iec
sh_T0(i,j,k) = flux_sh(i,j,k) - dshdt(i,j,k) * Rad%t_skin(i,j,k)
evap_T0(i,j,k) = evap(i,j,k) - devapdt(i,j,k) * Rad%t_skin(i,j,k)
lw_T0(i,j,k) = flux_lw(i,j,k) - dlwdt(i,j,k) * Rad%t_skin(i,j,k)
enddo ; enddo
enddo
if (CS%debug) then
call hchksum(flux_u(:,:,1:), "Mid do_fast flux_u", G%HI)
call hchksum(flux_v(:,:,1:), "Mid do_fast flux_v", G%HI)
call hchksum(flux_sh(:,:,1:), "Mid do_fast flux_sh", G%HI)
call hchksum(evap(:,:,1:), "Mid do_fast evap", G%HI)
call hchksum(flux_lw(:,:,1:), "Mid do_fast flux_lw", G%HI)
call hchksum(flux_sw(:,:,1:,nir_dir), "Mid do_fast flux_sw_nir_dir", G%HI)
call hchksum(flux_sw(:,:,1:,nir_dif), "Mid do_fast flux_sw_nir_dif", G%HI)
call hchksum(flux_sw(:,:,1:,vis_dir), "Mid do_fast flux_sw_vis_dir", G%HI)
call hchksum(flux_sw(:,:,1:,vis_dif), "Mid do_fast flux_sw_vis_dif", G%HI)
call hchksum(lprec(:,:,1:), "Mid do_fast lprec", G%HI)
call hchksum(fprec(:,:,1:), "Mid do_fast fprec", G%HI)
call hchksum(dshdt(:,:,1:), "Mid do_fast dshdt", G%HI)
call hchksum(devapdt(:,:,1:), "Mid do_fast devapdt", G%HI)
call hchksum(dlwdt(:,:,1:), "Mid do_fast dlwdt", G%HI)
endif
call get_SIS2_thermo_coefs(IST%ITV, ice_salinity=S_col, enthalpy_units=enth_units, &
Latent_fusion=LatHtFus, Latent_vapor=LatHtVap, slab_ice=slab_ice)
do j=jsc,jec ; do i=isc,iec
flux_lh(i,j,0) = LatHtVap * evap(i,j,0)
enddo ; enddo
!
! implicit update of ice surface temperature
!
dt_fast = time_type_to_real(Time_step)
enth_liq_0 = Enth_from_TS(0.0, 0.0, IST%ITV) ; I_enth_unit = 1.0 / enth_units
!$OMP parallel do default(none) shared(isc,iec,jsc,jec,ncat,NkIce,IST,dshdt,devapdt,dlwdt, &
!$OMP flux_sw,flux_sh,evap,flux_lw,enth_liq_0,&
!$OMP dt_fast,flux_lh,I_enth_unit,G,S_col,kg_H_Nk,slab_ice,&
!$OMP enth_units,LatHtFus,LatHtVap,IG,sOSS,FIA,Rad,CS) &
!$OMP private(latent,enth_col,sw_tot,dhf_dt, &
!$OMP hf_0,ts_new,dts,SW_abs_col,SW_absorbed,enth_here,&
!$OMP tot_heat_in,enth_imb,norm_enth_imb )
do j=jsc,jec ; do k=1,ncat ; do i=isc,iec
if (IST%part_size(i,j,k) > 0.0) then
enth_col(0) = IST%enth_snow(i,j,k,1)
do m=1,NkIce ; enth_col(m) = IST%enth_ice(i,j,k,m) ; enddo
! In the case of sublimation of either snow or ice, the vapor is at 0 C.
! If the vapor should be at a different temperature, a correction would be
! made here.
if (slab_ice) then
latent = LatHtVap + LatHtFus
elseif (IST%mH_snow(i,j,k)>0.0) then
latent = LatHtVap + (enth_liq_0 - IST%enth_snow(i,j,k,1)) * I_enth_unit
else
latent = LatHtVap + (enth_liq_0 - IST%enth_ice(i,j,k,1)) * I_enth_unit
endif
sw_tot = (flux_sw(i,j,k,vis_dir) + flux_sw(i,j,k,vis_dif)) + &
(flux_sw(i,j,k,nir_dir) + flux_sw(i,j,k,nir_dif))
dhf_dt = (dshdt(i,j,k) + devapdt(i,j,k)*latent) - dlwdt(i,j,k)
if (CS%Reorder_0C_heatflux) then
! This form seperately projects each contribution to the surface heat
! flux back to its value when T=0, so that a bitwise identical result
! can be obtained if the heating is redone. The two forms are
! mathematically equivalent.
hf_0 = ((flux_sh(i,j,k) - dshdt(i,j,k) * Rad%t_skin(i,j,k)) + &
(evap(i,j,k) - devapdt(i,j,k) * Rad%t_skin(i,j,k)) * latent) - &
((flux_lw(i,j,k) - dlwdt(i,j,k) * Rad%t_skin(i,j,k)) + &
Rad%sw_abs_sfc(i,j,k) * sw_tot)
else
hf_0 = ((flux_sh(i,j,k) + evap(i,j,k)*latent) - &
(flux_lw(i,j,k) + Rad%sw_abs_sfc(i,j,k)* sw_tot)) - &
dhf_dt * Rad%t_skin(i,j,k)
endif
SW_abs_col(0) = Rad%sw_abs_snow(i,j,k)*sw_tot
do m=1,NkIce ; SW_abs_col(m) = Rad%sw_abs_ice(i,j,k,m)*sw_tot ; enddo
! This call updates the snow and ice temperatures and accumulates the
! surface and bottom melting/freezing energy. The ice and snow do not
! actually lose or gain any mass from freezing or melting.
! mw/new - pass melt pond (surface temp fixed at freezing when present)
call ice_temp_SIS2(IST%mH_pond(i,j,k)*IG%H_to_kg_m2, &
IST%mH_snow(i,j,k)*IG%H_to_kg_m2, &
IST%mH_ice(i,j,k)*IG%H_to_kg_m2, &
enth_col, S_col, hf_0, dhf_dt, SW_abs_col, &
sOSS%T_fr_ocn(i,j), FIA%bheat(i,j), ts_new, &
dt_fast, NkIce, FIA%tmelt(i,j,k), FIA%bmelt(i,j,k), &
CS%ice_thm_CSp, IST%ITV, CS%column_check)
IST%enth_snow(i,j,k,1) = enth_col(0)
do m=1,NkIce ; IST%enth_ice(i,j,k,m) = enth_col(m) ; enddo
if (CS%Reorder_0C_heatflux) then
! These extended expressions for the new fluxes will reproduce answers
! with redo_update_ice_model_fast. They are mathematically equivalent
! to the next set of expressions.
flux_sh(i,j,k) = (flux_sh(i,j,k) - dshdt(i,j,k) * Rad%t_skin(i,j,k)) + ts_new * dshdt(i,j,k)
evap(i,j,k) = (evap(i,j,k) - devapdt(i,j,k) * Rad%t_skin(i,j,k)) + ts_new * devapdt(i,j,k)
flux_lw(i,j,k) = (flux_lw(i,j,k) - dlwdt(i,j,k) * Rad%t_skin(i,j,k)) + ts_new * dlwdt(i,j,k)
else
dts = ts_new - Rad%t_skin(i,j,k)
flux_sh(i,j,k) = flux_sh(i,j,k) + dts * dshdt(i,j,k)
evap(i,j,k) = evap(i,j,k) + dts * devapdt(i,j,k)
flux_lw(i,j,k) = flux_lw(i,j,k) + dts * dlwdt(i,j,k)
endif
flux_lh(i,j,k) = latent * evap(i,j,k)
Rad%t_skin(i,j,k) = ts_new
if (CS%column_check) then
SW_absorbed = SW_abs_col(0)
do m=1,NkIce ; SW_absorbed = SW_absorbed + SW_abs_col(m) ; enddo
CS%heat_in(i,j,k) = CS%heat_in(i,j,k) + dt_fast * &
((flux_lw(i,j,k) + Rad%sw_abs_sfc(i,j,k)*sw_tot) + SW_absorbed + &
FIA%bheat(i,j) - (flux_sh(i,j,k) + flux_lh(i,j,k)))
enth_here = (IG%H_to_kg_m2*IST%mH_snow(i,j,k)) * enth_col(0)
do m=1,NkIce
enth_here = enth_here + (IST%mH_ice(i,j,k)*kg_H_Nk) * enth_col(m)
enddo
tot_heat_in = enth_units * (CS%heat_in(i,j,k) - &
(FIA%bmelt(i,j,k) + FIA%tmelt(i,j,k)))
enth_imb = enth_here - (CS%enth_prev(i,j,k) + tot_heat_in)
if (abs(enth_imb) > CS%imb_tol * (abs(enth_here) + &
abs(CS%enth_prev(i,j,k)) + abs(tot_heat_in)) ) then
norm_enth_imb = enth_imb / (abs(enth_here) + &
abs(CS%enth_prev(i,j,k)) + abs(tot_heat_in))
enth_imb = enth_here - (CS%enth_prev(i,j,k) + tot_heat_in)
endif
endif
else ! IST%mH_ice <= 0
flux_lh(i,j,k) = LatHtVap * evap(i,j,k)
endif
enddo ; enddo ; enddo
call sum_top_quantities(FIA, Atmos_boundary, flux_u, flux_v, flux_sh, evap, &
flux_sw, flux_lw, lprec, fprec, flux_lh, Rad%t_skin, sOSS%SST_C, &
sh_T0, evap_T0, lw_T0, dshdt, devapdt, dlwdt, G, IG )
if (CS%debug) &
call IST_chksum("End do_update_ice_model_fast", IST, G, IG)
if (CS%bounds_check) &
call IST_bounds_check(IST, G, IG, "End of update_ice_fast", Rad=Rad) !, OSS=sOSS)
end subroutine do_update_ice_model_fast
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> accumulate_deposition_fluxes accumulates any fluxes that are labeled as having
!! the type 'air_sea_deposition'. These fluxes were not available when
!! update_ice_model_fast was called and the other fluxes were accumulated.
subroutine accumulate_deposition_fluxes(ABT, FIA, G, IG)
type(atmos_ice_boundary_type), intent(in) :: ABT
type(fast_ice_avg_type), intent(inout) :: FIA
type(SIS_hor_grid_type), intent(in) :: G
type(ice_grid_type), intent(in) :: IG
integer :: i, j, k, m, n, i_off, j_off, i2, j2, k2, isc, iec, jsc, jec, ncat, ind
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
i_off = LBOUND(ABT%t_flux,1) - G%isc
j_off = LBOUND(ABT%t_flux,2) - G%jsc
ind = 0
do n=1,ABT%fluxes%num_bcs ; do m=1,ABT%fluxes%bc(n)%num_fields
ind = ind + 1
if (ABT%fluxes%bc(n)%flux_type == 'air_sea_deposition') then
do k=0,ncat ; do j=jsc,jec ; do i=isc,iec
i2 = i+i_off ; j2 = j+j_off ; k2 = k+1
FIA%tr_flux_top(i,j,k,ind) = FIA%tr_flux_top(i,j,k,ind) + &
ABT%fluxes%bc(n)%field(m)%values(i2,j2,k2)
enddo ; enddo ; enddo
endif
enddo; enddo
end subroutine accumulate_deposition_fluxes
!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
!> redo_update_ice_model_fast applies the surface heat fluxes, shortwave radiation
!! diffusion of heat to the sea-ice to implicitly determine a new temperature
!! profile, subject to the constraint that ice and snow temperatures are never
!! above freezing, using fluxes that have been determined during previous calls
!! to do_update_ice_model_fast and stored in the fast_ice_avg_type.
subroutine redo_update_ice_model_fast(IST, sOSS, Rad, FIA, &
Time_step, CS, G, IG )
type(ice_state_type), intent(inout) :: IST
type(simple_OSS_type), intent(in) :: sOSS
type(ice_rad_type), intent(inout) :: Rad
type(fast_ice_avg_type), intent(inout) :: FIA
type(time_type), intent(in) :: Time_step ! The amount of time over which to advance the ice.
type(fast_thermo_CS), pointer :: CS
type(SIS_hor_grid_type), intent(inout) :: G
type(ice_grid_type), intent(in) :: IG
real, dimension(0:IG%NkIce) :: T_col ! The temperature of a column of ice and snow in degC.
real, dimension(IG%NkIce) :: S_col ! The thermodynamic salinity of a column of ice, in g/kg.
real, dimension(0:IG%NkIce) :: enth_col ! The enthalpy of a column of snow and ice, in enth_unit (J/kg?).
real, dimension(0:IG%NkIce) :: SW_abs_col
real :: dt_here, ts_new, dts, hf, hfd, latent
real :: hf_0 ! The positive upward surface heat flux when T_sfc = 0 C, in W m-2.
real :: dhf_dt ! The deriviative of the upward surface heat flux with Ts, in W m-2 C-1.
real :: sw_tot ! sum over dir/dif vis/nir components
real :: LatHtFus ! The latent heat of fusion of ice in J/kg.
real :: LatHtVap ! The latent heat of vaporization of water at 0C in J/kg.
real :: H_to_m_ice ! The specific volumes of ice and snow times the
real :: H_to_m_snow ! conversion factor from thickness units, in m H-1.
logical :: slab_ice ! If true, use the very old slab ice thermodynamics,
! with effectively zero heat capacity of ice and snow.
type(time_type) :: Dt_ice
logical :: sent
integer :: i, j, k, m, i2, j2, k2, isc, iec, jsc, jec, ncat, i_off, j_off, NkIce
real :: tot_heat_in, enth_here, enth_imb, norm_enth_imb, SW_absorbed
real :: enth_liq_0 ! The value of enthalpy for liquid fresh water at 0 C, in
! enthalpy units (sometimes J kg-1).
real :: enth_units ! A conversion factor from Joules kg-1 to enthalpy units.
real :: I_enth_unit ! The inverse of enth_units, in J kg-1 enth_unit-1.
real :: I_Nk
real :: kg_H_Nk ! The conversion factor from units of H to kg/m2 over Nk.
if (.not.associated(CS)) call SIS_error(FATAL, &
"SIS_fast_thermo: Module must be initialized before it is used.")
isc = G%isc ; iec = G%iec ; jsc = G%jsc ; jec = G%jec ; ncat = IG%CatIce
NkIce = IG%NkIce ; I_Nk = 1.0 / NkIce ; kg_H_Nk = IG%H_to_kg_m2 * I_Nk
if (CS%debug) &
call IST_chksum("Start redo_update_ice_model_fast", IST, G, IG)
call get_SIS2_thermo_coefs(IST%ITV, ice_salinity=S_col, enthalpy_units=enth_units, &
Latent_fusion=LatHtFus, Latent_vapor=LatHtVap, slab_ice=slab_ice)
!
! implicit update of ice surface temperature
!
dt_here = time_type_to_real(Time_step)
enth_liq_0 = Enth_from_TS(0.0, 0.0, IST%ITV) ; I_enth_unit = 1.0 / enth_units
!### RESCALE THE RADIATION
!### IS Rad A VALID TYPE TO USE?
!### DO WE NEED TO UPDATE RAD%T_SKIN? (PROBABLY ONLY FOR DIAGNOSTICS?)